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/tools/time_converter.hh"
#include "JROOT/JStyle.hh"
#include "JROOT/JCanvas.hh"
#include "JROOT/JRootToolkit.hh"
#include "JDAQ/JDAQEventIO.hh"
#include "JDAQ/JDAQSummarysliceIO.hh"
#include "JDetector/JDetector.hh"
#include "JDetector/JDetectorToolkit.hh"
#include "JDetector/JModuleRouter.hh"
#include "JGeometry3D/JCylinder3D.hh"
#include "JDynamics/JDynamics.hh"
#include "JTrigger/JHitL0.hh"
#include "JTrigger/JBuildL0.hh"
#include "JSupport/JSupport.hh"
#include "JSupport/JParallelFileScanner.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/JEvt.hh"
#include "JReconstruction/JEvtToolkit.hh"
#include "JReconstruction/JMuonGandalfParameters_t.hh"
#include "JReconstruction/JEventSelector.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.

Functions

int main (int argc, char **argv)
 

Detailed Description

Program to display hit probabilities.

Author
mdejong

Definition in file JReconstruction/JEvD.cc.

Function Documentation

◆ main()

int main ( int argc,
char ** argv )

Definition at line 108 of file JReconstruction/JEvD.cc.

109{
110 using namespace std;
111 using namespace JPP;
112 using namespace KM3NETDAQ;
113
115 typedef JParallelFileScanner_t::multi_pointer_type multi_pointer_type;
117 JACOUSTICS::JEvt>::typelist calibration_types;
118 typedef JMultipleFileScanner<calibration_types> JCalibration_t;
119
120 JParallelFileScanner_t inputFile;
121 JLimit_t& numberOfEvents = inputFile.getLimit();
122 string detectorFile;
123 JCalibration_t calibrationFile;
124 double Tmax_s;
125 string pdfFile;
126 string outputFile;
127 JMuonGandalfParameters_t parameters;
128 int application;
129 JEventSelector event_selector;
130 JCanvas canvas;
131 bool batch;
132 struct : JStyle::JParameters {
133 double arrowSize = 0.003;
134 string arrowType = "|->";
135 double arrowScale = 250.0;
136 Width_t lineWidth = 2;
137 Style_t lineStyle = 1;
138 int nbinsX = 50;
139 int nbinsY = 250;
140 double T_ns = 0.0;
141 } graphics;
142 string option;
143 int debug;
144
145
146 try {
147
148 JProperties properties = graphics.getProperties();
149
150 properties.insert(gmake_property(graphics.arrowSize));
151 properties.insert(gmake_property(graphics.arrowType));
152 properties.insert(gmake_property(graphics.arrowScale));
153 properties.insert(gmake_property(graphics.lineWidth));
154 properties.insert(gmake_property(graphics.lineStyle));
155 properties.insert(gmake_property(graphics.T_ns));
156
157 parameters.numberOfPrefits = 1;
158
159 JParser<> zap("Program to display hit probabilities.");
160
161 zap['w'] = make_field(canvas, "size of canvas <nx>x<ny> [pixels]") = JCanvas(1200, 600);
162 zap['f'] = make_field(inputFile, "input file (output of JXXXMuonReconstruction.sh)");
163 zap['a'] = make_field(detectorFile);
164 zap['+'] = make_field(calibrationFile) = JPARSER::initialised();
165 zap['T'] = make_field(Tmax_s) = 100.0;
166 zap['n'] = make_field(numberOfEvents) = JLimit::max();
167 zap['P'] = make_field(pdfFile);
168 zap['o'] = make_field(outputFile, "graphics output file name") = MAKE_STRING("display_" << WILDCARD << ".gif");
169 zap['@'] = make_field(parameters) = JPARSER::initialised();
170 zap['A'] = make_field(application) = JMUONGANDALF, JMUONENERGY, JMUONSTART;
171 zap['L'] = make_field(event_selector) = JPARSER::initialised();
172 zap['%'] = make_field(properties) = JPARSER::initialised();
173 zap['O'] = make_field(option, "draw option") = arrow_t, histogram_t;
174 zap['B'] = make_field(batch, "batch processing");
175 zap['d'] = make_field(debug) = 1;
176
177 zap(argc, argv);
178 }
179 catch(const exception& error) {
180 FATAL(error.what() << endl);
181 }
182
183 if (batch && outputFile == "") {
184 FATAL("Missing output file name " << outputFile << " in batch mode." << endl);
185 }
186
187 if (!batch && outputFile == "") {
188 outputFile = MAKE_STRING(WILDCARD << ".gif");
189 }
190
191 if (outputFile.find(WILDCARD) == string::npos) {
192 FATAL("Output file name " << outputFile << " has no wild card '" << WILDCARD << "'" << endl);
193 }
194
195
197
198 try {
199 load(detectorFile, detector);
200 }
201 catch(const JException& error) {
202 FATAL(error);
203 }
204 unique_ptr<JDynamics> dynamics;
205
206 try {
207
208 dynamics.reset(new JDynamics(detector, Tmax_s));
209
210 dynamics->load(calibrationFile);
211 }
212 catch(const exception& error) {
213 if (!calibrationFile.empty()) {
214 FATAL(error.what());
215 }
216 }
217
218 const double Zbed = 0.0;
219
220 JCylinder3D cylinder(detector.begin(), detector.end());
221
222 cylinder.addMargin(parameters.roadWidth_m);
223
224 if (cylinder.getZmin() < Zbed) {
225 cylinder.setZmin(Zbed);
226 }
227
228 const JModuleRouter router(dynamics ? dynamics->getDetector() : detector);
229
230 JSummaryFileRouter summary(inputFile);
231
232 const JMuonPDF_t pdf(pdfFile, parameters.TTS_ns);
233
234 const JTimeRange T_ns(parameters.TMin_ns, parameters.TMax_ns);
235
236 typedef vector<JHitL0> JDataL0_t;
237 typedef vector<JHitW0> JDataW0_t;
238
239 const JBuildL0<JHitL0> buildL0;
240
241
242 Vec offset(0.0, 0.0, 0.0);
243
244 try {
245 offset = getOffset(getHeader(inputFile));
246 } catch(const exception& error) {}
247
248
249 // ROOT
250
251 gROOT->SetBatch(batch);
252
253 TApplication* tp = new TApplication("user", NULL, NULL);
254 TCanvas* cv = new TCanvas("display", "", canvas.x, canvas.y);
255
256 if (!batch) {
257 ((TRootCanvas *) cv->GetCanvasImp())->Connect("CloseWindow()", "TApplication", tp, "Terminate()");
258 }
259
260 unique_ptr<TStyle> gStyle(new JStyle("gplot", cv->GetWw(), cv->GetWh(), graphics));
261
262 gROOT->SetStyle("gplot");
263 gROOT->ForceStyle();
264
265 const size_t NUMBER_OF_PADS = 3;
266
267 cv->SetFillStyle(4000);
268 cv->SetFillColor(kWhite);
269
270 TPad* p1 = new TPad("p1", NULL, 0.0, 0.00, 1.0, 0.95);
271 TPad* p2 = new TPad("p2", NULL, 0.0, 0.95, 1.0, 1.00);
272
273 p1->Divide(NUMBER_OF_PADS, 1);
274
275 p1->Draw();
276 p2->Draw();
277
278 const double Dmax = getMaximalDistance(detector);
279 const double Rmin = 0.0;
280 const double Rmax = min(parameters.roadWidth_m, 0.4 * Dmax);
281 const double Tmin = min(parameters.TMin_ns, -10.0);
282 const double Tmax = max(parameters.TMax_ns, +100.0);
283 const double Amin = 0.002 * (Tmax - Tmin); // minimal arrow length [ns]
284 const double Amax = 0.8 * (Tmax - Tmin); // maximal arrow length [ns]
285 const double ymin = Tmin - (option == arrow_t ? 0.2 * Amax : 0.0);
286 const double ymax = Tmax + (option == arrow_t ? 0.5 * Amax : 0.0);
287
288 const string Xlabel[NUMBER_OF_PADS] = { "R [m]", "#phi [rad]", "z [m]" };
289 const double Xmin [NUMBER_OF_PADS] = { Rmin, -PI, -0.4 * Dmax };
290 const double Xmax [NUMBER_OF_PADS] = { Rmax, +PI, +0.4 * Dmax };
291
292 double Xs[NUMBER_OF_PADS];
293
294 for (size_t i = 0; i != NUMBER_OF_PADS; ++i) {
295 Xs[i] = 0.003 * (Xmax[i] - Xmin[i]) * (0.5 * NUMBER_OF_PMTS); // x-offset arrow as function of PMT number
296 }
297
298 TH2D H2[NUMBER_OF_PADS];
299 TGraph G2[NUMBER_OF_PADS];
300
301 for (size_t i = 0; i != NUMBER_OF_PADS; ++i) {
302
303 H2[i] = TH2D(MAKE_CSTRING("h" << i), NULL, graphics.nbinsX, Xmin[i] - Xs[i], Xmax[i] + Xs[i], graphics.nbinsY, ymin, ymax);
304
305 H2[i].GetXaxis()->SetTitle(Xlabel[i].c_str());
306 H2[i].GetYaxis()->SetTitle("#Deltat [ns]");
307
308 H2[i].GetXaxis()->CenterTitle(true);
309 H2[i].GetYaxis()->CenterTitle(true);
310
311 H2[i].SetStats(kFALSE);
312
313 G2[i].Set(2);
314
315 G2[i].SetPoint(0, H2[i].GetXaxis()->GetXmin(), 0.0);
316 G2[i].SetPoint(1, H2[i].GetXaxis()->GetXmax(), 0.0);
317
318 p1->cd(i+1);
319
320 H2[i].Draw("AXIS");
321 G2[i].Draw("SAME");
322 }
323
324
325 for (JTreeScanner<Evt> mc(inputFile); inputFile.hasNext(); ) {
326
327 cout << "event: " << setw(8) << inputFile.getCounter() << endl;
328
329 multi_pointer_type ps = inputFile.next();
330
331 JDAQEvent* tev = ps;
332 JFIT::JEvt* in = ps;
333 Evt* event = NULL;
334
335 if (dynamics) {
336 dynamics->update(*tev);
337 }
338
339 if (mc.getEntries() != 0) {
340 event = mc.getEntry(tev->getCounter()); // Monte Carlo true information
341 }
342
343 in->select(JHistory::is_application(application));
344
345 if (!in->empty()) {
346
347 sort(in->begin(), in->end(), qualitySorter);
348
349 if (!event_selector(*tev, *in, event)) {
350 continue;
351 }
352
353
354 JDataL0_t dataL0;
355
356 buildL0(*tev, router, true, back_inserter(dataL0));
357
358 summary.update(*tev);
359
360 JFIT::JFit muon; // Monte Carlo true muon
361
362 if (event != NULL) {
363
364 const time_converter converter = time_converter(*event, *tev);
365
366 for (const auto& t1 : event->mc_trks) {
367 if (is_muon(t1)) {
368 if (t1.E > muon.getE()) {
369
370 JTrack3E ta = getTrack(t1);
371
372 ta.add(getPosition(offset));
373 ta.add(converter.putTime());
374
375 muon = getFit(0, ta, 0.0, 0, t1.E, 1);
376
377 muon.setW(JSTART_LENGTH_METRES, fabs(t1.len));
378 }
379 }
380 }
381 }
382
383 bool monte_carlo = false; // show Monte Carlo true muon
384 size_t index = 0; // index of fit
385
386 for (bool next = false; !next; ) {
387
388 for (size_t i = 0; i != NUMBER_OF_PADS; ++i) {
389 H2[i].Reset();
390 }
391
392 JFIT::JFit fit;
393
394 if (!monte_carlo)
395 fit = (*in)[index];
396 else
397 fit = muon;
398
399 JRotation3D R (getDirection(fit));
400 JLine1Z tz(getPosition (fit).rotate(R), fit.getT());
401 JRange<double> Z_m;
402 /*
403 if (fit.hasW(JSTART_LENGTH_METRES)) {
404 Z_m = JRange<double>(0.0, fit.getW(JSTART_LENGTH_METRES)) + JRange<double>(parameters.ZMin_m, parameters.ZMax_m);
405 }
406 */
407 const JFIT::JModel<JLine1Z> match(tz, parameters.roadWidth_m, T_ns, Z_m);
408
409 // hit selection based on fit result
410
411 JDataW0_t data;
412
413 for (JDataL0_t::const_iterator i = dataL0.begin(); i != dataL0.end(); ++i) {
414
415 JHitW0 hit(*i, summary.getRate(i->getPMTIdentifier(), parameters.R_Hz));
416
417 hit.rotate(R);
418
419 if (match(hit)) {
420 data.push_back(hit);
421 }
422 }
423
424 // select first hit in PMT
425
426 sort(data.begin(), data.end(), JHitW0::compare);
427
428 JDataW0_t::iterator __end = unique(data.begin(), data.end(), equal_to<JDAQPMTIdentifier>());
429
430 double E_GeV = parameters.E_GeV;
431 /*
432 if (fit.getE() > 0.1) {
433 E_GeV = fit.getE();
434 }
435 */
436
437 // move fit to geometrical center of hits
438
440
441 for (JDataW0_t::iterator hit = data.begin(); hit != __end; ++hit) {
442
443 const double x = hit->getX() - tz.getX();
444 const double y = hit->getY() - tz.getY();
445 const double z = hit->getZ();
446 const double R = sqrt(x*x + y*y);
447
448 zs.include(z - R/getTanThetaC());
449 }
450
451 const double z0 = tz.getZ();
452 const double z1 = 0.5 * (zs.getLowerLimit() + zs.getUpperLimit());
453
454 tz.setZ(z1, getSpeedOfLight());
455
456 // graphics
457
458 ostringstream os;
459 vector<TArrow> arrow [NUMBER_OF_PADS];
460 vector<TMarker> marker[NUMBER_OF_PADS];
461
462 if (fit.hasW(JSTART_LENGTH_METRES) && fit.getW(JSTART_LENGTH_METRES) > 0.0) {
463
464 marker[2].push_back(TMarker(z0 - tz.getZ(), 0.0, kFullCircle));
465 marker[2].push_back(TMarker(z0 - tz.getZ() + fit.getW(JSTART_LENGTH_METRES), 0.0, kFullCircle));
466
467 static_cast<TAttMarker&>(marker[2][0]) = TAttMarker(kRed, kFullCircle, 0.7);
468 static_cast<TAttMarker&>(marker[2][1]) = TAttMarker(kRed, kFullCircle, 0.7);
469 }
470
471 DEBUG("trk: "
472 << FIXED(7,2) << tz.getX() << ' '
473 << FIXED(7,2) << tz.getY() << ' '
474 << FIXED(7,2) << tz.getZ() << ' '
475 << FIXED(12,2) << tz.getT() << endl);
476
477 double chi2 = 0;
478
479 for (JDataW0_t::const_iterator hit = data.begin(); hit != __end; ++hit) {
480
481 const double x = hit->getX() - tz.getX();
482 const double y = hit->getY() - tz.getY();
483 const double z = hit->getZ() - tz.getZ();
484 const double R = sqrt(x*x + y*y);
485
486 const double t1 = tz.getT() + (z + R * getTanThetaC()) * getInverseSpeedOfLight();
487
488 JDirection3D dir(hit->getDX(), hit->getDY(), hit->getDZ()); // PMT orientation
489
490 dir.rotate(JRotation3Z(-atan2(y,x))); // rotate PMT axis to x-z plane
491
492 const double theta = dir.getTheta();
493 const double phi = fabs(dir.getPhi()); // rotational symmetry of Cherenkov cone
494
495 //const double E = gWater.getE(E_GeV, z); // correct for energy loss
496 const double E = E_GeV;
497 const double dt = T_ns.constrain(hit->getT() - t1);
498
499 JMuonPDF_t::result_type H1 = pdf.calculate(E, R, theta, phi, dt);
500 JMuonPDF_t::result_type H0(hit->getR() * 1e-9, 0.0, T_ns);
501
502 if (H1.V >= parameters.VMax_npe) {
503 H1 *= parameters.VMax_npe / H1.V;
504 }
505
506 H1 += H0; // signal + background
507
508 chi2 += H1.getChi2() - H0.getChi2();
509
510 DEBUG("hit: "
511 << setw(8) << hit->getModuleID() << '.' << FILL(2,'0') << (int) hit->getPMTAddress() << FILL() << ' '
512 << SCIENTIFIC(8,2) << E << ' '
513 << FIXED(7,2) << R << ' '
514 << FIXED(7,4) << theta << ' '
515 << FIXED(7,4) << phi << ' '
516 << FIXED(7,3) << dt << ' '
517 << FIXED(7,3) << H1.getChi2() << ' '
518 << FIXED(7,3) << H0.getChi2() << endl);
519
520 const double derivative = H1.getDerivativeOfChi2() - H0.getDerivativeOfChi2();
521
522 double size = derivative * graphics.arrowScale; // size of arrow
523
524 if (fabs(size) < Amin) {
525 size = (size > 0.0 ? +Amin : -Amin);
526 } else if (fabs(size) > Amax) {
527 size = (size > 0.0 ? +Amax : -Amax);
528 }
529
530 const double X[NUMBER_OF_PADS] = { R, atan2(y,x), z - R/getTanThetaC() };
531
532 const double xs = (double) (NUMBER_OF_PMTS - 2 * hit->getPMTAddress()) / (double) NUMBER_OF_PMTS;
533
534 for (size_t i = 0; i != NUMBER_OF_PADS; ++i) {
535
536 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());
537
538 a1.SetLineWidth(graphics.lineWidth);
539 a1.SetLineStyle(graphics.lineStyle);
540
541 arrow[i].push_back(a1);
542
543 H2[i].Fill(X[i], dt + graphics.T_ns);
544 }
545 }
546
547 double zmax = 0.0;
548
549 for (size_t i = 0; i != NUMBER_OF_PADS; ++i) {
550 if (H2[i].GetMaximum() > zmax) {
551 zmax = H2[i].GetMaximum();
552 }
553 }
554
555 zmax *= 1.2;
556
557 for (size_t i = 0; i != NUMBER_OF_PADS; ++i) {
558 H2[i].SetMaximum(zmax);
559 }
560
561 os << FILL(6,'0') << tev->getRunNumber() << ":" << tev->getFrameIndex() << "/" << tev->getCounter() << FILL();
562 os << " Q = " << FIXED(4,0) << fit.getQ()
563 << '/' << FIXED(4,0) << -chi2;
564 os << " E = " << SCIENTIFIC(7,1) << fit.getE() << " [GeV]";
565 os << " cos(#theta) = " << FIXED(6,3) << fit.getDZ();
566 os << " L = " << FIXED(6,2) << fit.getW(JSTART_LENGTH_METRES, 0.0) << " [m]";
567
568 if (monte_carlo)
569 os << " Monte Carlo";
570 else if (muon.getStatus() >= 0)
571 os << " #Delta#alpha = " << FIXED(6,2) << getAngle(getDirection(muon), getDirection(fit)) << " [deg]";
572
573
574 // draw
575
576 TLatex title(0.05, 0.5, os.str().c_str());
577
578 title.SetTextAlign(12);
579 title.SetTextFont(42);
580 title.SetTextSize(0.6);
581
582 p2->cd();
583
584 title.Draw();
585
586 for (int i = 0; i != NUMBER_OF_PADS; ++i) {
587
588 p1->cd(i+1);
589
590 if (option == arrow_t) {
591
592 for (auto& a1 : arrow[i]) {
593 a1.Draw();
594 }
595
596 for (auto& m1 : marker[i]) {
597 m1.Draw();
598 }
599 }
600
601 if (option == histogram_t) {
602 H2[i].Draw("SAME");
603 }
604 }
605
606 cv->Update();
607
608
609 // action
610
611 if (batch) {
612
613 cv->SaveAs(replace(outputFile, WILDCARD, MAKE_STRING(inputFile.getCounter())).c_str());
614
615 next = true;
616
617 } else {
618
619 static int count = 0;
620
621 if (count++ == 0) {
622 cout << endl << "Type '?' for possible options." << endl;
623 }
624
625 for (bool user = true; user; ) {
626
627 JAxis3D ts(getAxis(fit));
628
629 const JCylinder3D::intersection_type intersection = cylinder.getIntersection(ts);
630
631 ts.move(intersection.first);
632
633 cout << "\n> " << flush;
634
635 switch (JKeypress(true).get()) {
636
637 case '?':
638 cout << endl;
639 cout << "possible options: " << endl;
640 cout << 'p' << " -> " << "print information" << endl;
641 cout << 'q' << " -> " << "exit application" << endl;
642 cout << 'u' << " -> " << "update canvas" << endl;
643 cout << 's' << " -> " << "save graphics to file" << endl;
644 cout << '+' << " -> " << "next fit" << endl;
645 cout << '-' << " -> " << "previous fit" << endl;
646 cout << 'M' << " -> " << "Monte Carlo true muon information" << endl;
647 cout << 'F' << " -> " << "fit information" << endl;
648 if (event_selector.is_valid()) {
649 cout << 'L' << " -> " << "reload event selector" << endl;
650 }
651 cout << 'r' << " -> " << "rewind input file" << endl;
652 cout << 'R' << " -> " << "switch to ROOT mode (quit ROOT to continue)" << endl;
653 cout << ' ' << " -> " << "next event (as well as any other key)" << endl;
654 break;
655
656 case 'p':
657
658 cout << endl;
659 cout << "intersection: " << FIXED(6,1) << intersection.first << ' '<< FIXED(6,1) << intersection.second << endl;
660 cout << "entry point: "
661 << FIXED(6,1) << ts.getX() - cylinder.getX() << ' '
662 << FIXED(6,1) << ts.getY() - cylinder.getY() << ' '
663 << FIXED(6,1) << ts.getZ() << endl;
664 for (const auto& i : getWeight) {
665 cout << LEFT(32) << i.first << RIGHT(1) << ' ' << FIXED(12,5) << getWeight(fit, i.first, 0.0) << endl;
666 }
667 break;
668
669 case 'q':
670 cout << endl;
671 return 0;
672
673 case 'u':
674 cv->Update();
675 break;
676
677 case 's':
678 cv->SaveAs(replace(outputFile, WILDCARD, MAKE_STRING(inputFile.getCounter())).c_str());
679 break;
680
681 case '+':
682 monte_carlo = false;
683 index = (index != in->size() - 1 ? index + 1 : 0);
684 user = false;
685 break;
686
687 case '-':
688 monte_carlo = false;
689 index = (index != 0 ? index - 1 : in->size() - 1);
690 user = false;
691 break;
692
693 case 'M':
694 if (muon.getStatus() >= 0)
695 monte_carlo = true;
696 else
697 ERROR(endl << "No Monte Carlo muon available." << endl);
698 user = false;
699 break;
700
701 case 'F':
702 monte_carlo = false;
703 user = false;
704 break;
705
706 case 'L':
707 if (event_selector.is_valid()) {
708 execute(MAKE_STRING("make -f " << getPath(argv[0]) << "/JMakeEventSelector libs"), 3);
709 event_selector.reload();
710 }
711 break;
712
713 case 'R':
714 tp->Run(kTRUE);
715 break;
716
717 case 'r':
718 inputFile.rewind();
719
720 default:
721 next = true;
722 user = false;
723 break;
724 }
725 }
726 }
727 }
728 }
729 }
730 cout << endl;
731}
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
Detector data structure.
Definition JDetector.hh:96
Router for direct addressing of module data in detector data structure.
Utility class to parse parameter values.
Data structure for set of track fit results.
void select(const JSelector_t &selector)
Select fits.
Data structure for track fit results with history and optional associated values.
void setW(const std::vector< double > &W)
Set associated values.
double getDZ() const
Get Z-slope.
double getE() const
Get energy.
int getStatus() const
Get status of the fit; negative values should refer to a bad fit.
double getQ() const
Get quality.
const std::vector< double > & getW() const
Get associated values.
double getT() const
Get time.
bool hasW(const int i) const
Check availability of value.
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.
JTime & add(const JTime &value)
Addition operator.
3D track with energy.
Definition JTrack3E.hh:32
General exception.
Definition JException.hh:24
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
General purpose class for object reading from a list of file names.
General purpose class for parallel reading of objects from a single file or multiple files.
Object reading from a list of files.
File router for fast addressing of summary data.
Template definition for direct access of elements in ROOT TChain.
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
Template L0 hit builder.
Definition JBuildL0.hh:38
int getRunNumber() const
Get run number.
int getFrameIndex() const
Get frame index.
JTriggerCounter_t getCounter() const
Get trigger counter.
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 JMUONENERGY
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 ...
JAxis3D getAxis(const Trk &track)
Get axis.
JDirection3D getDirection(const Vec &dir)
Get direction.
JTrack3E getTrack(const Trk &track)
Get track.
JPosition3D getPosition(const Vec &pos)
Get position.
bool is_muon(const Trk &track)
Test whether given track is a (anti-)muon.
Vec getOffset(const JHead &header)
Get offset.
JFit getFit(const int id, const JMODEL::JString &string)
Get fit parameters of string.
void load(const std::string &file_name, JDetector &detector)
Load detector from input file.
double getMaximalDistance(const JDetector &detector, const bool option=false)
Get maximal distance between modules in detector.
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.
@ LEFT
Definition JTwosome.hh:18
@ RIGHT
Definition JTwosome.hh:18
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).
bool qualitySorter(const JFit &first, const JFit &second)
Comparison of fit results.
JRECONSTRUCTION::JWeight getWeight
Head getHeader(const JMultipleFileScanner_t &file_list)
Get Monte Carlo header.
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
Auxiliary data structure for sequence of same character.
Definition JManip.hh:330
Auxiliary data structure for floating point format specification.
Definition JManip.hh:448
Detector file.
Definition JHead.hh:227
Acoustic event fit.
Orientation of module.
Dynamic detector calibration.
Definition JDynamics.hh:86
bool is_valid() const
Check validity of function.
void reload()
Reload function from shared library.
Auxiliary class to test history.
Definition JHistory.hh:167
Auxiliary class to match data points with given model.
Auxiliary class for recursive type list generation.
Definition JTypeList.hh:351
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
const JLimit & getLimit() const
Get limit.
Definition JLimit.hh:84
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 Vec class is a straightforward 3-d vector, which also works in pyroot.
Definition Vec.hh:13