Jpp in_tag_pdf_generation
the software that should make you happy
Loading...
Searching...
No Matches
JEvD.cc File Reference

Program to display hit probabilities. More...

#include <string>
#include <iostream>
#include <iomanip>
#include <vector>
#include <algorithm>
#include <memory>
#include "TROOT.h"
#include "TApplication.h"
#include "TCanvas.h"
#include "TRootCanvas.h"
#include "TStyle.h"
#include "TH2D.h"
#include "TArrow.h"
#include "TLatex.h"
#include "TMarker.h"
#include "km3net-dataformat/offline/Head.hh"
#include "km3net-dataformat/offline/Evt.hh"
#include "km3net-dataformat/online/JDAQ.hh"
#include "km3net-dataformat/online/JDAQHeader.hh"
#include "km3net-dataformat/tools/time_converter.hh"
#include "JROOT/JStyle.hh"
#include "JROOT/JCanvas.hh"
#include "JROOT/JRootToolkit.hh"
#include "JDAQ/JDAQSummarysliceIO.hh"
#include "JTrigger/JHitL0.hh"
#include "JSupport/JSupport.hh"
#include "JSupport/JSingleFileScanner.hh"
#include "JSupport/JSummaryFileRouter.hh"
#include "JSupport/JMonteCarloFileSupportkit.hh"
#include "JAAnet/JAAnetToolkit.hh"
#include "JPhysics/JPDF_t.hh"
#include "JFit/JLine1Z.hh"
#include "JFit/JModel.hh"
#include "JReconstruction/JHitW0.hh"
#include "JReconstruction/JMuonGandalfParameters_t.hh"
#include "JLang/JPredicate.hh"
#include "JLang/JComparator.hh"
#include "JMath/JMathToolkit.hh"
#include "JSystem/JKeypress.hh"
#include "JSystem/JProcess.hh"
#include "Jeep/JFunctionAdaptor.hh"
#include "Jeep/JProperties.hh"
#include "Jeep/JPrint.hh"
#include "Jeep/JParser.hh"
#include "Jeep/JMessage.hh"

Go to the source code of this file.

Classes

struct  JAANET::JEventSelector
 Event selector. More...
 

Namespaces

namespace  JAANET
 Extensions to Evt data format.
 

Functions

bool JAANET::hasW (const Trk &trk, const int i)
 Check availability of value.
 
double JAANET::getW (const Trk &trk, const int i)
 Get associated value.
 
double JAANET::getW (const Trk &trk, const int i, const double value)
 Get associated value.
 
void JAANET::setW (Trk &trk, const int i, const double value)
 Set associated value.
 
int main (int argc, char **argv)
 

Detailed Description

Program to display hit probabilities.

Author
mdejong

Definition in file JAAnet/JEvD.cc.

Function Documentation

◆ main()

int main ( int argc,
char ** argv )

Definition at line 195 of file JAAnet/JEvD.cc.

