Jpp 20.0.0-195-g190c9e876
the software that should make you happy
Loading...
Searching...
No Matches
JAAnet/JEvD.cc
Go to the documentation of this file.
1#include <string>
2#include <iostream>
3#include <iomanip>
4#include <vector>
5#include <algorithm>
6#include <memory>
7
8#include "TROOT.h"
9#include "TApplication.h"
10#include "TCanvas.h"
11#include "TRootCanvas.h"
12#include "TStyle.h"
13#include "TH2D.h"
14#include "TArrow.h"
15#include "TLatex.h"
16#include "TMarker.h"
17
23
24#include "JROOT/JStyle.hh"
25#include "JROOT/JCanvas.hh"
26#include "JROOT/JRootToolkit.hh"
27
29
30#include "JTrigger/JHitL0.hh"
31
32#include "JSupport/JSupport.hh"
37
38#include "JPhysics/JPDF_t.hh"
39
40#include "JFit/JLine1Z.hh"
41#include "JFit/JModel.hh"
44
45#include "JLang/JPredicate.hh"
46#include "JLang/JComparator.hh"
47
48#include "JMath/JMathToolkit.hh"
49#include "JSystem/JKeypress.hh"
50#include "JSystem/JProcess.hh"
51
53#include "Jeep/JProperties.hh"
54#include "Jeep/JPrint.hh"
55#include "Jeep/JParser.hh"
56#include "Jeep/JMessage.hh"
57
58
59namespace JAANET {
60
62
63 /**
64 * Event selector.
65 *
66 * The default constructor will accept all events.\n
67 * A different method can dynamically be loaded from a shared library using class JEEP::JFunctionAdaptor.
68 */
70 public JFunctionAdaptor<bool, const Trk&, const Evt&>
71 {
72 /**
73 * Default event selection.
74 *
75 * \param trk track
76 * \param evt event
77 * \return true
78 */
79 static inline bool select(const Trk& trk, const Evt& evt)
80 {
81 return true;
82 }
83
84
85 /**
86 * Default constructor.
87 */
89 {
90 this->function = select;
91 this->symbol = "select";
92 }
93 };
94
95
96 /**
97 * Check availability of value.
98 *
99 * \param trk track
100 * \param i index
101 * \return true if available; else false
102 */
103 inline bool hasW(const Trk& trk, const int i)
104 {
105 return (i >= 0 && i < (int) trk.fitinf.size());
106 }
107
108
109 /**
110 * Get associated value.
111 *
112 * \param trk track
113 * \param i index
114 * \return value
115 */
116 inline double getW(const Trk& trk, const int i)
117 {
118 return trk.fitinf.at(i);
119 }
120
121
122 /**
123 * Get associated value.
124 *
125 * \param trk track
126 * \param i index
127 * \param value default value
128 * \return value
129 */
130 inline double getW(const Trk& trk, const int i, const double value)
131 {
132 if (hasW(trk,i))
133 return trk.fitinf.at(i);
134 else
135 return value;
136 }
137
138
139 /**
140 * Set associated value.
141 *
142 * \param trk track
143 * \param i index
144 * \param value value
145 */
146 void setW(Trk& trk, const int i, const double value)
147 {
148 if (i >= (int) trk.fitinf.size()) {
149 trk.fitinf.resize(i + 1, 0.0);
150 }
151
152 trk.fitinf[i] = value;
153 }
154}
155
156namespace {
157
158 /**
159 * Wild card character for file name substition.
160 */
161 const char WILDCARD = '%';
162
163
164 /**
165 * Execute command in shell.
166 *
167 * \param command command
168 */
169 inline void execute(const std::string& command, int debug)
170 {
171 using namespace std;
172 using namespace JPP;
173
174 JProcess process(command);
175
176 istream in(process.getInputStreamBuffer());
177
178 for (string buffer; getline(in, buffer); ) {
179 DEBUG(buffer << endl);
180 }
181 }
182
183 const char* const histogram_t = "histogram"; //!< draw histogram
184 const char* const arrow_t = "arrow"; //!< draw arrow
185}
186
187
188/**
189 * \file
190 *
191 * Program to display hit probabilities.
192 *
193 * \author mdejong
194 */
195int main(int argc, char **argv)
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 bool equalize = false;
220 } graphics;
221 string option;
222 int debug;
223
224
225 try {
226
227 JProperties properties = graphics.getProperties();
228
229 properties.insert(gmake_property(graphics.arrowSize));
230 properties.insert(gmake_property(graphics.arrowType));
231 properties.insert(gmake_property(graphics.arrowScale));
232 properties.insert(gmake_property(graphics.lineWidth));
233 properties.insert(gmake_property(graphics.lineStyle));
234 properties.insert(gmake_property(graphics.T_ns));
235 properties.insert(gmake_property(graphics.equalize));
236
237 parameters.numberOfPrefits = 1;
238
239 JParser<> zap("Program to display hit probabilities.");
240
241 zap['w'] = make_field(canvas, "size of canvas <nx>x<ny> [pixels]") = JCanvas(1200, 600);
242 zap['f'] = make_field(inputFile, "input file (output of JXXXMuonReconstruction.sh)");
243 zap['n'] = make_field(numberOfEvents) = JLimit::max();
244 zap['F'] = make_field(pdfFile);
245 zap['o'] = make_field(outputFile, "graphics output file name") = MAKE_STRING("display_" << WILDCARD << ".gif");
246 zap['@'] = make_field(parameters) = JPARSER::initialised();
248 zap['L'] = make_field(event_selector) = JPARSER::initialised();
249 zap['%'] = make_field(properties) = JPARSER::initialised();
250 zap['O'] = make_field(option, "draw option") = arrow_t, histogram_t;
251 zap['B'] = make_field(batch, "batch processing");
252 zap['d'] = make_field(debug) = 1;
253
254 zap(argc, argv);
255 }
256 catch(const exception& error) {
257 FATAL(error.what() << endl);
258 }
259
260 if (batch && outputFile == "") {
261 FATAL("Missing output file name " << outputFile << " in batch mode." << endl);
262 }
263
264 if (!batch && outputFile == "") {
265 outputFile = MAKE_STRING(WILDCARD << ".gif");
266 }
267
268 if (outputFile.find(WILDCARD) == string::npos) {
269 FATAL("Output file name " << outputFile << " has no wild card '" << WILDCARD << "'" << endl);
270 }
271
272 JHit::setSlewing(false);
273
274 JSummaryFileRouter summary(inputFile);
275
276 const JMuonPDF_t pdf(pdfFile, parameters.TTS_ns);
277
278 const JTimeRange T_ns(parameters.TMin_ns, parameters.TMax_ns);
279
280 typedef vector<JHitL0> JDataL0_t;
281 typedef vector<JHitW0> JDataW0_t;
282
283
284 // ROOT
285
286 gROOT->SetBatch(batch);
287
288 TApplication* tp = new TApplication("user", NULL, NULL);
289 TCanvas* cv = new TCanvas("display", "", canvas.x, canvas.y);
290
291 if (!batch) {
292 ((TRootCanvas *) cv->GetCanvasImp())->Connect("CloseWindow()", "TApplication", tp, "Terminate()");
293 }
294
295 unique_ptr<TStyle> gStyle(new JStyle("gplot", cv->GetWw(), cv->GetWh(), graphics));
296
297 gROOT->SetStyle("gplot");
298 gROOT->ForceStyle();
299
300 const size_t NUMBER_OF_PADS = 3;
301
302 cv->SetFillStyle(4000);
303 cv->SetFillColor(kWhite);
304
305 TPad* p1 = new TPad("p1", NULL, 0.0, 0.00, 1.0, 0.95);
306 TPad* p2 = new TPad("p2", NULL, 0.0, 0.95, 1.0, 1.00);
307
308 p1->Divide(NUMBER_OF_PADS, 1);
309
310 p1->Draw();
311 p2->Draw();
312
313 const double Dmax = 1000.0;
314 const double Rmin = 0.0;
315 const double Rmax = min(parameters.roadWidth_m, 0.4 * Dmax);
316 const double Tmin = min(parameters.TMin_ns, -10.0);
317 const double Tmax = max(parameters.TMax_ns, +100.0);
318 const double Amin = 0.002 * (Tmax - Tmin); // minimal arrow length [ns]
319 const double Amax = 0.8 * (Tmax - Tmin); // maximal arrow length [ns]
320 const double ymin = Tmin - (option == arrow_t ? 0.2 * Amax : 0.0);
321 const double ymax = Tmax + (option == arrow_t ? 0.5 * Amax : 0.0);
322
323 const string Xlabel[NUMBER_OF_PADS] = { "R [m]", "#phi [rad]", "z [m]" };
324 const double Xmin [NUMBER_OF_PADS] = { Rmin, -PI, -0.4 * Dmax };
325 const double Xmax [NUMBER_OF_PADS] = { Rmax, +PI, +0.4 * Dmax };
326
327 double Xs[NUMBER_OF_PADS];
328
329 for (size_t i = 0; i != NUMBER_OF_PADS; ++i) {
330 Xs[i] = 0.003 * (Xmax[i] - Xmin[i]) * (0.5 * NUMBER_OF_PMTS); // x-offset arrow as function of PMT number
331 }
332
333 TH2D H2[NUMBER_OF_PADS];
334 TGraph G2[NUMBER_OF_PADS];
335
336 for (size_t i = 0; i != NUMBER_OF_PADS; ++i) {
337
338 H2[i] = TH2D(MAKE_CSTRING("h" << i), NULL, graphics.nbinsX, Xmin[i] - Xs[i], Xmax[i] + Xs[i], graphics.nbinsY, ymin, ymax);
339
340 H2[i].GetXaxis()->SetTitle(Xlabel[i].c_str());
341 H2[i].GetYaxis()->SetTitle("#Deltat [ns]");
342
343 H2[i].GetXaxis()->CenterTitle(true);
344 H2[i].GetYaxis()->CenterTitle(true);
345
346 H2[i].SetStats(kFALSE);
347
348 G2[i].Set(2);
349
350 G2[i].SetPoint(0, H2[i].GetXaxis()->GetXmin(), 0.0);
351 G2[i].SetPoint(1, H2[i].GetXaxis()->GetXmax(), 0.0);
352
353 p1->cd(i+1);
354
355 H2[i].Draw("AXIS");
356 G2[i].Draw("SAME");
357 }
358
359
360 while (inputFile.hasNext()) {
361
362 cout << "event: " << setw(8) << inputFile.getCounter() << endl;
363
364 const Evt* evt = inputFile.next();
365
367
369
370 if (!event_selector(fit, *evt)) {
371 continue;
372 }
373
374 JDataL0_t dataL0;
375
376 for (const Hit& hit : evt->hits) {
377 dataL0.push_back(JHitL0(JDAQPMTIdentifier(hit.dom_id, hit.channel_id),
379 JHit(hit.t, hit.tot)));
380 }
381
382 summary.update(JDAQChronometer(evt->det_id,
383 evt->run_id,
384 evt->frame_index,
385 JDAQUTCExtended(evt->t.GetSec(), evt->t.GetNanoSec() / 16)));
386
387 const time_converter converter = time_converter(*evt);
388
389 Trk muon; // Monte Carlo true muon
390
391 if (has_muon(*evt)) {
392
393 for (const auto& trk : evt->mc_trks) {
394 if (is_muon(trk)) {
395 if (trk.E > muon.E) {
396
397 muon = trk;
398 muon.t += converter.putTime();
399
400 setW(muon, JSTART_LENGTH_METRES, fabs(muon.len));
401 }
402 }
403 }
404 }
405
406
407 bool monte_carlo = false; // show Monte Carlo true muon
408
409 for (bool next = false; !next; ) {
410
411 for (size_t i = 0; i != NUMBER_OF_PADS; ++i) {
412 H2[i].Reset();
413 }
414
415 Trk trk;
416
417 if (!monte_carlo)
418 trk = fit;
419 else
420 trk = muon;
421
422 JRotation3D R (getDirection(trk));
423 JLine1Z tz(getPosition (trk).rotate(R), trk.t);
424 JRange<double> Z_m;
425 /*
426 if (hasW(trk, JSTART_LENGTH_METRES)) {
427 Z_m = JRange<double>(0.0, getW(fit,JSTART_LENGTH_METRES)) + JRange<double>(parameters.ZMin_m, parameters.ZMax_m);
428 }
429 */
430 const JModel<JLine1Z> match(tz, parameters.roadWidth_m, T_ns, Z_m);
431
432 // hit selection based on fit result
433
434 JDataW0_t data;
435
436 for (JDataL0_t::const_iterator i = dataL0.begin(); i != dataL0.end(); ++i) {
437
438 const int type = 0;
439 const double QE = 1.0;
440 const double R_Hz = summary.getRate(i->getPMTIdentifier(), parameters.R_Hz);
441
442 JHitW0 hit(*i, type, QE, R_Hz);
443
444 hit.rotate(R);
445
446 if (match(hit)) {
447 data.push_back(hit);
448 }
449 }
450
451 // select first hit in PMT
452
453 sort(data.begin(), data.end(), JHitW0::compare);
454
455 JDataW0_t::iterator __end = unique(data.begin(), data.end(), equal_to<JDAQPMTIdentifier>());
456
457 double E_GeV = parameters.E_GeV;
458 /*
459 if (trk.E > 0.1) {
460 E_GeV = trk.E;
461 }
462 */
463
464 // move fit to geometrical center of hits
465
467
468 for (JDataW0_t::iterator hit = data.begin(); hit != __end; ++hit) {
469
470 const double x = hit->getX() - tz.getX();
471 const double y = hit->getY() - tz.getY();
472 const double z = hit->getZ();
473 const double R = sqrt(x*x + y*y);
474
475 zs.include(z - R/getTanThetaC());
476 }
477
478 const double z0 = tz.getZ();
479 const double z1 = 0.5 * (zs.getLowerLimit() + zs.getUpperLimit());
480
481 tz.setZ(z1, getSpeedOfLight());
482
483
484 // graphics
485
486 ostringstream os;
487 vector<TArrow> arrow [NUMBER_OF_PADS];
488 vector<TMarker> marker[NUMBER_OF_PADS];
489
490 if (hasW(trk, JSTART_LENGTH_METRES) && getW(trk, JSTART_LENGTH_METRES) > 0.0) {
491
492 marker[2].push_back(TMarker(z0 - tz.getZ(), 0.0, kFullCircle));
493 marker[2].push_back(TMarker(z0 - tz.getZ() + getW(trk, JSTART_LENGTH_METRES), 0.0, kFullCircle));
494
495 static_cast<TAttMarker&>(marker[2][0]) = TAttMarker(kRed, kFullCircle, 0.7);
496 static_cast<TAttMarker&>(marker[2][1]) = TAttMarker(kRed, kFullCircle, 0.7);
497 }
498
499 DEBUG("trk: "
500 << FIXED(7,2) << tz.getX() << ' '
501 << FIXED(7,2) << tz.getY() << ' '
502 << FIXED(7,2) << tz.getZ() << ' '
503 << FIXED(12,2) << tz.getT() << endl);
504
505 double chi2 = 0;
506
507 for (JDataW0_t::const_iterator hit = data.begin(); hit != __end; ++hit) {
508
509 const double x = hit->getX() - tz.getX();
510 const double y = hit->getY() - tz.getY();
511 const double z = hit->getZ() - tz.getZ();
512 const double R = sqrt(x*x + y*y);
513
514 const double t1 = tz.getT() + (z + R * getTanThetaC()) * getInverseSpeedOfLight();
515
516 JDirection3D dir(hit->getDX(), hit->getDY(), hit->getDZ()); // PMT orientation
517
518 dir.rotate(JRotation3Z(-atan2(y,x))); // rotate PMT axis to x-z plane
519
520 const double theta = dir.getTheta();
521 const double phi = fabs(dir.getPhi()); // rotational symmetry of Cherenkov cone
522
523 //const double E = gWater.getE(E_GeV, z); // correct for energy loss
524 const double E = E_GeV;
525 const double dt = T_ns.constrain(hit->getT() - t1);
526
527 JMuonPDF_t::result_type H1 = pdf.calculate(E, R, theta, phi, dt);
528 JMuonPDF_t::result_type H0(hit->getR() * 1e-9, 0.0, T_ns);
529
530 if (H1.V >= parameters.VMax_npe) {
531 H1 *= parameters.VMax_npe / H1.V;
532 }
533
534 H1 += H0; // signal + background
535
536 chi2 += H1.getChi2() - H0.getChi2();
537
538 DEBUG("hit: "
539 << setw(8) << hit->getModuleID() << '.' << FILL(2,'0') << (int) hit->getPMTAddress() << FILL() << ' '
540 << SCIENTIFIC(8,2) << E << ' '
541 << FIXED(7,2) << R << ' '
542 << FIXED(7,4) << theta << ' '
543 << FIXED(7,4) << phi << ' '
544 << FIXED(7,3) << dt << ' '
545 << FIXED(7,3) << H1.getChi2() << ' '
546 << FIXED(7,3) << H0.getChi2() << endl);
547
548 const double derivative = H1.getDerivativeOfChi2() - H0.getDerivativeOfChi2();
549
550 double size = derivative * graphics.arrowScale; // size of arrow
551
552 if (fabs(size) < Amin) {
553 size = (size > 0.0 ? +Amin : -Amin);
554 } else if (fabs(size) > Amax) {
555 size = (size > 0.0 ? +Amax : -Amax);
556 }
557
558 const double X[NUMBER_OF_PADS] = { R, atan2(y,x), z - R/getTanThetaC() };
559
560 const double xs = (double) (NUMBER_OF_PMTS - 2 * hit->getPMTAddress()) / (double) NUMBER_OF_PMTS;
561
562 for (size_t i = 0; i != NUMBER_OF_PADS; ++i) {
563
564 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());
565
566 a1.SetLineWidth(graphics.lineWidth);
567 a1.SetLineStyle(graphics.lineStyle);
568
569 arrow[i].push_back(a1);
570
571 H2[i].Fill(X[i], dt + graphics.T_ns);
572 }
573 }
574
575 if (graphics.equalize) {
576
577 double zmax = 0.0;
578
579 for (size_t i = 0; i != NUMBER_OF_PADS; ++i) {
580 if (H2[i].GetMaximum() > zmax) {
581 zmax = H2[i].GetMaximum();
582 }
583 }
584
585 zmax *= 1.2;
586
587 for (size_t i = 0; i != NUMBER_OF_PADS; ++i) {
588 H2[i].SetMaximum(zmax);
589 }
590 }
591
592 os << FILL(6,'0') << evt->run_id << ":" << evt->frame_index << "/" << evt->trigger_counter << FILL();
593 os << " Q = " << FIXED(4,0) << trk.lik
594 << '/' << FIXED(4,0) << -chi2;
595 os << " E = " << SCIENTIFIC(7,1) << trk.E << " [GeV]";
596 os << " cos(#theta) = " << FIXED(6,3) << trk.dir.z;
597 os << " L = " << FIXED(6,2) << getW(trk, JSTART_LENGTH_METRES, 0.0) << " [m]";
598
599 if (monte_carlo)
600 os << " Monte Carlo";
601 else if (is_muon(muon))
602 os << " #Delta#alpha = " << FIXED(6,2) << getAngle(getDirection(muon), getDirection(trk)) << " [deg]";
603
604
605 // draw
606
607 TLatex title(0.05, 0.5, os.str().c_str());
608
609 title.SetTextAlign(12);
610 title.SetTextFont(42);
611 title.SetTextSize(0.6);
612
613 p2->cd();
614
615 title.Draw();
616
617 for (int i = 0; i != NUMBER_OF_PADS; ++i) {
618
619 p1->cd(i+1);
620
621 if (option == arrow_t) {
622
623 for (auto& a1 : arrow[i]) {
624 a1.Draw();
625 }
626
627 for (auto& m1 : marker[i]) {
628 m1.Draw();
629 }
630 }
631
632 if (option == histogram_t) {
633 H2[i].Draw("SAME");
634 }
635 }
636
637 cv->Update();
638
639
640 // action
641
642 if (batch) {
643
644 cv->SaveAs(replace(outputFile, WILDCARD, MAKE_STRING(inputFile.getCounter())).c_str());
645
646 next = true;
647
648 } else {
649
650 static int count = 0;
651
652 if (count++ == 0) {
653 cout << endl << "Type '?' for possible options." << endl;
654 }
655
656 for (bool user = true; user; ) {
657
658 cout << "\n> " << flush;
659
660 switch (JKeypress(true).get()) {
661
662 case '?':
663 cout << endl;
664 cout << "possible options: " << endl;
665 cout << 'q' << " -> " << "exit application" << endl;
666 cout << 'u' << " -> " << "update canvas" << endl;
667 cout << 's' << " -> " << "save graphics to file" << endl;
668 cout << 'M' << " -> " << "Monte Carlo true muon information" << endl;
669 cout << 'F' << " -> " << "fit information" << endl;
670 if (event_selector.is_valid()) {
671 cout << 'L' << " -> " << "reload event selector" << endl;
672 }
673 cout << 'r' << " -> " << "rewind input file" << endl;
674 cout << 'R' << " -> " << "switch to ROOT mode (quit ROOT to continue)" << endl;
675 cout << 'p' << " -> " << "print event information" << endl;
676 cout << ' ' << " -> " << "next event (as well as any other key)" << endl;
677 break;
678
679 case 'q':
680 cout << endl;
681 return 0;
682
683 case 'u':
684 cv->Update();
685 break;
686
687 case 's':
688 cv->SaveAs(replace(outputFile, WILDCARD, MAKE_STRING(inputFile.getCounter())).c_str());
689 break;
690
691 case 'M':
692 if (is_muon(muon))
693 monte_carlo = true;
694 else
695 ERROR(endl << "No Monte Carlo muon available." << endl);
696 user = false;
697 break;
698
699 case 'F':
700 monte_carlo = false;
701 user = false;
702 break;
703
704 case 'L':
705 if (event_selector.is_valid()) {
706 execute(MAKE_STRING("make -f " << getPath(argv[0]) << "/JMakeEventSelector libs"), 3);
707 event_selector.reload();
708 }
709 break;
710
711 case 'R':
712 tp->Run(kTRUE);
713 break;
714
715 case 'p':
716 cout << endl;
717 evt->print(cout);
718 if (debug >= debug_t) {
719 for (const auto& trk : evt->mc_trks) {
720 cout << "MC "; trk.print(cout); cout << endl;
721 }
722 for (const auto& trk : evt->trks) {
723 cout << "fit "; trk.print(cout); cout << endl;
724 }
725 for (const auto& hit : evt->hits) {
726 cout << "hit "; hit.print(cout); cout << endl;
727 }
728 }
729 break;
730
731 case 'r':
732 inputFile.rewind();
733
734 default:
735 next = true;
736 user = false;
737 break;
738 }
739 }
740 }
741 }
742 }
743 }
744 cout << endl;
745}
Definition of hit and track types and auxiliary methods for handling Monte Carlo data.
int main(int argc, char **argv)
string outputFile
KM3NeT DAQ constants, bit handling, etc.
TPaveText * p1
Basic data structure for L0 hit.
Keyboard settings for unbuffered input.
Binary methods for member methods.
General purpose messaging.
#define DEBUG(A)
Message macros.
Definition JMessage.hh:62
#define FATAL(A)
Definition JMessage.hh:67
int debug
debug level
Definition JSirene.cc:74
Auxiliary data structure for muon PDF.
Utility class to parse command line options.
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition JParser.hh:2142
I/O formatting auxiliaries.
#define MAKE_CSTRING(A)
Make C-string.
Definition JPrint.hh:72
#define MAKE_STRING(A)
Make string.
Definition JPrint.hh:63
Utility class to parse parameter values.
#define gmake_property(A)
macros to convert (template) parameter to JPropertiesElement object
Scanning of objects from a single file according a format that follows from the extension of each fil...
ROOT TTree parameter settings of various packages.
Utility class to parse parameter values.
Data structure for fit of straight line paralel to z-axis.
Definition JLine1Z.hh:29
double getT(const JVector3D &pos) const
Get arrival time of Cherenkov light at given position.
Definition JLine1Z.hh:114
double getZ(const JPosition3D &pos) const
Get point of emission of Cherenkov light along muon path.
Definition JLine1Z.hh:134
void setZ(const double z, const double velocity)
Set z-position of vertex.
Definition JLine1Z.hh:75
Axis object.
Definition JAxis3D.hh:41
JAxis3D & rotate(const JRotation3D &R)
Rotate axis.
Definition JAxis3D.hh:225
Data structure for direction in three dimensions.
JDirection3D & rotate(const JRotation3D &R)
Rotate.
Rotation around Z-axis.
double getY() const
Get y position.
Definition JVector3D.hh:104
double getX() const
Get x position.
Definition JVector3D.hh:94
double getTheta() const
Get theta angle.
Definition JVersor3D.hh:128
double getPhi() const
Get phi angle.
Definition JVersor3D.hh:144
Utility class to parse command line options.
Definition JParser.hh:1698
Auxiliary class for a hit with background rate value.
Definition JHitW0.hh:25
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.
void update(const JDAQHeader &header)
Update router.
double getRate(const JDAQPMTIdentifier &id) const
Get rate.
Enable unbuffered terminal input.
Definition JKeypress.hh:32
Streaming of input and output from Linux command.
Definition JProcess.hh:30
T constrain(argument_type x) const
Constrain value to range.
Definition JRange.hh:350
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 JSTART_LENGTH_METRES
distance between projected positions on the track of optical modules for which the response does not ...
Extensions to Evt data format.
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::istream & getline(std::istream &in, JString &object)
Read string from input stream until end of line.
Definition JString.hh:478
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
Event selector.
JEventSelector()
Default constructor.
static bool select(const Trk &trk, const Evt &evt)
Default event selection.
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
result_type calculate(const double E, const double R, const double theta, const double phi, const double t1) const
Get PDF.
Definition JPDF_t.hh:235
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
std::vector< double > fitinf
place to store additional fit info, see km3net-dataformat/definitions/fitparameters....
Definition Trk.hh:32
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.
Auxiliary include file for time conversion between DAQ/trigger hit and Monte Carlo hit.