Jpp test-rotations-old
the software that should make you happy
Loading...
Searching...
No Matches
JReconstruction/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 "TStyle.h"
12#include "TH2D.h"
13#include "TArrow.h"
14#include "TLatex.h"
15#include "TMarker.h"
16
21
22#include "JROOT/JStyle.hh"
23#include "JROOT/JCanvas.hh"
24#include "JROOT/JRootToolkit.hh"
25
26#include "JDAQ/JDAQEventIO.hh"
28
32
34
35#include "JTrigger/JHitL0.hh"
36#include "JTrigger/JBuildL0.hh"
37
38#include "JSupport/JSupport.hh"
43
44#include "JPhysics/JPDF_t.hh"
45
46#include "JFit/JLine1Z.hh"
47#include "JFit/JModel.hh"
53
54#include "JLang/JPredicate.hh"
55#include "JLang/JComparator.hh"
56
57#include "JMath/JMathToolkit.hh"
58#include "JSystem/JKeypress.hh"
59#include "JSystem/JProcess.hh"
60
62#include "Jeep/JProperties.hh"
63#include "Jeep/JPrint.hh"
64#include "Jeep/JParser.hh"
65#include "Jeep/JMessage.hh"
66
67
68namespace {
69
70 /**
71 * Wild card character for file name substition.
72 */
73 const char WILDCARD = '%';
74
75 /**
76 * Execute command in shell.
77 *
78 * \param command command
79 */
80 inline void execute(const std::string& command, int debug)
81 {
82 using namespace std;
83 using namespace JPP;
84
85 JProcess process(command);
86
87 istream in(process.getInputStreamBuffer());
88
89 for (string buffer; getline(in, buffer); ) {
90 DEBUG(buffer << endl);
91 }
92 }
93
94 const char* const histogram_t = "histogram"; //!< draw histogram
95 const char* const arrow_t = "arrow"; //!< draw arrow
96}
97
98
99/**
100 * \file
101 *
102 * Program to display hit probabilities.
103 *
104 * \author mdejong
105 */
106int main(int argc, char **argv)
107{
108 using namespace std;
109 using namespace JPP;
110 using namespace KM3NETDAQ;
111
112 typedef JParallelFileScanner< JTypeList<JDAQEvent, JFIT::JEvt> > JParallelFileScanner_t;
113 typedef JParallelFileScanner_t::multi_pointer_type multi_pointer_type;
115 JACOUSTICS::JEvt>::typelist calibration_types;
116 typedef JMultipleFileScanner<calibration_types> JCalibration_t;
117
118 JParallelFileScanner_t inputFile;
119 JLimit_t& numberOfEvents = inputFile.getLimit();
120 string detectorFile;
121 JCalibration_t calibrationFile;
122 double Tmax_s;
123 string pdfFile;
124 string outputFile;
125 JMuonGandalfParameters_t parameters;
126 int application;
127 JEventSelector event_selector;
128 JCanvas canvas;
129 bool batch;
130 struct : JStyle::JParameters {
131 double arrowSize = 0.003;
132 string arrowType = "|->";
133 double arrowScale = 250.0;
134 Width_t lineWidth = 2;
135 Style_t lineStyle = 1;
136 int nbinsX = 50;
137 int nbinsY = 250;
138 double T_ns = 0.0;
139 } graphics;
140 string option;
141 int debug;
142
143
144 try {
145
146 JProperties properties = graphics.getProperties();
147
148 properties.insert(gmake_property(graphics.arrowSize));
149 properties.insert(gmake_property(graphics.arrowType));
150 properties.insert(gmake_property(graphics.arrowScale));
151 properties.insert(gmake_property(graphics.lineWidth));
152 properties.insert(gmake_property(graphics.lineStyle));
153 properties.insert(gmake_property(graphics.T_ns));
154
155 parameters.numberOfPrefits = 1;
156
157 JParser<> zap("Program to display hit probabilities.");
158
159 zap['w'] = make_field(canvas, "size of canvas <nx>x<ny> [pixels]") = JCanvas(1200, 600);
160 zap['f'] = make_field(inputFile, "input file (output of JXXXMuonReconstruction.sh)");
161 zap['a'] = make_field(detectorFile);
162 zap['+'] = make_field(calibrationFile) = JPARSER::initialised();
163 zap['T'] = make_field(Tmax_s) = 100.0;
164 zap['n'] = make_field(numberOfEvents) = JLimit::max();
165 zap['P'] = make_field(pdfFile);
166 zap['o'] = make_field(outputFile, "graphics output file name") = MAKE_STRING("display_" << WILDCARD << ".gif");
167 zap['@'] = make_field(parameters) = JPARSER::initialised();
168 zap['A'] = make_field(application) = JMUONGANDALF, JMUONENERGY, JMUONSTART;
169 zap['L'] = make_field(event_selector) = JPARSER::initialised();
170 zap['%'] = make_field(properties) = JPARSER::initialised();
171 zap['O'] = make_field(option, "draw option") = arrow_t, histogram_t;
172 zap['B'] = make_field(batch, "batch processing");
173 zap['d'] = make_field(debug) = 1;
174
175 zap(argc, argv);
176 }
177 catch(const exception& error) {
178 FATAL(error.what() << endl);
179 }
180
181 if (batch && outputFile == "") {
182 FATAL("Missing output file name " << outputFile << " in batch mode." << endl);
183 }
184
185 if (!batch && outputFile == "") {
186 outputFile = MAKE_STRING(WILDCARD << ".gif");
187 }
188
189 if (outputFile.find(WILDCARD) == string::npos) {
190 FATAL("Output file name " << outputFile << " has no wild card '" << WILDCARD << "'" << endl);
191 }
192
193
195
196 try {
197 load(detectorFile, detector);
198 }
199 catch(const JException& error) {
200 FATAL(error);
201 }
202 unique_ptr<JDynamics> dynamics;
203
204 try {
205
206 dynamics.reset(new JDynamics(detector, Tmax_s));
207
208 dynamics->load(calibrationFile);
209 }
210 catch(const exception& error) {
211 if (!calibrationFile.empty()) {
212 FATAL(error.what());
213 }
214 }
215
216 const JModuleRouter router(dynamics ? dynamics->getDetector() : detector);
217
218 JSummaryFileRouter summary(inputFile);
219
220 const JMuonPDF_t pdf(pdfFile, parameters.TTS_ns);
221
222 const JTimeRange T_ns(parameters.TMin_ns, parameters.TMax_ns);
223
224 typedef vector<JHitL0> JDataL0_t;
225 typedef vector<JHitW0> JDataW0_t;
226
227 const JBuildL0<JHitL0> buildL0;
228
229
230 Vec offset(0.0, 0.0, 0.0);
231
232 try {
233 offset = getOffset(getHeader(inputFile));
234 } catch(const exception& error) {}
235
236
237 // ROOT
238
239 gROOT->SetBatch(batch);
240
241 TApplication* tp = new TApplication("user", NULL, NULL);
242 TCanvas* cv = new TCanvas("display", "", canvas.x, canvas.y);
243
244 unique_ptr<TStyle> gStyle(new JStyle("gplot", cv->GetWw(), cv->GetWh(), graphics));
245
246 gROOT->SetStyle("gplot");
247 gROOT->ForceStyle();
248
249 const size_t NUMBER_OF_PADS = 3;
250
251 cv->SetFillStyle(4000);
252 cv->SetFillColor(kWhite);
253
254 TPad* p1 = new TPad("p1", NULL, 0.0, 0.00, 1.0, 0.95);
255 TPad* p2 = new TPad("p2", NULL, 0.0, 0.95, 1.0, 1.00);
256
257 p1->Divide(NUMBER_OF_PADS, 1);
258
259 p1->Draw();
260 p2->Draw();
261
262 const double Dmax = getMaximalDistance(detector);
263 const double Rmin = 0.0;
264 const double Rmax = min(parameters.roadWidth_m, 0.4 * Dmax);
265 const double Tmin = min(parameters.TMin_ns, -10.0);
266 const double Tmax = max(parameters.TMax_ns, +100.0);
267 const double Amin = 0.002 * (Tmax - Tmin); // minimal arrow length [ns]
268 const double Amax = 0.8 * (Tmax - Tmin); // maximal arrow length [ns]
269 const double ymin = Tmin - (option == arrow_t ? 0.2 * Amax : 0.0);
270 const double ymax = Tmax + (option == arrow_t ? 0.5 * Amax : 0.0);
271
272 const string Xlabel[NUMBER_OF_PADS] = { "R [m]", "#phi [rad]", "z [m]" };
273 const double Xmin [NUMBER_OF_PADS] = { Rmin, -PI, -0.3 * Dmax };
274 const double Xmax [NUMBER_OF_PADS] = { Rmax, +PI, +0.3 * Dmax };
275
276 double Xs[NUMBER_OF_PADS];
277
278 for (size_t i = 0; i != NUMBER_OF_PADS; ++i) {
279 Xs[i] = 0.003 * (Xmax[i] - Xmin[i]) * (0.5 * NUMBER_OF_PMTS); // x-offset arrow as function of PMT number
280 }
281
282 TH2D H2[NUMBER_OF_PADS];
283 TGraph G2[NUMBER_OF_PADS];
284
285 for (size_t i = 0; i != NUMBER_OF_PADS; ++i) {
286
287 H2[i] = TH2D(MAKE_CSTRING("h" << i), NULL, graphics.nbinsX, Xmin[i] - Xs[i], Xmax[i] + Xs[i], graphics.nbinsY, ymin, ymax);
288
289 H2[i].GetXaxis()->SetTitle(Xlabel[i].c_str());
290 H2[i].GetYaxis()->SetTitle("#Deltat [ns]");
291
292 H2[i].GetXaxis()->CenterTitle(true);
293 H2[i].GetYaxis()->CenterTitle(true);
294
295 H2[i].SetStats(kFALSE);
296
297 G2[i].Set(2);
298
299 G2[i].SetPoint(0, H2[i].GetXaxis()->GetXmin(), 0.0);
300 G2[i].SetPoint(1, H2[i].GetXaxis()->GetXmax(), 0.0);
301
302 p1->cd(i+1);
303
304 H2[i].Draw("AXIS");
305 G2[i].Draw("SAME");
306 }
307
308
309 for (JTreeScanner<Evt> mc(inputFile); inputFile.hasNext(); ) {
310
311 cout << "event: " << setw(8) << inputFile.getCounter() << endl;
312
313 multi_pointer_type ps = inputFile.next();
314
315 JDAQEvent* tev = ps;
316 JFIT::JEvt* in = ps;
317 Evt* event = NULL;
318
319 if (dynamics) {
320 dynamics->update(*tev);
321 }
322
323 if (mc.getEntries() != 0) {
324 event = mc.getEntry(tev->getCounter()); // Monte Carlo true information
325 }
326
327 in->select(JHistory::is_application(application));
328
329 if (!in->empty()) {
330
331 sort(in->begin(), in->end(), qualitySorter);
332
333 if (!event_selector(*tev, *in, event)) {
334 continue;
335 }
336
337
338 JDataL0_t dataL0;
339
340 buildL0(*tev, router, true, back_inserter(dataL0));
341
342 summary.update(*tev);
343
344 JFIT::JFit muon; // Monte Carlo true muon
345
346 if (event != NULL) {
347
348 const time_converter converter = time_converter(*event, *tev);
349
350 for (const auto& t1 : event->mc_trks) {
351 if (is_muon(t1)) {
352 if (t1.E > muon.getE()) {
353
354 JTrack3E ta = getTrack(t1);
355
356 ta.add(getPosition(offset));
357 ta.add(converter.putTime());
358
359 muon = getFit(0, ta, 0.0, 0, t1.E, 1);
360
361 muon.setW(JSTART_LENGTH_METRES, fabs(t1.len));
362 }
363 }
364 }
365 }
366
367 bool monte_carlo = false; // show Monte Carlo true muon
368 size_t index = 0; // index of fit
369
370 for (bool next = false; !next; ) {
371
372 for (size_t i = 0; i != NUMBER_OF_PADS; ++i) {
373 H2[i].Reset();
374 }
375
376 JFIT::JFit fit;
377
378 if (!monte_carlo)
379 fit = (*in)[index];
380 else
381 fit = muon;
382
383 JRotation3D R (getDirection(fit));
384 JLine1Z tz(getPosition (fit).rotate(R), fit.getT());
385 JRange<double> Z_m;
386 /*
387 if (fit.hasW(JSTART_LENGTH_METRES)) {
388 Z_m = JRange<double>(0.0, fit.getW(JSTART_LENGTH_METRES)) + JRange<double>(parameters.ZMin_m, parameters.ZMax_m);
389 }
390 */
391 const JFIT::JModel<JLine1Z> match(tz, parameters.roadWidth_m, T_ns, Z_m);
392
393 // hit selection based on fit result
394
395 JDataW0_t data;
396
397 for (JDataL0_t::const_iterator i = dataL0.begin(); i != dataL0.end(); ++i) {
398
399 JHitW0 hit(*i, summary.getRate(i->getPMTIdentifier(), parameters.R_Hz));
400
401 hit.rotate(R);
402
403 if (match(hit)) {
404 data.push_back(hit);
405 }
406 }
407
408 // select first hit in PMT
409
410 sort(data.begin(), data.end(), JHitW0::compare);
411
412 JDataW0_t::iterator __end = unique(data.begin(), data.end(), equal_to<JDAQPMTIdentifier>());
413
414 double E_GeV = parameters.E_GeV;
415 /*
416 if (fit.getE() > 0.1) {
417 E_GeV = fit.getE();
418 }
419 */
420
421 // move fit to geometrical center of hits
422
424
425 for (JDataW0_t::iterator hit = data.begin(); hit != __end; ++hit) {
426
427 const double x = hit->getX() - tz.getX();
428 const double y = hit->getY() - tz.getY();
429 const double z = hit->getZ();
430 const double R = sqrt(x*x + y*y);
431
432 zs.include(z - R/getTanThetaC());
433 }
434
435 const double z0 = tz.getZ();
436 const double z1 = 0.5 * (zs.getLowerLimit() + zs.getUpperLimit());
437
438 tz.setZ(z1, getSpeedOfLight());
439
440 // graphics
441
442 ostringstream os;
443 vector<TArrow> arrow [NUMBER_OF_PADS];
444 vector<TMarker> marker[NUMBER_OF_PADS];
445
446 if (fit.hasW(JSTART_LENGTH_METRES) && fit.getW(JSTART_LENGTH_METRES) > 0.0) {
447
448 marker[2].push_back(TMarker(z0 - tz.getZ(), 0.0, kFullCircle));
449 marker[2].push_back(TMarker(z0 - tz.getZ() + fit.getW(JSTART_LENGTH_METRES), 0.0, kFullCircle));
450
451 static_cast<TAttMarker&>(marker[2][0]) = TAttMarker(kRed, kFullCircle, 0.7);
452 static_cast<TAttMarker&>(marker[2][1]) = TAttMarker(kRed, kFullCircle, 0.7);
453 }
454
455 DEBUG("trk: "
456 << FIXED(7,2) << tz.getX() << ' '
457 << FIXED(7,2) << tz.getY() << ' '
458 << FIXED(7,2) << tz.getZ() << ' '
459 << FIXED(12,2) << tz.getT() << endl);
460
461 double chi2 = 0;
462
463 for (JDataW0_t::const_iterator hit = data.begin(); hit != __end; ++hit) {
464
465 const double x = hit->getX() - tz.getX();
466 const double y = hit->getY() - tz.getY();
467 const double z = hit->getZ() - tz.getZ();
468 const double R = sqrt(x*x + y*y);
469
470 const double t1 = tz.getT() + (z + R * getTanThetaC()) * getInverseSpeedOfLight();
471
472 JDirection3D dir(hit->getDX(), hit->getDY(), hit->getDZ()); // PMT orientation
473
474 dir.rotate(JRotation3Z(-atan2(y,x))); // rotate PMT axis to x-z plane
475
476 const double theta = dir.getTheta();
477 const double phi = fabs(dir.getPhi()); // rotational symmetry of Cherenkov cone
478
479 //const double E = gWater.getE(E_GeV, z); // correct for energy loss
480 const double E = E_GeV;
481 const double dt = T_ns.constrain(hit->getT() - t1);
482
483 JMuonPDF_t::result_type H1 = pdf.calculate(E, R, theta, phi, dt);
484 JMuonPDF_t::result_type H0(hit->getR() * 1e-9, 0.0, T_ns);
485
486 if (H1.V >= parameters.VMax_npe) {
487 H1 *= parameters.VMax_npe / H1.V;
488 }
489
490 H1 += H0; // signal + background
491
492 chi2 += H1.getChi2() - H0.getChi2();
493
494 DEBUG("hit: "
495 << setw(8) << hit->getModuleID() << '.' << FILL(2,'0') << (int) hit->getPMTAddress() << FILL() << ' '
496 << SCIENTIFIC(8,2) << E << ' '
497 << FIXED(7,2) << R << ' '
498 << FIXED(7,4) << theta << ' '
499 << FIXED(7,4) << phi << ' '
500 << FIXED(7,3) << dt << ' '
501 << FIXED(7,3) << H1.getChi2() << ' '
502 << FIXED(7,3) << H0.getChi2() << endl);
503
504 const double derivative = H1.getDerivativeOfChi2() - H0.getDerivativeOfChi2();
505
506 double size = derivative * graphics.arrowScale; // size of arrow
507
508 if (fabs(size) < Amin) {
509 size = (size > 0.0 ? +Amin : -Amin);
510 } else if (fabs(size) > Amax) {
511 size = (size > 0.0 ? +Amax : -Amax);
512 }
513
514 const double X[NUMBER_OF_PADS] = { R, atan2(y,x), z - R/getTanThetaC() };
515
516 const double xs = (double) (NUMBER_OF_PMTS - 2 * hit->getPMTAddress()) / (double) NUMBER_OF_PMTS;
517
518 for (size_t i = 0; i != NUMBER_OF_PADS; ++i) {
519
520 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());
521
522 a1.SetLineWidth(graphics.lineWidth);
523 a1.SetLineStyle(graphics.lineStyle);
524
525 arrow[i].push_back(a1);
526
527 H2[i].Fill(X[i], dt + graphics.T_ns);
528 }
529 }
530
531 os << FILL(6,'0') << tev->getRunNumber() << ":" << tev->getFrameIndex() << "/" << tev->getCounter() << FILL();
532 os << " Q = " << FIXED(4,0) << fit.getQ()
533 << '/' << FIXED(4,0) << -chi2;
534 os << " E = " << SCIENTIFIC(7,1) << fit.getE() << " [GeV]";
535 os << " cos(#theta) = " << FIXED(6,3) << fit.getDZ();
536
537 if (monte_carlo)
538 os << " Monte Carlo";
539 else if (muon.getStatus() >= 0)
540 os << " #Delta#alpha = " << FIXED(6,2) << getAngle(getDirection(muon), getDirection(fit)) << " [deg]";
541
542
543 // draw
544
545 TLatex title(0.05, 0.5, os.str().c_str());
546
547 title.SetTextAlign(12);
548 title.SetTextFont(42);
549 title.SetTextSize(0.6);
550
551 p2->cd();
552
553 title.Draw();
554
555 for (int i = 0; i != NUMBER_OF_PADS; ++i) {
556
557 p1->cd(i+1);
558
559 if (option == arrow_t) {
560
561 for (auto& a1 : arrow[i]) {
562 a1.Draw();
563 }
564
565 for (auto& m1 : marker[i]) {
566 m1.Draw();
567 }
568 }
569
570 if (option == histogram_t) {
571 H2[i].Draw("SAME");
572 }
573 }
574
575 cv->Update();
576
577
578 // action
579
580 if (batch) {
581
582 cv->SaveAs(replace(outputFile, WILDCARD, MAKE_STRING(inputFile.getCounter())).c_str());
583
584 next = true;
585
586 } else {
587
588 static int count = 0;
589
590 if (count++ == 0) {
591 cout << endl << "Type '?' for possible options." << endl;
592 }
593
594 for (bool user = true; user; ) {
595
596 cout << "\n> " << flush;
597
598 switch (JKeypress(true).get()) {
599
600 case '?':
601 cout << endl;
602 cout << "possible options: " << endl;
603 cout << 'q' << " -> " << "exit application" << endl;
604 cout << 'u' << " -> " << "update canvas" << endl;
605 cout << 's' << " -> " << "save graphics to file" << endl;
606 cout << '+' << " -> " << "next fit" << endl;
607 cout << '-' << " -> " << "previous fit" << endl;
608 cout << 'M' << " -> " << "Monte Carlo true muon information" << endl;
609 cout << 'F' << " -> " << "fit information" << endl;
610 if (event_selector.is_valid()) {
611 cout << 'L' << " -> " << "reload event selector" << endl;
612 }
613 cout << 'r' << " -> " << "rewind input file" << endl;
614 cout << 'R' << " -> " << "switch to ROOT mode (quit ROOT to continue)" << endl;
615 cout << ' ' << " -> " << "next event (as well as any other key)" << endl;
616 break;
617
618 case 'q':
619 cout << endl;
620 return 0;
621
622 case 'u':
623 cv->Update();
624 break;
625
626 case 's':
627 cv->SaveAs(replace(outputFile, WILDCARD, MAKE_STRING(inputFile.getCounter())).c_str());
628 break;
629
630 case '+':
631 monte_carlo = false;
632 index = (index != in->size() - 1 ? index + 1 : 0);
633 user = false;
634 break;
635
636 case '-':
637 monte_carlo = false;
638 index = (index != 0 ? index - 1 : in->size() - 1);
639 user = false;
640 break;
641
642 case 'M':
643 if (muon.getStatus() >= 0)
644 monte_carlo = true;
645 else
646 ERROR(endl << "No Monte Carlo muon available." << endl);
647 user = false;
648 break;
649
650 case 'F':
651 monte_carlo = false;
652 user = false;
653 break;
654
655 case 'L':
656 if (event_selector.is_valid()) {
657 execute(MAKE_STRING("make -f " << getPath(argv[0]) << "/JMakeEventSelector libs"), 3);
658 event_selector.reload();
659 }
660 break;
661
662 case 'R':
663 tp->Run(kTRUE);
664 break;
665
666 case 'r':
667 inputFile.rewind();
668
669 default:
670 next = true;
671 user = false;
672 break;
673 }
674 }
675 }
676 }
677 }
678 }
679 cout << endl;
680}
Definition of hit and track types and auxiliary methods for handling Monte Carlo data.
string outputFile
KM3NeT DAQ constants, bit handling, etc.
Data structure for detector geometry and calibration.
TPaveText * p1
Dynamic detector calibration.
Basic data structure for L0 hit.
Keyboard settings for unbuffered input.
Auxiliary methods for geometrical 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:72
Direct access to module in detector data structure.
Auxiliary data structure for muon PDF.
Parallel scanning of objects from a single file or multiple files according a format that follows fro...
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
int main(int argc, char **argv)
ROOT TTree parameter settings of various packages.
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
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
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.
JTime & add(const JTime &value)
Addition operator.
3D track with energy.
Definition JTrack3E.hh:32
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
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.
File router for fast addressing of summary data.
void update(const JDAQHeader &header)
Update router.
double getRate(const JDAQPMTIdentifier &id, const double rate_Hz) const
Get rate.
Template definition for direct access of elements in ROOT TChain.
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
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 ...
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::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).
bool qualitySorter(const JFit &first, const JFit &second)
Comparison of fit results.
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
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:233
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
Auxiliary include file for time conversion between DAQ/trigger hit and Monte Carlo hit.