196{
197 using namespace std;
198 using namespace JPP;
199 using namespace KM3NETDAQ;
200
201 JSingleFileScanner<Evt> inputFile;
202 JLimit_t& numberOfEvents = inputFile.getLimit();
203 string pdfFile;
204 string outputFile;
205 JMuonGandalfParameters_t parameters;
206 int application;
207 JEventSelector event_selector;
208 JCanvas canvas;
209 bool batch;
210 struct : JStyle::JParameters {
211 double arrowSize = 0.003;
212 string arrowType = "|->";
213 double arrowScale = 250.0;
214 Width_t lineWidth = 2;
215 Style_t lineStyle = 1;
216 int nbinsX = 50;
217 int nbinsY = 250;
218 double T_ns = 0.0;
219 } graphics;
220 string option;
221 int debug;
222
223
224 try {
225
226 JProperties properties = graphics.getProperties();
227
228 properties.insert(gmake_property(graphics.arrowSize));
229 properties.insert(gmake_property(graphics.arrowType));
230 properties.insert(gmake_property(graphics.arrowScale));
231 properties.insert(gmake_property(graphics.lineWidth));
232 properties.insert(gmake_property(graphics.lineStyle));
233 properties.insert(gmake_property(graphics.T_ns));
234
235 parameters.numberOfPrefits = 1;
236
237 JParser<> zap("Program to display hit probabilities.");
238
239 zap['w'] = make_field(canvas, "size of canvas <nx>x<ny> [pixels]") = JCanvas(1200, 600);
240 zap['f'] = make_field(inputFile, "input file (output of JXXXMuonReconstruction.sh)");
241 zap['n'] = make_field(numberOfEvents) = JLimit::max();
242 zap['P'] = make_field(pdfFile);
243 zap['o'] = make_field(outputFile, "graphics output file name") = MAKE_STRING("display_" << WILDCARD << ".gif");
244 zap['@'] = make_field(parameters) = JPARSER::initialised();
246 zap['L'] = make_field(event_selector) = JPARSER::initialised();
247 zap['%'] = make_field(properties) = JPARSER::initialised();
248 zap['O'] = make_field(option, "draw option") = arrow_t, histogram_t;
249 zap['B'] = make_field(batch, "batch processing");
250 zap['d'] = make_field(debug) = 1;
251
252 zap(argc, argv);
253 }
254 catch(const exception& error) {
255 FATAL(error.what() << endl);
256 }
257
258 if (batch && outputFile == "") {
259 FATAL("Missing output file name " << outputFile << " in batch mode." << endl);
260 }
261
262 if (!batch && outputFile == "") {
263 outputFile = MAKE_STRING(WILDCARD << ".gif");
264 }
265
266 if (outputFile.find(WILDCARD) == string::npos) {
267 FATAL("Output file name " << outputFile << " has no wild card '" << WILDCARD << "'" << endl);
268 }
269
270 JHit::setSlewing(false);
271
272 JSummaryFileRouter summary(inputFile);
273
274 const JMuonPDF_t pdf(pdfFile, parameters.TTS_ns);
275
276 const JTimeRange T_ns(parameters.TMin_ns, parameters.TMax_ns);
277
278 typedef vector<JHitL0> JDataL0_t;
279 typedef vector<JHitW0> JDataW0_t;
280
281
282 // ROOT
283
284 gROOT->SetBatch(batch);
285
286 TApplication* tp = new TApplication("user", NULL, NULL);
287 TCanvas* cv = new TCanvas("display", "", canvas.x, canvas.y);
288
289 if (!batch) {
290 ((TRootCanvas *) cv->GetCanvasImp())->Connect("CloseWindow()", "TApplication", tp, "Terminate()");
291 }
292
293 unique_ptr<TStyle> gStyle(new JStyle("gplot", cv->GetWw(), cv->GetWh(), graphics));
294
295 gROOT->SetStyle("gplot");
296 gROOT->ForceStyle();
297
298 const size_t NUMBER_OF_PADS = 3;
299
300 cv->SetFillStyle(4000);
301 cv->SetFillColor(kWhite);
302
303 TPad* p1 = new TPad("p1", NULL, 0.0, 0.00, 1.0, 0.95);
304 TPad* p2 = new TPad("p2", NULL, 0.0, 0.95, 1.0, 1.00);
305
306 p1->Divide(NUMBER_OF_PADS, 1);
307
308 p1->Draw();
309 p2->Draw();
310
311 const double Dmax = 1000.0;
312 const double Rmin = 0.0;
313 const double Rmax = min(parameters.roadWidth_m, 0.4 * Dmax);
314 const double Tmin = min(parameters.TMin_ns, -10.0);
315 const double Tmax = max(parameters.TMax_ns, +100.0);
316 const double Amin = 0.002 * (Tmax - Tmin); // minimal arrow length [ns]
317 const double Amax = 0.8 * (Tmax - Tmin); // maximal arrow length [ns]
318 const double ymin = Tmin - (option == arrow_t ? 0.2 * Amax : 0.0);
319 const double ymax = Tmax + (option == arrow_t ? 0.5 * Amax : 0.0);
320
321 const string Xlabel[NUMBER_OF_PADS] = { "R [m]", "#phi [rad]", "z [m]" };
322 const double Xmin [NUMBER_OF_PADS] = { Rmin, -PI, -0.4 * Dmax };
323 const double Xmax [NUMBER_OF_PADS] = { Rmax, +PI, +0.4 * Dmax };
324
325 double Xs[NUMBER_OF_PADS];
326
327 for (size_t i = 0; i != NUMBER_OF_PADS; ++i) {
328 Xs[i] = 0.003 * (Xmax[i] - Xmin[i]) * (0.5 * NUMBER_OF_PMTS); // x-offset arrow as function of PMT number
329 }
330
331 TH2D H2[NUMBER_OF_PADS];
332 TGraph G2[NUMBER_OF_PADS];
333
334 for (size_t i = 0; i != NUMBER_OF_PADS; ++i) {
335
336 H2[i] = TH2D(MAKE_CSTRING("h" << i), NULL, graphics.nbinsX, Xmin[i] - Xs[i], Xmax[i] + Xs[i], graphics.nbinsY, ymin, ymax);
337
338 H2[i].GetXaxis()->SetTitle(Xlabel[i].c_str());
339 H2[i].GetYaxis()->SetTitle("#Deltat [ns]");
340
341 H2[i].GetXaxis()->CenterTitle(true);
342 H2[i].GetYaxis()->CenterTitle(true);
343
344 H2[i].SetStats(kFALSE);
345
346 G2[i].Set(2);
347
348 G2[i].SetPoint(0, H2[i].GetXaxis()->GetXmin(), 0.0);
349 G2[i].SetPoint(1, H2[i].GetXaxis()->GetXmax(), 0.0);
350
351 p1->cd(i+1);
352
353 H2[i].Draw("AXIS");
354 G2[i].Draw("SAME");
355 }
356
357
358 while (inputFile.hasNext()) {
359
360 cout << "event: " << setw(8) << inputFile.getCounter() << endl;
361
362 const Evt* evt = inputFile.next();
363
365
367
368 if (!event_selector(fit, *evt)) {
369 continue;
370 }
371
372 JDataL0_t dataL0;
373
374 for (const Hit& hit : evt->hits) {
375 dataL0.push_back(JHitL0(JDAQPMTIdentifier(hit.dom_id, hit.channel_id),
377 JHit(hit.t, hit.tot)));
378 }
379
380 summary.update(JDAQChronometer(evt->det_id,
381 evt->run_id,
382 evt->frame_index,
383 JDAQUTCExtended(evt->t.GetSec(), evt->t.GetNanoSec() / 16)));
384
385 const time_converter converter = time_converter(*evt);
386
387 Trk muon; // Monte Carlo true muon
388
389 if (has_muon(*evt)) {
390
391 for (const auto& trk : evt->mc_trks) {
392 if (is_muon(trk)) {
393 if (trk.E > muon.E) {
394
395 muon = trk;
396 muon.t += converter.putTime();
397
398 setW(muon, JSTART_LENGTH_METRES, fabs(muon.len));
399 }
400 }
401 }
402 }
403
404
405 bool monte_carlo = false; // show Monte Carlo true muon
406
407 for (bool next = false; !next; ) {
408
409 for (size_t i = 0; i != NUMBER_OF_PADS; ++i) {
410 H2[i].Reset();
411 }
412
413 Trk trk;
414
415 if (!monte_carlo)
416 trk = fit;
417 else
418 trk = muon;
419
420 JRotation3D R (getDirection(trk));
421 JLine1Z tz(getPosition (trk).rotate(R), trk.t);
422 JRange<double> Z_m;
423 /*
424 if (hasW(trk, JSTART_LENGTH_METRES)) {
425 Z_m = JRange<double>(0.0, getW(fit,JSTART_LENGTH_METRES)) + JRange<double>(parameters.ZMin_m, parameters.ZMax_m);
426 }
427 */
428 const JModel<JLine1Z> match(tz, parameters.roadWidth_m, T_ns, Z_m);
429
430 // hit selection based on fit result
431
432 JDataW0_t data;
433
434 for (JDataL0_t::const_iterator i = dataL0.begin(); i != dataL0.end(); ++i) {
435
436 JHitW0 hit(*i, summary.getRate(i->getPMTIdentifier(), parameters.R_Hz));
437
438 hit.rotate(R);
439
440 if (match(hit)) {
441 data.push_back(hit);
442 }
443 }
444
445 // select first hit in PMT
446
447 sort(data.begin(), data.end(), JHitW0::compare);
448
449 JDataW0_t::iterator __end = unique(data.begin(), data.end(), equal_to<JDAQPMTIdentifier>());
450
451 double E_GeV = parameters.E_GeV;
452 /*
453 if (trk.E > 0.1) {
454 E_GeV = trk.E;
455 }
456 */
457
458 // move fit to geometrical center of hits
459
461
462 for (JDataW0_t::iterator hit = data.begin(); hit != __end; ++hit) {
463
464 const double x = hit->getX() - tz.getX();
465 const double y = hit->getY() - tz.getY();
466 const double z = hit->getZ();
467 const double R = sqrt(x*x + y*y);
468
469 zs.include(z - R/getTanThetaC());
470 }
471
472 const double z0 = tz.getZ();
473 const double z1 = 0.5 * (zs.getLowerLimit() + zs.getUpperLimit());
474
475 tz.setZ(z1, getSpeedOfLight());
476
477
478 // graphics
479
480 ostringstream os;
481 vector<TArrow> arrow [NUMBER_OF_PADS];
482 vector<TMarker> marker[NUMBER_OF_PADS];
483
484 if (hasW(trk, JSTART_LENGTH_METRES) && getW(trk, JSTART_LENGTH_METRES) > 0.0) {
485
486 marker[2].push_back(TMarker(z0 - tz.getZ(), 0.0, kFullCircle));
487 marker[2].push_back(TMarker(z0 - tz.getZ() + getW(trk, JSTART_LENGTH_METRES), 0.0, kFullCircle));
488
489 static_cast<TAttMarker&>(marker[2][0]) = TAttMarker(kRed, kFullCircle, 0.7);
490 static_cast<TAttMarker&>(marker[2][1]) = TAttMarker(kRed, kFullCircle, 0.7);
491 }
492
493 DEBUG("trk: "
494 << FIXED(7,2) << tz.getX() << ' '
495 << FIXED(7,2) << tz.getY() << ' '
496 << FIXED(7,2) << tz.getZ() << ' '
497 << FIXED(12,2) << tz.getT() << endl);
498
499 double chi2 = 0;
500
501 for (JDataW0_t::const_iterator hit = data.begin(); hit != __end; ++hit) {
502
503 const double x = hit->getX() - tz.getX();
504 const double y = hit->getY() - tz.getY();
505 const double z = hit->getZ() - tz.getZ();
506 const double R = sqrt(x*x + y*y);
507
508 const double t1 = tz.getT() + (z + R * getTanThetaC()) * getInverseSpeedOfLight();
509
510 JDirection3D dir(hit->getDX(), hit->getDY(), hit->getDZ()); // PMT orientation
511
512 dir.rotate(JRotation3Z(-atan2(y,x))); // rotate PMT axis to x-z plane
513
514 const double theta = dir.getTheta();
515 const double phi = fabs(dir.getPhi()); // rotational symmetry of Cherenkov cone
516
517 //const double E = gWater.getE(E_GeV, z); // correct for energy loss
518 const double E = E_GeV;
519 const double dt = T_ns.constrain(hit->getT() - t1);
520
521 JMuonPDF_t::result_type H1 = pdf.calculate(E, R, theta, phi, dt);
522 JMuonPDF_t::result_type H0(hit->getR() * 1e-9, 0.0, T_ns);
523
524 if (H1.V >= parameters.VMax_npe) {
525 H1 *= parameters.VMax_npe / H1.V;
526 }
527
528 H1 += H0; // signal + background
529
530 chi2 += H1.getChi2() - H0.getChi2();
531
532 DEBUG("hit: "
533 << setw(8) << hit->getModuleID() << '.' << FILL(2,'0') << (int) hit->getPMTAddress() << FILL() << ' '
534 << SCIENTIFIC(8,2) << E << ' '
535 << FIXED(7,2) << R << ' '
536 << FIXED(7,4) << theta << ' '
537 << FIXED(7,4) << phi << ' '
538 << FIXED(7,3) << dt << ' '
539 << FIXED(7,3) << H1.getChi2() << ' '
540 << FIXED(7,3) << H0.getChi2() << endl);
541
542 const double derivative = H1.getDerivativeOfChi2() - H0.getDerivativeOfChi2();
543
544 double size = derivative * graphics.arrowScale; // size of arrow
545
546 if (fabs(size) < Amin) {
547 size = (size > 0.0 ? +Amin : -Amin);
548 } else if (fabs(size) > Amax) {
549 size = (size > 0.0 ? +Amax : -Amax);
550 }
551
552 const double X[NUMBER_OF_PADS] = { R, atan2(y,x), z - R/getTanThetaC() };
553
554 const double xs = (double) (NUMBER_OF_PMTS - 2 * hit->getPMTAddress()) / (double) NUMBER_OF_PMTS;
555
556 for (size_t i = 0; i != NUMBER_OF_PADS; ++i) {
557
558 TArrow a1(X[i] + xs*Xs[i], dt + graphics.T_ns, X[i] + xs*Xs[i], dt + graphics.T_ns + size, graphics.arrowSize, graphics.arrowType.c_str());
559
560 a1.SetLineWidth(graphics.lineWidth);
561 a1.SetLineStyle(graphics.lineStyle);
562
563 arrow[i].push_back(a1);
564
565 H2[i].Fill(X[i], dt + graphics.T_ns);
566 }
567 }
568
569 double zmax = 0.0;
570
571 for (size_t i = 0; i != NUMBER_OF_PADS; ++i) {
572 if (H2[i].GetMaximum() > zmax) {
573 zmax = H2[i].GetMaximum();
574 }
575 }
576
577 zmax *= 1.2;
578
579 for (size_t i = 0; i != NUMBER_OF_PADS; ++i) {
580 H2[i].SetMaximum(zmax);
581 }
582
583 os << FILL(6,'0') << evt->run_id << ":" << evt->frame_index << "/" << evt->trigger_counter << FILL();
584 os << " Q = " << FIXED(4,0) << trk.lik
585 << '/' << FIXED(4,0) << -chi2;
586 os << " E = " << SCIENTIFIC(7,1) << trk.E << " [GeV]";
587 os << " cos(#theta) = " << FIXED(6,3) << trk.dir.z;
588 os << " L = " << FIXED(6,2) << getW(trk, JSTART_LENGTH_METRES, 0.0) << " [m]";
589
590 if (monte_carlo)
591 os << " Monte Carlo";
592 else if (is_muon(muon))
593 os << " #Delta#alpha = " << FIXED(6,2) << getAngle(getDirection(muon), getDirection(trk)) << " [deg]";
594
595
596 // draw
597
598 TLatex title(0.05, 0.5, os.str().c_str());
599
600 title.SetTextAlign(12);
601 title.SetTextFont(42);
602 title.SetTextSize(0.6);
603
604 p2->cd();
605
606 title.Draw();
607
608 for (int i = 0; i != NUMBER_OF_PADS; ++i) {
609
610 p1->cd(i+1);
611
612 if (option == arrow_t) {
613
614 for (auto& a1 : arrow[i]) {
615 a1.Draw();
616 }
617
618 for (auto& m1 : marker[i]) {
619 m1.Draw();
620 }
621 }
622
623 if (option == histogram_t) {
624 H2[i].Draw("SAME");
625 }
626 }
627
628 cv->Update();
629
630
631 // action
632
633 if (batch) {
634
635 cv->SaveAs(replace(outputFile, WILDCARD, MAKE_STRING(inputFile.getCounter())).c_str());
636
637 next = true;
638
639 } else {
640
641 static int count = 0;
642
643 if (count++ == 0) {
644 cout << endl << "Type '?' for possible options." << endl;
645 }
646
647 for (bool user = true; user; ) {
648
649 cout << "\n> " << flush;
650
651 switch (JKeypress(true).get()) {
652
653 case '?':
654 cout << endl;
655 cout << "possible options: " << endl;
656 cout << 'q' << " -> " << "exit application" << endl;
657 cout << 'u' << " -> " << "update canvas" << endl;
658 cout << 's' << " -> " << "save graphics to file" << endl;
659 cout << 'M' << " -> " << "Monte Carlo true muon information" << endl;
660 cout << 'F' << " -> " << "fit information" << endl;
661 if (event_selector.is_valid()) {
662 cout << 'L' << " -> " << "reload event selector" << endl;
663 }
664 cout << 'r' << " -> " << "rewind input file" << endl;
665 cout << 'R' << " -> " << "switch to ROOT mode (quit ROOT to continue)" << endl;
666 cout << 'p' << " -> " << "print event information" << endl;
667 cout << ' ' << " -> " << "next event (as well as any other key)" << endl;
668 break;
669
670 case 'q':
671 cout << endl;
672 return 0;
673
674 case 'u':
675 cv->Update();
676 break;
677
678 case 's':
679 cv->SaveAs(replace(outputFile, WILDCARD, MAKE_STRING(inputFile.getCounter())).c_str());
680 break;
681
682 case 'M':
683 if (is_muon(muon))
684 monte_carlo = true;
685 else
686 ERROR(endl << "No Monte Carlo muon available." << endl);
687 user = false;
688 break;
689
690 case 'F':
691 monte_carlo = false;
692 user = false;
693 break;
694
695 case 'L':
696 if (event_selector.is_valid()) {
697 execute(MAKE_STRING("make -f " << getPath(argv[0]) << "/JMakeEventSelector libs"), 3);
698 event_selector.reload();
699 }
700 break;
701
702 case 'R':
703 tp->Run(kTRUE);
704 break;
705
706 case 'p':
707 cout << endl;
708 evt->print(cout);
709 if (debug >= debug_t) {
710 for (const auto& trk : evt->mc_trks) {
711 cout << "MC "; trk.print(cout); cout << endl;
712 }
713 for (const auto& trk : evt->trks) {
714 cout << "fit "; trk.print(cout); cout << endl;
715 }
716 for (const auto& hit : evt->hits) {
717 cout << "hit "; hit.print(cout); cout << endl;
718 }
719 }
720 break;
721
722 case 'r':
723 inputFile.rewind();
724
725 default:
726 next = true;
727 user = false;
728 break;
729 }
730 }
731 }
732 }
733 }
734 }
735 cout << endl;
736}
string outputFile
TPaveText * p1
#define DEBUG(A)
Message macros.
Definition JMessage.hh:62
#define FATAL(A)
Definition JMessage.hh:67
int debug
debug level
Definition JSirene.cc:72
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition JParser.hh:2142
#define MAKE_CSTRING(A)
Make C-string.
Definition JPrint.hh:72
#define MAKE_STRING(A)
Make string.
Definition JPrint.hh:63
#define gmake_property(A)
macros to convert (template) parameter to JPropertiesElement object
Utility class to parse parameter values.
Data structure for fit of straight line paralel to z-axis.
Definition JLine1Z.hh:29
Axis object.
Definition JAxis3D.hh:41
Data structure for direction in three dimensions.
Rotation around Z-axis.
Utility class to parse command line options.
Definition JParser.hh:1698
Auxiliary class for a hit with background rate value.
Definition JHitW0.hh:23
Data structure for size of TCanvas.
Definition JCanvas.hh:26
int y
number of pixels in Y
Definition JCanvas.hh:99
int x
number of pixels in X
Definition JCanvas.hh:98
Wrapper class around ROOT TStyle.
Definition JStyle.hh:24
Object reading from a list of files.
virtual void rewind() override
Rewind.
virtual bool hasNext() override
Check availability of next element.
counter_type getCounter() const
Get counter.
virtual const pointer_type & next() override
Get next element.
File router for fast addressing of summary data.
Enable unbuffered terminal input.
Definition JKeypress.hh:32
static JRange< T, JComparator_t > DEFAULT_RANGE()
Default range.
Definition JRange.hh:555
range_type & include(argument_type x)
Include given value to range.
Definition JRange.hh:397
T getLowerLimit() const
Get lower limit.
Definition JRange.hh:202
T getUpperLimit() const
Get upper limit.
Definition JRange.hh:213
Data structure for L0 hit.
Definition JHitL0.hh:31
static void setSlewing(const bool slewing)
Set slewing option.
Data structure for UTC time.
Auxiliary class to convert DAQ hit time to/from Monte Carlo hit time.
double putTime() const
Get Monte Carlo time minus DAQ/trigger time.
static const int JMUONGANDALF
static const int JMUONPREFIT
static const int JMUONENERGY
static const int JMUONSIMPLEX
static const int JMUONSTART
static const int JSTART_LENGTH_METRES
distance between projected positions on the track of optical modules for which the response does not ...
JDirection3D getDirection(const Vec &dir)
Get direction.
bool hasW(const Trk &trk, const int i)
Check availability of value.
bool has_muon(const Evt &evt)
Test whether given event has a muon.
void setW(Trk &trk, const int i, const double value)
Set associated value.
double getW(const Trk &track, const size_t index, const double value)
Get track information.
JPosition3D getPosition(const Vec &pos)
Get position.
bool is_muon(const Trk &track)
Test whether given track is a (anti-)muon.
@ debug_t
debug
Definition JMessage.hh:29
std::string getPath(const std::string &file_name)
Get path, i.e. part before last JEEP::PATHNAME_SEPARATOR if any.
double getAngle(const JQuaternion3D &first, const JQuaternion3D &second)
Get space angle between quanternions.
std::string replace(const std::string &input, const std::string &target, const std::string &replacement)
Replace tokens in string.
static const double PI
Mathematical constants.
const double getInverseSpeedOfLight()
Get inverse speed of light.
double getTanThetaC()
Get average tangent of Cherenkov angle of water corresponding to group velocity.
const double getSpeedOfLight()
Get speed of light.
This name space includes all other name spaces (except KM3NETDAQ, KM3NET and ANTARES).
const JFit & get_best_reconstructed_track(const JEvt &evt, JTrackSelector_t selector, JQualitySorter_t comparator)
Get best reconstructed track.
bool has_reconstructed_track(const JEvt &evt, JTrackSelector_t selector)
Test whether given event has a track according selection.
KM3NeT DAQ data structures and auxiliaries.
Definition DataQueue.cc:39
static const char WILDCARD
Definition JDAQTags.hh:56
The Evt class respresent a Monte Carlo (MC) event as well as an offline event.
Definition Evt.hh:21
int frame_index
from the raw data
Definition Evt.hh:29
int run_id
DAQ run identifier.
Definition Evt.hh:26
std::vector< Hit > hits
list of hits
Definition Evt.hh:38
std::vector< Trk > mc_trks
MC: list of MC truth tracks.
Definition Evt.hh:49
int det_id
detector identifier from DAQ
Definition Evt.hh:23
ULong64_t trigger_counter
trigger counter
Definition Evt.hh:31
void print(std::ostream &out=std::cout) const
Print event.
Definition Evt.hh:74
std::vector< Trk > trks
list of reconstructed tracks (can be several because of prefits,showers, etc).
Definition Evt.hh:39
TTimeStamp t
UTC time of the timeslice, or the event_time for MC. (default: 01 Jan 1970 00:00:00)
Definition Evt.hh:33
Auxiliary data structure for sequence of same character.
Definition JManip.hh:330
Auxiliary data structure for floating point format specification.
Definition JManip.hh:448
Definition Hit.hh:10
int dom_id
module identifier from the data (unique in the detector).
Definition Hit.hh:14
Vec pos
hit position
Definition Hit.hh:25
void print(std::ostream &out=std::cout) const
Print hit.
Definition Hit.hh:60
Vec dir
hit direction; i.e. direction of the PMT
Definition Hit.hh:26
unsigned int channel_id
PMT channel id {0,1, .., 30} local to moduke.
Definition Hit.hh:15
unsigned int tot
tot value as stored in raw data (int for pyroot)
Definition Hit.hh:17
double t
hit time (from tdc+calibration or MC truth)
Definition Hit.hh:23
Acoustics hit.
Model for fit to acoustics data.
bool is_valid() const
Check validity of function.
void reload()
Reload function from shared library.
Auxiliary data structure for muon PDF.
Definition JPDF_t.hh:135
JFunction1D_t::result_type result_type
Definition JPDF_t.hh:145
Empty structure for specification of parser element that is initialised (i.e. does not require input)...
Definition JParser.hh:68
double TMin_ns
minimal time w.r.t. Cherenkov hypothesis [ns]
double TMax_ns
maximal time w.r.t. Cherenkov hypothesis [ns]
double VMax_npe
maximum number of of photo-electrons
Auxiliary class for defining the range of iterations of objects.
Definition JLimit.hh:45
static counter_type max()
Get maximum counter value.
Definition JLimit.hh:128
Auxiliary data structure for floating point format specification.
Definition JManip.hh:488
The Trk class represents a Monte Carlo (MC) particle as well as a reconstructed track/shower.
Definition Trk.hh:15
void print(std::ostream &out=std::cout) const
Print track.
Definition Trk.hh:182
Vec dir
track direction
Definition Trk.hh:18
double E
Energy [GeV] (either MC truth or reconstructed)
Definition Trk.hh:20
double t
track time [ns] (when the particle is at pos )
Definition Trk.hh:19
double len
length, if applicable [m]
Definition Trk.hh:22
double lik
likelihood or lambda value (for aafit, lambda)
Definition Trk.hh:23
double z
Definition Vec.hh:14
Range of reconstruction stages.