Jpp  15.0.0
the software that should make you happy
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 
7 #include "TROOT.h"
8 #include "TApplication.h"
9 #include "TCanvas.h"
10 #include "TStyle.h"
11 #include "TH2D.h"
12 #include "TArrow.h"
13 #include "TLatex.h"
14 
19 
20 #include "JROOT/JStyle.hh"
21 #include "JROOT/JCanvas.hh"
22 #include "JROOT/JRootToolkit.hh"
23 
24 #include "JDAQ/JDAQEventIO.hh"
26 
27 #include "JDetector/JDetector.hh"
30 
31 #include "JTrigger/JHitL0.hh"
32 #include "JTrigger/JBuildL0.hh"
33 
34 #include "JSupport/JSupport.hh"
38 #include "JAAnet/JAAnetToolkit.hh"
39 
40 #include "JPhysics/JPDF_t.hh"
41 
42 #include "JFit/JLine1Z.hh"
43 #include "JFit/JModel.hh"
45 #include "JReconstruction/JEvt.hh"
48 
49 #include "JLang/JPredicate.hh"
50 #include "JLang/JComparator.hh"
51 
52 #include "JMath/JMathToolkit.hh"
53 #include "JSystem/JKeypress.hh"
54 #include "JSystem/JProcess.hh"
55 
56 #include "Jeep/JFunctionAdaptor.hh"
57 #include "Jeep/JPrint.hh"
58 #include "Jeep/JParser.hh"
59 #include "Jeep/JMessage.hh"
60 
61 
62 namespace {
63 
64  /**
65  * Wild card character for file name substition.
66  */
67  const char WILDCARD = '%';
68 
69 
72 
73 
74  /**
75  * Event selector.
76  *
77  * The default constructor will accept all events.\n
78  * A different method can dynamically be loaded from a (shared) library using class JEEP::JFunctionAdaptor.
79  * For the definition of an alternative method, see e.g.\ event_selector.cc
80  */
81  struct JEventSelector :
82  public JFunctionAdaptor<bool, const JEvt&, const Evt*>
83  {
84  typedef JFunctionAdaptor<bool, const JEvt&, const Evt*> function_adaptor_type;
85 
86  /**
87  * Default event selection.
88  *
89  * \param in input event
90  * \param event pointer to Monte Carlo event
91  * \return true
92  */
93  static inline bool accept(const JEvt& in, const Evt* event)
94  {
95  return true;
96  }
97 
98 
99  /**
100  * Default constructor.
101  */
102  JEventSelector() :
103  JFunctionAdaptor<bool, const JEvt&, const Evt*>(accept)
104  {}
105 
106 
107  /**
108  * Check validity of function.
109  *
110  * \return true if valid function; else false
111  */
112  bool is_valid() const
113  {
114  return (libso != "" && function_adaptor_type::is_valid());
115  }
116 
117 
118  /**
119  * Re-load function from shared library.
120  */
121  void reload()
122  {
123  reset();
124 
125  load(libso);
126  }
127 
128 
129  /**
130  * Read event selector from input stream.
131  *
132  * \param in input stream
133  * \param object event selector
134  * \return input stream
135  */
136  friend std::istream& operator>>(std::istream& in, JEventSelector& object)
137  {
138  using namespace std;
139 
140  object.reset();
141 
142  if (in >> object.libso) {
143 
144  if (object.libso.find(function_adaptor_type::SEPARATOR) == string::npos) {
145  object.libso += MAKE_CSTRING(function_adaptor_type::SEPARATOR << "accept");
146  }
147 
148  object.load(object.libso);
149  }
150 
151  return in;
152  }
153 
154  private:
155  std::string libso; //!< shared library : symbol
156  };
157 
158 
159  /**
160  * Execute command in shell.
161  *
162  * \param command command
163  */
164  inline void execute(const std::string& command, int debug)
165  {
166  using namespace std;
167  using namespace JPP;
168 
169  JProcess process(command);
170 
171  istream in(process.getInputStreamBuffer());
172 
173  for (string buffer; getline(in, buffer); ) {
174  DEBUG(buffer << endl);
175  }
176  }
177 }
178 
179 
180 /**
181  * \file
182  *
183  * Program to display hit probabilities.
184  *
185  * \author mdejong
186  */
187 int main(int argc, char **argv)
188 {
189  using namespace std;
190  using namespace JPP;
191  using namespace KM3NETDAQ;
192 
193  typedef JParallelFileScanner< JTypeList<JDAQEvent, JEvt> > JParallelFileScanner_t;
194  typedef JParallelFileScanner_t::multi_pointer_type multi_pointer_type;
195 
196  JParallelFileScanner_t inputFile;
197  JLimit_t& numberOfEvents = inputFile.getLimit();
198  string detectorFile;
199  string pdfFile;
200  string outputFile;
202  int application;
203  JEventSelector event_selector;
204  JCanvas canvas;
205  bool batch;
206  double arrowSize;
207  string arrowType;
208  double arrowScale;
209  int debug;
210 
211 
212  try {
213 
214  parameters.numberOfPrefits = 1;
215 
216  JParser<> zap("Program to display hit probabilities.");
217 
218  zap['w'] = make_field(canvas, "size of canvas <nx>x<ny> [pixels]") = JCanvas(1200, 600);
219  zap['f'] = make_field(inputFile, "input file (output of JXXXMuonReconstruction.sh)");
220  zap['a'] = make_field(detectorFile);
221  zap['n'] = make_field(numberOfEvents) = JLimit::max();
222  zap['P'] = make_field(pdfFile);
223  zap['o'] = make_field(outputFile, "graphics output file name") = MAKE_STRING("display_" << WILDCARD << ".gif");
225  zap['A'] = make_field(application) = JMUONGANDALF, JMUONENERGY;
226  zap['L'] = make_field(event_selector) = JPARSER::initialised();
227  zap['S'] = make_field(arrowSize) = 0.003;
228  zap['T'] = make_field(arrowType) = "|->";
229  zap['F'] = make_field(arrowScale) = 250.0;
230  zap['B'] = make_field(batch, "batch processing");
231  zap['d'] = make_field(debug) = 1;
232 
233  zap(argc, argv);
234  }
235  catch(const exception& error) {
236  FATAL(error.what() << endl);
237  }
238 
239  if (batch && outputFile == "") {
240  FATAL("Missing output file name " << outputFile << " in batch mode." << endl);
241  }
242 
243  if (!batch && outputFile == "") {
244  outputFile = MAKE_STRING(WILDCARD << ".gif");
245  }
246 
247  if (outputFile.find(WILDCARD) == string::npos) {
248  FATAL("Output file name " << outputFile << " has no wild card '" << WILDCARD << "'" << endl);
249  }
250 
251 
253 
254  try {
255  load(detectorFile, detector);
256  }
257  catch(const JException& error) {
258  FATAL(error);
259  }
260 
261  const JModuleRouter router(detector);
262 
263  JSummaryFileRouter summary(inputFile, parameters.R_Hz);
264 
265  const JMuonPDF_t pdf(pdfFile, parameters.TTS_ns);
266 
267  const JTimeRange T_ns(parameters.TMin_ns, parameters.TMax_ns);
268 
269  typedef vector<JHitL0> JDataL0_t;
270  typedef vector<JHitW0> JDataW0_t;
271 
272  const JBuildL0<JHitL0> buildL0;
273 
274 
275  // Monte Carlo
276 
277  JVector3D center(0.0, 0.0, 0.0);
278 
279  try {
280  center = get<JPosition3D>(getHeader(inputFile));
281  } catch(const exception& error) {}
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  JSinglePointer<TStyle> gStyle(new JStyle("gplot", cv->GetWw(), cv->GetWh()));
292 
293  gROOT->SetStyle("gplot");
294  gROOT->ForceStyle();
295 
296  const size_t NUMBER_OF_PADS = 3;
297 
298  cv->SetFillStyle(4000);
299  cv->SetFillColor(kWhite);
300  cv->Divide(NUMBER_OF_PADS, 1);
301 
302  const double Dmax = getMaximalDistance(detector);
303  const double Rmin = 0.0;
304  const double Rmax = min(parameters.roadWidth_m, 0.4 * Dmax);
305  const double Tmin = min(parameters.TMin_ns, -10.0);
306  const double Tmax = max(parameters.TMax_ns, +100.0);
307  const double Amin = 0.002 * (Tmax - Tmin); // minimal arrow length [ns]
308  const double Amax = 0.8 * (Tmax - Tmin); // maximal arrow length [ns]
309  const double ymin = min(-Amax, Tmin - 0.3 * Amax);
310  const double ymax = max(+Amax, Tmax + 0.5 * Amax);
311 
312  const string Xlabel[NUMBER_OF_PADS] = { "R [m]", "#phi [rad]", "z [m]" };
313  const double Xmin [NUMBER_OF_PADS] = { Rmin, -PI, -0.5 * Dmax };
314  const double Xmax [NUMBER_OF_PADS] = { Rmax, +PI, +0.5 * Dmax };
315 
316  double Xs[NUMBER_OF_PADS];
317 
318  for (size_t i = 0; i != NUMBER_OF_PADS; ++i) {
319  Xs[i] = 0.003 * (Xmax[i] - Xmin[i]) * (0.5 * NUMBER_OF_PMTS); // x-offset arrow as function of PMT number
320  }
321 
322  TH2D H2[NUMBER_OF_PADS];
323  TGraph G2[NUMBER_OF_PADS];
324 
325  for (size_t i = 0; i != NUMBER_OF_PADS; ++i) {
326 
327  H2[i] = TH2D(MAKE_CSTRING("h" << i), NULL, 100, Xmin[i] - Xs[i], Xmax[i] + Xs[i], 100, ymin, ymax);
328 
329  H2[i].GetXaxis()->SetTitle(Xlabel[i].c_str());
330  H2[i].GetYaxis()->SetTitle("#Deltat [ns]");
331 
332  H2[i].GetXaxis()->CenterTitle(true);
333  H2[i].GetYaxis()->CenterTitle(true);
334 
335  H2[i].SetStats(kFALSE);
336 
337  G2[i].Set(2);
338 
339  G2[i].SetPoint(0, H2[i].GetXaxis()->GetXmin(), 0.0);
340  G2[i].SetPoint(1, H2[i].GetXaxis()->GetXmax(), 0.0);
341 
342  cv->cd(i+1);
343 
344  H2[i].Draw("AXIS");
345  G2[i].Draw("SAME");
346  }
347 
348 
349  for (JTreeScanner<Evt> mc(inputFile); inputFile.hasNext(); ) {
350 
351  cout << "\revent: " << setw(8) << inputFile.getCounter() << flush;
352 
353  multi_pointer_type ps = inputFile.next();
354 
355  JDAQEvent* tev = ps;
356  JEvt* in = ps;
357  Evt* event = NULL;
358 
359  if (mc.getEntries() != 0) {
360  event = mc.getEntry(tev->getCounter()); // Monte Carlo true information
361  }
362 
363 
364  in->select(JHistory::is_event(application));
365 
366  if (!in->empty()) {
367 
368  sort(in->begin(), in->end(), qualitySorter);
369 
370  if (!event_selector(*in, event)) {
371  continue;
372  }
373 
374 
375  JDataL0_t dataL0;
376 
377  buildL0(*tev, router, true, back_inserter(dataL0));
378 
379  summary.update(*tev);
380 
381 
382  JFit muon; // Monte Carlo true muon
383 
384  if (event != NULL) {
385 
386  const time_converter converter = time_converter(*event, *tev);
387 
388  for (const auto& t1 : event->mc_trks) {
389  if (is_muon(t1)) {
390  if (t1.E > muon.getE()) {
391 
392  JTrack3E ta = getTrack(t1);
393 
394  ta.add(center);
395  ta.add(converter.putTime());
396 
397  muon = getFit(0, ta, 0.0, 0, t1.E, 1);
398  }
399  }
400  }
401  }
402 
403 
404  bool next_event = false; // goto next event
405  bool monte_carlo = false; // show Monte Carlo true muon
406  size_t index = 0; // index of fit
407 
408  while (!next_event) {
409 
410  JFit fit;
411 
412  if (!monte_carlo)
413  fit = (*in)[index];
414  else
415  fit = muon;
416 
417 
418  JRotation3D R (getDirection(fit));
419  JLine1Z tz(getPosition (fit).rotate(R), fit.getT());
420  JRange<double> Z_m;
421  /*
422  if (fit.hasW(JSTART_LENGTH_METRES)) {
423  Z_m = JRange<double>(0.0, fit.getW(JSTART_LENGTH_METRES)) + JRange<double>(parameters.ZMin_m, parameters.ZMax_m);
424  }
425  */
426  const JModel<JLine1Z> match(tz, parameters.roadWidth_m, T_ns, Z_m);
427 
428  // hit selection based on fit result
429 
430  JDataW0_t data;
431 
432  for (JDataL0_t::const_iterator i = dataL0.begin(); i != dataL0.end(); ++i) {
433 
434  JHitW0 hit(*i, summary.getRate(i->getPMTIdentifier()));
435 
436  hit.rotate(R);
437 
438  if (match(hit)) {
439  data.push_back(hit);
440  }
441  }
442 
443  // select first hit in PMT
444 
445  sort(data.begin(), data.end(), JHitW0::compare);
446 
447  JDataW0_t::iterator __end = unique(data.begin(), data.end(), equal_to<JDAQPMTIdentifier>());
448 
449  double E_GeV = parameters.E_GeV;
450  /*
451  if (fit.getE() > 0.1) {
452  E_GeV = fit.getE();
453  }
454  */
455 
456  // move fit to geometrical center of hits
457 
458  JRange<double> zs(make_array(data.begin(), __end, &JHitW0::getZ));
459 
460  tz.setZ(0.5 * (zs.getLowerLimit() + zs.getUpperLimit()), getSpeedOfLight());
461 
462 
463  // graphics
464 
465  TLatex latex[NUMBER_OF_PADS];
466  vector<TArrow> arrow[NUMBER_OF_PADS];
467 
468  for (int i = 0; i != NUMBER_OF_PADS; ++i) {
469 
470  latex[i].SetTextAlign(12);
471  latex[i].SetTextFont(42);
472  latex[i].SetTextSize(0.06);
473 
474  latex[i].SetX(H2[i].GetXaxis()->GetXmin() + 0.05 * (H2[i].GetXaxis()->GetXmax() - H2[i].GetXaxis()->GetXmin()));
475  latex[i].SetY(H2[i].GetYaxis()->GetXmax() - 0.05 * (H2[i].GetYaxis()->GetXmax() - H2[i].GetYaxis()->GetXmin()));
476  }
477 
478  double chi2 = 0;
479 
480  for (JDataW0_t::const_iterator hit = data.begin(); hit != __end; ++hit) {
481 
482  const double x = hit->getX() - tz.getX();
483  const double y = hit->getY() - tz.getY();
484  const double z = hit->getZ() - tz.getZ();
485  const double R = sqrt(x*x + y*y);
486 
487  const double t1 = tz.getT() + (z + R * getTanThetaC()) * getInverseSpeedOfLight();
488 
489  JDirection3D dir(hit->getDX(), hit->getDY(), hit->getDZ()); // PMT orientation
490 
491  dir.rotate(JRotation3Z(-atan2(y,x))); // rotate PMT axis to x-z plane
492 
493  const double theta = dir.getTheta();
494  const double phi = fabs(dir.getPhi()); // rotational symmetry of Cherenkov cone
495 
496  //const double E = gWater.getE(E_GeV, z); // correct for energy loss
497  const double E = E_GeV;
498  const double dt = T_ns.constrain(hit->getT() - t1);
499 
500  JMuonPDF_t::result_type H1 = pdf.calculate(E, R, theta, phi, dt);
501  JMuonPDF_t::result_type H0(hit->getR() * 1e-9, 0.0, T_ns);
502 
503  if (H1.V >= parameters.VMax_npe) {
504  H1 *= parameters.VMax_npe / H1.V;
505  }
506 
507  H1 += H0; // signal + background
508 
509  chi2 += H1.getChi2() - H0.getChi2();
510 
511  const double derivative = H1.getDerivativeOfChi2() - H0.getDerivativeOfChi2();
512 
513  double size = derivative * arrowScale; // size of arrow
514 
515  if (fabs(size) < Amin) {
516  size = (size > 0.0 ? +Amin : -Amin);
517  } else if (fabs(size) > Amax) {
518  size = (size > 0.0 ? +Amax : -Amax);
519  }
520 
521  const double X[NUMBER_OF_PADS] = { R, atan2(y,x), z };
522 
523  const double xs = (double) (NUMBER_OF_PMTS - 2 * hit->getPMTAddress()) / (double) NUMBER_OF_PMTS;
524 
525  for (size_t i = 0; i != NUMBER_OF_PADS; ++i) {
526  arrow[i].push_back(TArrow(X[i] + xs*Xs[i], dt, X[i] + xs*Xs[i], dt + size, arrowSize, arrowType.c_str()));
527  }
528  }
529 
530  latex[0].SetTitle(MAKE_CSTRING("" << FILL(6,'0') << tev->getRunNumber() << FILL() <<
531  ":" << tev->getFrameIndex() <<
532  " " << tev->getCounter()));
533 
534  if (!monte_carlo)
535  latex[1].SetTitle(MAKE_CSTRING("[" << index << "]" << " "
536  "Q = " << FIXED(4,0) << getQuality(chi2) << " "
537  "#beta = " << FIXED(6,2) << fit.getW(JGANDALF_BETA0_RAD, 0.0) * 180.0/PI << " [deg]"));
538  else
539  latex[1].SetTitle(MAKE_CSTRING("Q = " << FIXED(4,0) << getQuality(chi2)));
540 
541  if (muon.getStatus() >= 0) {
542 
543  if (!monte_carlo)
544  latex[2].SetTitle(MAKE_CSTRING("#Delta#alpha = " << FIXED(6,2) << getAngle(getDirection(muon), getDirection(fit)) << " [deg]"));
545  else
546  latex[2].SetTitle("Monte Carlo");
547  }
548 
549 
550  // draw
551 
552  for (int i = 0; i != NUMBER_OF_PADS; ++i) {
553 
554  cv->cd(i+1);
555 
556  latex[i].Draw();
557 
558  for (auto& a1 : arrow[i]) {
559  a1.Draw();
560  }
561  }
562 
563  cv->Update();
564 
565 
566  // action
567 
568  if (batch) {
569 
570  cv->SaveAs(replace(outputFile, WILDCARD, MAKE_STRING(inputFile.getCounter())).c_str());
571 
572  next_event = true;
573 
574  } else {
575 
576  for (int c = 0; c == 0; ) {
577 
578  if (!batch) {
579 
580  static bool first = true;
581 
582  if (first) {
583  cout << endl << "Type '?' for possible options." << endl;
584  }
585 
586  first = false;
587  }
588 
589  cout << "\n> " << flush;
590 
591  switch((c = JKeypress(true).get())) {
592 
593  case '?':
594  cout << endl
595  << "possible options: " << endl
596  << 'q' << " -> " << "exit application" << endl
597  << 'u' << " -> " << "update canvas" << endl
598  << 's' << " -> " << "save graphics to file" << endl
599  << '+' << " -> " << "next fit" << endl
600  << '-' << " -> " << "previous fit" << endl
601  << 'M' << " -> " << "Monte Carlo true muon information" << endl
602  << 'F' << " -> " << "fit information" << endl
603  << 'x' << " -> " << "remake standard event selector" << endl
604  << 'L' << " -> " << "reload event selector" << endl
605  << 'r' << " -> " << "rewind input file" << endl
606  << 'R' << " -> " << "switch to ROOT mode (quit ROOT to continue)" << endl
607  << ' ' << " -> " << "next event (and any other key)" << endl;
608  c = 0;
609  break;
610 
611  case 'q':
612  cout << endl;
613  return 0;
614 
615  case 'u':
616  cv->Update();
617  c = 0;
618  break;
619 
620  case 's':
621  cv->SaveAs(replace(outputFile, WILDCARD, MAKE_STRING(inputFile.getCounter())).c_str());
622  c = 0;
623  break;
624 
625  case '+':
626  monte_carlo = false;
627  index = (index != in->size() - 1 ? index + 1 : 0);
628  break;
629 
630  case '-':
631  monte_carlo = false;
632  index = (index != 0 ? index - 1 : in->size() - 1);
633  break;
634 
635  case 'M':
636  if (muon.getStatus() >= 0)
637  monte_carlo = true;
638  else
639  ERROR(endl << "No Monte Carlo muon available." << endl);
640  break;
641 
642  case 'F':
643  monte_carlo = false;
644  break;
645 
646  case 'r':
647  inputFile.rewind();
648  next_event = true;
649  break;
650 
651  case 'x':
652  execute(MAKE_STRING("make -C " << getPath(argv[0]) << "./" << " libs"), debug);
653  c = 0;
654  break;
655 
656  case 'L':
657  if (event_selector.is_valid())
658  event_selector.reload();
659  else
660  ERROR(endl << "No event selector (use command line option -L)." << endl);
661  break;
662 
663  case 'R':
664  tp->Run(kTRUE);
665  c = 0;
666  break;
667 
668  default:
669  next_event = true;
670  break;
671  }
672  }
673  }
674  }
675  }
676  }
677  cout << endl;
678 }
Utility class to parse command line options.
Definition: JParser.hh:1500
double getAngle(const JQuaternion3D &first, const JQuaternion3D &second)
Get space angle between quanternions.
General exception.
Definition: JException.hh:23
T getLowerLimit() const
Get lower limit.
Definition: JRange.hh:202
int main(int argc, char *argv[])
Definition: Main.cc:15
ROOT TTree parameter settings of various packages.
bool accept(const JEvt &in, const Evt *event)
Event selection.
Auxiliary methods for geometrical methods.
TString replace(const TString &target, const TRegexp &regexp, const T &replacement)
Replace regular expression in input by given replacement.
Definition: JPrintResult.cc:63
Data structure for direction in three dimensions.
Definition: JDirection3D.hh:33
double getQuality(const double chi2, const int NDF)
Get quality of fit.
JTrack3E getTrack(const Trk &track)
Get track.
Template specialisation of L0 builder for JHitL0 data type.
Definition: JBuildL0.hh:102
do rm f tmp H1
std::string getPath(const std::string &file_name)
Get path, i.e. part before last JEEP::PATHNAME_SEPARATOR if any.
Definition: JeepToolkit.hh:108
void update(const JDAQHeader &header)
Update router.
Detector data structure.
Definition: JDetector.hh:89
bool is_muon(const Trk &track)
Test whether given track is a (anti-)muon.
Router for direct addressing of module data in detector data structure.
Rotation matrix.
Definition: JRotation3D.hh:111
double getRate() const
Get default rate.
*fatal Wrong number of arguments esac JCookie sh typeset Z DETECTOR typeset Z SOURCE_RUN typeset Z TARGET_RUN set_variable PARAMETERS_FILE $WORKDIR parameters
Definition: diff-Tuna.sh:38
General purpose class for parallel reading of objects from a single file or multiple files...
#define MAKE_CSTRING(A)
Make C-string.
Definition: JPrint.hh:151
Template specialisation of class JModel to match hit with muon trajectory along z-axis.
Definition: JFit/JModel.hh:34
Empty structure for specification of parser element that is initialised (i.e. does not require input)...
Definition: JParser.hh:66
Auxiliary class to convert DAQ hit time to/from Monte Carlo hit time.
Auxiliary data structure for floating point format specification.
Definition: JManip.hh:446
string outputFile
int getRunNumber() const
Get run number.
3D track with energy.
Definition: JTrack3E.hh:30
Data structure for detector geometry and calibration.
then usage E
Definition: JMuonPostfit.sh:35
void load(const std::string &file_name, const std::string &symbol)
Load function from shared library.
Acoustic single fit.
bool qualitySorter(const JFit &first, const JFit &second)
Comparison of fit results.
then echo The file $DIR KM3NeT_00000001_00000000 root already please rename or remove it first
#define MAKE_STRING(A)
Make string.
Definition: JPrint.hh:142
int getFrameIndex() const
Get frame index.
JTime & add(const JTime &value)
Addition operator.
Head getHeader(const JMultipleFileScanner_t &file_list)
Get Monte Carlo header.
Basic data structure for L0 hit.
Rotation around Z-axis.
Definition: JRotation3D.hh:85
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
The template JSinglePointer class can be used to hold a pointer to an object.
Auxiliary class for defining the range of iterations of objects.
Definition: JLimit.hh:41
double getMaximalDistance(const JDetector &detector)
Get maximal distance between modules in detector.
I/O formatting auxiliaries.
JAxis3D & rotate(const JRotation3D &R)
Rotate axis.
Definition: JAxis3D.hh:225
Detector file.
Definition: JHead.hh:196
JDirection3D getDirection(const Vec &dir)
Get direction.
JFunction1D_t::result_type result_type
Definition: JPDF_t.hh:145
Data structure for vector in three dimensions.
Definition: JVector3D.hh:34
JDirection3D & rotate(const JRotation3D &R)
Rotate.
Acoustic event fit.
Keyboard settings for unbuffered input.
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:1961
Enable unbuffered terminal input.
Definition: JKeypress.hh:32
const array_type< JValue_t > & make_array(const JValue_t(&array)[N])
Method to create array of values.
Definition: JVectorize.hh:54
Auxiliary class to test history.
Definition: JHistory.hh:104
static const int JMUONGANDALF
double getTheta() const
Get theta angle.
Definition: JVersor3D.hh:128
Parallel scanning of objects from a single file or multiple files according a format that follows fro...
bool is_valid(const json &js)
Check validity of JSon data.
JPosition3D getPosition(const Vec &pos)
Get position.
T getUpperLimit() const
Get upper limit.
Definition: JRange.hh:213
Definition of hit and track types and auxiliary methods for handling Monte Carlo data.
#define ERROR(A)
Definition: JMessage.hh:66
then break fi done getCenter read X Y Z let X
std::istream & getline(std::istream &in, JString &object)
Read string from input stream until end of line.
Definition: JString.hh:478
static const int JGANDALF_BETA0_RAD
KM3NeT Data Definitions v2.1.0-7-g97845ea https://git.km3net.de/common/km3net-dataformat.
double putTime() const
Get Monte Carlo time minus DAQ/trigger time.
static const double PI
Mathematical constants.
T constrain(argument_type x) const
Constrain value to range.
Definition: JRange.hh:350
File router for fast addressing of summary data.
int debug
debug level
Definition: JSirene.cc:63
Auxiliary data structure for muon PDF.
then usage $script[distance] fi case set_variable R
Definition: JDrawLED.sh:43
General purpose messaging.
Auxiliary data structure for sequence of same character.
Definition: JManip.hh:328
Auxiliary include file for time conversion between DAQ/trigger hit and Monte Carlo hit...
Auxiliary data structure for muon PDF.
Definition: JPDF_t.hh:135
#define FATAL(A)
Definition: JMessage.hh:67
Direct access to module in detector data structure.
Streaming of input and output from Linux command.
Definition: JProcess.hh:29
JFit getFit(const JHistory &history, const JTrack3D &track, const double Q, const int NDF, const double energy=0.0, const int status=0)
Get fit.
Auxiliary class for a hit with background rate value.
Definition: JHitW0.hh:21
const double getSpeedOfLight()
Get speed of light.
void load(const std::string &file_name, JDetector &detector)
Load detector from input file.
Utility class to parse command line options.
const double getInverseSpeedOfLight()
Get inverse speed of light.
then usage $script< input_file >< detector_file > fi set_variable OUTPUT_DIR set_variable SELECTOR JDAQTimesliceL1 set_variable DEBUG case set_variable DEBUG
Data structure for fit of straight line paralel to z-axis.
Definition: JLine1Z.hh:27
double getTanThetaC()
Get average tangent of Cherenkov angle of water corresponding to group velocity.
const JLimit & getLimit() const
Get limit.
Definition: JLimit.hh:73
do set_variable DETECTOR_TXT $WORKDIR detector
KM3NeT DAQ constants, bit handling, etc.
static const int NUMBER_OF_PMTS
Total number of PMTs in module.
Definition: JDAQ.hh:26
Function adaptor.
then fatal Wrong number of arguments fi set_variable DETECTOR $argv[1] set_variable INPUT_FILE $argv[2] eval JPrintDetector a $DETECTOR O IDENTIFIER eval JPrintDetector a $DETECTOR O SUMMARY source JAcoustics sh $DETECTOR_ID CHECK_EXIT_CODE typeset A TRIPODS get_tripods $WORKDIR tripod txt TRIPODS for EMITTER in
Definition: JCanberra.sh:41
JTriggerCounter_t getCounter() const
Get trigger counter.
Wrapper class around ROOT TStyle.
Definition: JStyle.hh:20
static const int JMUONENERGY
The Evt class respresent a Monte Carlo (MC) event as well as an offline event.
Definition: Evt.hh:19
Data structure for size of TCanvas.
Definition: JCanvas.hh:26
friend std::istream & operator>>(std::istream &in, JFunctionAdaptorHelper &object)
Read function adaptor helper from input stream.