Jpp  18.5.2
the software that should make you happy
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 "TStyle.h"
12 #include "TH2D.h"
13 #include "TArrow.h"
14 #include "TLatex.h"
15 #include "TMarker.h"
16 
22 
23 #include "JROOT/JStyle.hh"
24 #include "JROOT/JCanvas.hh"
25 #include "JROOT/JRootToolkit.hh"
26 
28 
29 #include "JTrigger/JHitL0.hh"
30 
31 #include "JSupport/JSupport.hh"
35 #include "JAAnet/JAAnetToolkit.hh"
36 
37 #include "JPhysics/JPDF_t.hh"
38 
39 #include "JFit/JLine1Z.hh"
40 #include "JFit/JModel.hh"
43 
44 #include "JLang/JPredicate.hh"
45 #include "JLang/JComparator.hh"
46 
47 #include "JMath/JMathToolkit.hh"
48 #include "JSystem/JKeypress.hh"
49 #include "JSystem/JProcess.hh"
50 
51 #include "Jeep/JFunctionAdaptor.hh"
52 #include "Jeep/JPrint.hh"
53 #include "Jeep/JParser.hh"
54 #include "Jeep/JMessage.hh"
55 
56 
57 namespace JAANET {
58 
60 
61  /**
62  * Event selector.
63  *
64  * The default constructor will accept all events.\n
65  * A different method can dynamically be loaded from a shared library using class JEEP::JFunctionAdaptor.
66  */
67  struct JEventSelector :
68  public JFunctionAdaptor<bool, const Trk&, const Evt&>
69  {
70  /**
71  * Default event selection.
72  *
73  * \param trk track
74  * \param evt event
75  * \return true
76  */
77  static inline bool select(const Trk& trk, const Evt& evt)
78  {
79  return true;
80  }
81 
82 
83  /**
84  * Default constructor.
85  */
87  {
88  this->function = select;
89  this->symbol = "select";
90  }
91  };
92 
93 
94  /**
95  * Check availability of value.
96  *
97  * \param trk track
98  * \param i index
99  * \return true if available; else false
100  */
101  inline bool hasW(const Trk& trk, const int i)
102  {
103  return (i >= 0 && i < (int) trk.fitinf.size());
104  }
105 
106 
107  /**
108  * Get associated value.
109  *
110  * \param trk track
111  * \param i index
112  * \return value
113  */
114  inline double getW(const Trk& trk, const int i)
115  {
116  return trk.fitinf.at(i);
117  }
118 
119 
120  /**
121  * Get associated value.
122  *
123  * \param trk track
124  * \param i index
125  * \param value default value
126  * \return value
127  */
128  inline double getW(const Trk& trk, const int i, const double value)
129  {
130  if (hasW(trk,i))
131  return trk.fitinf.at(i);
132  else
133  return value;
134  }
135 
136 
137  /**
138  * Set associated value.
139  *
140  * \param trk track
141  * \param i index
142  * \param value value
143  */
144  void setW(Trk& trk, const int i, const double value)
145  {
146  if (i >= (int) trk.fitinf.size()) {
147  trk.fitinf.resize(i + 1, 0.0);
148  }
149 
150  trk.fitinf[i] = value;
151  }
152 }
153 
154 namespace {
155 
156  /**
157  * Wild card character for file name substition.
158  */
159  const char WILDCARD = '%';
160 
161 
162  /**
163  * Execute command in shell.
164  *
165  * \param command command
166  */
167  inline void execute(const std::string& command, int debug)
168  {
169  using namespace std;
170  using namespace JPP;
171 
172  JProcess process(command);
173 
174  istream in(process.getInputStreamBuffer());
175 
176  for (string buffer; getline(in, buffer); ) {
177  DEBUG(buffer << endl);
178  }
179  }
180 }
181 
182 
183 /**
184  * \file
185  *
186  * Program to display hit probabilities.
187  *
188  * \author mdejong
189  */
190 int main(int argc, char **argv)
191 {
192  using namespace std;
193  using namespace JPP;
194  using namespace KM3NETDAQ;
195 
196  JSingleFileScanner<Evt> inputFile;
197  JLimit_t& numberOfEvents = inputFile.getLimit();
198  string pdfFile;
199  string outputFile;
201  int application;
202  JEventSelector event_selector;
203  JCanvas canvas;
204  bool batch;
205  double arrowSize;
206  string arrowType;
207  double arrowScale;
208  int debug;
209 
210 
211  try {
212 
213  parameters.numberOfPrefits = 1;
214 
215  JParser<> zap("Program to display hit probabilities.");
216 
217  zap['w'] = make_field(canvas, "size of canvas <nx>x<ny> [pixels]") = JCanvas(1200, 600);
218  zap['f'] = make_field(inputFile, "input file (output of JXXXMuonReconstruction.sh)");
219  zap['n'] = make_field(numberOfEvents) = JLimit::max();
220  zap['P'] = make_field(pdfFile);
221  zap['o'] = make_field(outputFile, "graphics output file name") = MAKE_STRING("display_" << WILDCARD << ".gif");
224  zap['L'] = make_field(event_selector) = JPARSER::initialised();
225  zap['S'] = make_field(arrowSize) = 0.003;
226  zap['T'] = make_field(arrowType) = "|->";
227  zap['F'] = make_field(arrowScale) = 250.0;
228  zap['B'] = make_field(batch, "batch processing");
229  zap['d'] = make_field(debug) = 1;
230 
231  zap(argc, argv);
232  }
233  catch(const exception& error) {
234  FATAL(error.what() << endl);
235  }
236 
237  if (batch && outputFile == "") {
238  FATAL("Missing output file name " << outputFile << " in batch mode." << endl);
239  }
240 
241  if (!batch && outputFile == "") {
242  outputFile = MAKE_STRING(WILDCARD << ".gif");
243  }
244 
245  if (outputFile.find(WILDCARD) == string::npos) {
246  FATAL("Output file name " << outputFile << " has no wild card '" << WILDCARD << "'" << endl);
247  }
248 
249 
250  JSummaryFileRouter summary(inputFile, parameters.R_Hz);
251 
252  const JMuonPDF_t pdf(pdfFile, parameters.TTS_ns);
253 
254  const JTimeRange T_ns(parameters.TMin_ns, parameters.TMax_ns);
255 
256  typedef vector<JHitL0> JDataL0_t;
257  typedef vector<JHitW0> JDataW0_t;
258 
259 
260  // ROOT
261 
262  gROOT->SetBatch(batch);
263 
264  TApplication* tp = new TApplication("user", NULL, NULL);
265  TCanvas* cv = new TCanvas("display", "", canvas.x, canvas.y);
266 
267  unique_ptr<TStyle> gStyle(new JStyle("gplot", cv->GetWw(), cv->GetWh()));
268 
269  gROOT->SetStyle("gplot");
270  gROOT->ForceStyle();
271 
272  const size_t NUMBER_OF_PADS = 3;
273 
274  cv->SetFillStyle(4000);
275  cv->SetFillColor(kWhite);
276 
277  TPad* p1 = new TPad("p1", NULL, 0.0, 0.00, 1.0, 0.95);
278  TPad* p2 = new TPad("p2", NULL, 0.0, 0.95, 1.0, 1.00);
279 
280  p1->Divide(NUMBER_OF_PADS, 1);
281 
282  p1->Draw();
283  p2->Draw();
284 
285  const double Dmax = 1000.0;
286  const double Rmin = 0.0;
287  const double Rmax = min(parameters.roadWidth_m, 0.4 * Dmax);
288  const double Tmin = min(parameters.TMin_ns, -10.0);
289  const double Tmax = max(parameters.TMax_ns, +100.0);
290  const double Amin = 0.002 * (Tmax - Tmin); // minimal arrow length [ns]
291  const double Amax = 0.8 * (Tmax - Tmin); // maximal arrow length [ns]
292  const double ymin = min(-Amax, Tmin - 0.3 * Amax);
293  const double ymax = max(+Amax, Tmax + 0.5 * Amax);
294 
295  const string Xlabel[NUMBER_OF_PADS] = { "R [m]", "#phi [rad]", "z [m]" };
296  const double Xmin [NUMBER_OF_PADS] = { Rmin, -PI, -0.5 * Dmax };
297  const double Xmax [NUMBER_OF_PADS] = { Rmax, +PI, +0.5 * Dmax };
298 
299  double Xs[NUMBER_OF_PADS];
300 
301  for (size_t i = 0; i != NUMBER_OF_PADS; ++i) {
302  Xs[i] = 0.003 * (Xmax[i] - Xmin[i]) * (0.5 * NUMBER_OF_PMTS); // x-offset arrow as function of PMT number
303  }
304 
305  TH2D H2[NUMBER_OF_PADS];
306  TGraph G2[NUMBER_OF_PADS];
307 
308  for (size_t i = 0; i != NUMBER_OF_PADS; ++i) {
309 
310  H2[i] = TH2D(MAKE_CSTRING("h" << i), NULL, 100, Xmin[i] - Xs[i], Xmax[i] + Xs[i], 100, ymin, ymax);
311 
312  H2[i].GetXaxis()->SetTitle(Xlabel[i].c_str());
313  H2[i].GetYaxis()->SetTitle("#Deltat [ns]");
314 
315  H2[i].GetXaxis()->CenterTitle(true);
316  H2[i].GetYaxis()->CenterTitle(true);
317 
318  H2[i].SetStats(kFALSE);
319 
320  G2[i].Set(2);
321 
322  G2[i].SetPoint(0, H2[i].GetXaxis()->GetXmin(), 0.0);
323  G2[i].SetPoint(1, H2[i].GetXaxis()->GetXmax(), 0.0);
324 
325  p1->cd(i+1);
326 
327  H2[i].Draw("AXIS");
328  G2[i].Draw("SAME");
329  }
330 
331 
332  while (inputFile.hasNext()) {
333 
334  cout << "\revent: " << setw(8) << inputFile.getCounter() << flush;
335 
336  Evt* evt = inputFile.next();
337 
338  if (has_reconstructed_track<JPP_RECONSTRUCTION_TYPE>(*evt, rec_stages_range(application))) {
339 
340  Trk fit = get_best_reconstructed_track<JPP_RECONSTRUCTION_TYPE>(*evt, rec_stages_range(application));
341 
342  if (!event_selector(fit, *evt)) {
343  continue;
344  }
345 
346 
347  JDataL0_t dataL0;
348 
349  for (const Hit& hit : evt->hits) {
350  dataL0.push_back(JHitL0(JDAQPMTIdentifier(hit.dom_id, hit.channel_id),
351  JAxis3D(getPosition(hit.pos), getDirection(hit.dir)),
352  JHit(hit.t, hit.tot)));
353  }
354 
355  summary.update(JDAQChronometer(evt->det_id,
356  evt->run_id,
357  evt->frame_index,
358  JDAQUTCExtended(evt->t.GetSec(), evt->t.GetNanoSec() / 16)));
359 
360  const time_converter converter = time_converter(*evt);
361 
362  Trk muon; // Monte Carlo true muon
363 
364  if (has_muon(*evt)) {
365 
366  for (const auto& trk : evt->mc_trks) {
367  if (is_muon(trk)) {
368  if (trk.E > muon.E) {
369 
370  muon = trk;
371  muon.t += converter.putTime();
372 
373  setW(muon, JSTART_LENGTH_METRES, fabs(muon.len));
374  }
375  }
376  }
377  }
378 
379 
380  bool monte_carlo = false; // show Monte Carlo true muon
381 
382  for (bool next = false; !next; ) {
383 
384  Trk trk;
385 
386  if (!monte_carlo)
387  trk = fit;
388  else
389  trk = muon;
390 
391  JRotation3D R (getDirection(trk));
392  JLine1Z tz(getPosition (trk).rotate(R), trk.t);
393  JRange<double> Z_m;
394  /*
395  if (hasW(trk, JSTART_LENGTH_METRES)) {
396  Z_m = JRange<double>(0.0, getW(fit,JSTART_LENGTH_METRES)) + JRange<double>(parameters.ZMin_m, parameters.ZMax_m);
397  }
398  */
399  const JModel<JLine1Z> match(tz, parameters.roadWidth_m, T_ns, Z_m);
400 
401  // hit selection based on fit result
402 
403  JDataW0_t data;
404 
405  for (JDataL0_t::const_iterator i = dataL0.begin(); i != dataL0.end(); ++i) {
406 
407  JHitW0 hit(*i, summary.getRate(i->getPMTIdentifier()));
408 
409  hit.rotate(R);
410 
411  if (match(hit)) {
412  data.push_back(hit);
413  }
414  }
415 
416  // select first hit in PMT
417 
418  sort(data.begin(), data.end(), JHitW0::compare);
419 
420  JDataW0_t::iterator __end = unique(data.begin(), data.end(), equal_to<JDAQPMTIdentifier>());
421 
422  double E_GeV = parameters.E_GeV;
423  /*
424  if (trk.E > 0.1) {
425  E_GeV = trk.E;
426  }
427  */
428 
429  // move fit to geometrical center of hits
430 
431  JRange<double> zs(make_array(data.begin(), __end, &JHitW0::getZ));
432 
433  const double z0 = tz.getZ();
434  const double z1 = 0.5 * (zs.getLowerLimit() + zs.getUpperLimit());
435 
436  tz.setZ(z1, getSpeedOfLight());
437 
438 
439  // graphics
440 
441  ostringstream os;
442  vector<TArrow> arrow [NUMBER_OF_PADS];
443  vector<TMarker> marker[NUMBER_OF_PADS];
444 
445  if (hasW(trk, JSTART_LENGTH_METRES) && getW(trk, JSTART_LENGTH_METRES) > 0.0) {
446 
447  marker[2].push_back(TMarker(z0 - tz.getZ(), 0.0, kFullCircle));
448  marker[2].push_back(TMarker(z0 - tz.getZ() + getW(trk, JSTART_LENGTH_METRES), 0.0, kFullCircle));
449 
450  static_cast<TAttMarker&>(marker[2][0]) = TAttMarker(kRed, kFullCircle, 0.7);
451  static_cast<TAttMarker&>(marker[2][1]) = TAttMarker(kRed, kFullCircle, 0.7);
452  }
453 
454  double chi2 = 0;
455 
456  for (JDataW0_t::const_iterator hit = data.begin(); hit != __end; ++hit) {
457 
458  const double x = hit->getX() - tz.getX();
459  const double y = hit->getY() - tz.getY();
460  const double z = hit->getZ() - tz.getZ();
461  const double R = sqrt(x*x + y*y);
462 
463  const double t1 = tz.getT() + (z + R * getTanThetaC()) * getInverseSpeedOfLight();
464 
465  JDirection3D dir(hit->getDX(), hit->getDY(), hit->getDZ()); // PMT orientation
466 
467  dir.rotate(JRotation3Z(-atan2(y,x))); // rotate PMT axis to x-z plane
468 
469  const double theta = dir.getTheta();
470  const double phi = fabs(dir.getPhi()); // rotational symmetry of Cherenkov cone
471 
472  //const double E = gWater.getE(E_GeV, z); // correct for energy loss
473  const double E = E_GeV;
474  const double dt = T_ns.constrain(hit->getT() - t1);
475 
476  JMuonPDF_t::result_type H1 = pdf.calculate(E, R, theta, phi, dt);
477  JMuonPDF_t::result_type H0(hit->getR() * 1e-9, 0.0, T_ns);
478 
479  if (H1.V >= parameters.VMax_npe) {
480  H1 *= parameters.VMax_npe / H1.V;
481  }
482 
483  H1 += H0; // signal + background
484 
485  chi2 += H1.getChi2() - H0.getChi2();
486 
487  const double derivative = H1.getDerivativeOfChi2() - H0.getDerivativeOfChi2();
488 
489  double size = derivative * arrowScale; // size of arrow
490 
491  if (fabs(size) < Amin) {
492  size = (size > 0.0 ? +Amin : -Amin);
493  } else if (fabs(size) > Amax) {
494  size = (size > 0.0 ? +Amax : -Amax);
495  }
496 
497  const double X[NUMBER_OF_PADS] = { R, atan2(y,x), z - R/getTanThetaC() };
498 
499  const double xs = (double) (NUMBER_OF_PMTS - 2 * hit->getPMTAddress()) / (double) NUMBER_OF_PMTS;
500 
501  for (size_t i = 0; i != NUMBER_OF_PADS; ++i) {
502  arrow[i].push_back(TArrow(X[i] + xs*Xs[i], dt, X[i] + xs*Xs[i], dt + size, arrowSize, arrowType.c_str()));
503  }
504  }
505 
506  os << FILL(6,'0') << evt->run_id << ":" << evt->frame_index << "/" << evt->trigger_counter << FILL();
507  os << " Q = " << FIXED(4,0) << -chi2;
508  os << " E = " << SCIENTIFIC(7,1) << trk.E << " [GeV]";
509  os << " cos(#theta) = " << FIXED(6,3) << trk.dir.z;
510 
511  if (monte_carlo)
512  os << " Monte Carlo";
513  else if (is_muon(muon))
514  os << " #Delta#alpha = " << FIXED(6,2) << getAngle(getDirection(muon), getDirection(trk)) << " [deg]";
515 
516 
517  // draw
518 
519  TLatex title(0.05, 0.5, os.str().c_str());
520 
521  title.SetTextAlign(12);
522  title.SetTextFont(42);
523  title.SetTextSize(0.6);
524 
525  p2->cd();
526 
527  title.Draw();
528 
529  for (int i = 0; i != NUMBER_OF_PADS; ++i) {
530 
531  p1->cd(i+1);
532 
533  for (auto& a1 : arrow[i]) {
534  a1.Draw();
535  }
536 
537  for (auto& m1 : marker[i]) {
538  m1.Draw();
539  }
540  }
541 
542  cv->Update();
543 
544 
545  // action
546 
547  if (batch) {
548 
549  cv->SaveAs(replace(outputFile, WILDCARD, MAKE_STRING(inputFile.getCounter())).c_str());
550 
551  next = true;
552 
553  } else {
554 
555  static int count = 0;
556 
557  if (count++ == 0) {
558  cout << endl << "Type '?' for possible options." << endl;
559  }
560 
561  for (bool user = true; user; ) {
562 
563  cout << "\n> " << flush;
564 
565  switch (JKeypress(true).get()) {
566 
567  case '?':
568  cout << endl;
569  cout << "possible options: " << endl;
570  cout << 'q' << " -> " << "exit application" << endl;
571  cout << 'u' << " -> " << "update canvas" << endl;
572  cout << 's' << " -> " << "save graphics to file" << endl;
573  cout << 'M' << " -> " << "Monte Carlo true muon information" << endl;
574  cout << 'F' << " -> " << "fit information" << endl;
575  if (event_selector.is_valid()) {
576  cout << 'L' << " -> " << "reload event selector" << endl;
577  }
578  cout << 'r' << " -> " << "rewind input file" << endl;
579  cout << 'R' << " -> " << "switch to ROOT mode (quit ROOT to continue)" << endl;
580  cout << 'p' << " -> " << "print event information" << endl;
581  cout << ' ' << " -> " << "next event (as well as any other key)" << endl;
582  break;
583 
584  case 'q':
585  cout << endl;
586  return 0;
587 
588  case 'u':
589  cv->Update();
590  break;
591 
592  case 's':
593  cv->SaveAs(replace(outputFile, WILDCARD, MAKE_STRING(inputFile.getCounter())).c_str());
594  break;
595 
596  case 'M':
597  if (is_muon(muon))
598  monte_carlo = true;
599  else
600  ERROR(endl << "No Monte Carlo muon available." << endl);
601  user = false;
602  break;
603 
604  case 'F':
605  monte_carlo = false;
606  user = false;
607  break;
608 
609  case 'L':
610  if (event_selector.is_valid()) {
611  execute(MAKE_STRING("make -f " << getPath(argv[0]) << "/JMakeEventSelector libs"), 3);
612  event_selector.reload();
613  }
614  break;
615 
616  case 'R':
617  tp->Run(kTRUE);
618  break;
619 
620  case 'p':
621  cout << endl;
622  evt->print(cout);
623  if (debug >= debug_t) {
624  for (const auto& trk : evt->mc_trks) {
625  cout << "MC "; trk.print(cout); cout << endl;
626  }
627  for (const auto& trk : evt->trks) {
628  cout << "fit "; trk.print(cout); cout << endl;
629  }
630  for (const auto& hit : evt->hits) {
631  cout << "hit "; hit.print(cout); cout << endl;
632  }
633  }
634  break;
635 
636  case 'r':
637  inputFile.rewind();
638 
639  default:
640  next = true;
641  user = false;
642  break;
643  }
644  }
645  }
646  }
647  }
648  }
649  cout << endl;
650 }
static const int JMUONSTART
Utility class to parse command line options.
Definition: JParser.hh:1514
double getAngle(const JQuaternion3D &first, const JQuaternion3D &second)
Get space angle between quanternions.
debug
Definition: JMessage.hh:29
T getLowerLimit() const
Get lower limit.
Definition: JRange.hh:202
int main(int argc, char *argv[])
Definition: Main.cc:15
Vec pos
hit position
Definition: Hit.hh:25
ROOT TTree parameter settings of various packages.
TPaveText * p1
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
then usage $script< input file >[option[primary[working directory]]] nWhere option can be E
Definition: JMuonPostfit.sh:40
double t
track time [ns] (when the particle is at pos )
Definition: Trk.hh:19
double z
Definition: Vec.hh:14
JEventSelector()
Default constructor.
Definition: JAAnet/JEvD.cc:86
then wget no check certificate user
std::string getPath(const std::string &file_name)
Get path, i.e. part before last JEEP::PATHNAME_SEPARATOR if any.
Definition: JeepToolkit.hh:148
void update(const JDAQHeader &header)
Update router.
bool is_muon(const Trk &track)
Test whether given track is a (anti-)muon.
void setW(Trk &trk, const int i, const double value)
Set associated value.
Definition: JAAnet/JEvD.cc:144
Rotation matrix.
Definition: JRotation3D.hh:111
Vec dir
track direction
Definition: Trk.hh:18
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
#define MAKE_CSTRING(A)
Make C-string.
Definition: JPrint.hh:136
Template specialisation of class JModel to match hit with muon trajectory along z-axis.
Definition: JFit/JModel.hh:34
double getZ(const JPosition3D &pos) const
Get point of emission of Cherenkov light along muon path.
Definition: JLine1Z.hh:134
Empty structure for specification of parser element that is initialised (i.e. does not require input)...
Definition: JParser.hh:83
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
Data structure for UTC time.
double E
Energy [GeV] (either MC truth or reconstructed)
Definition: Trk.hh:20
Range of reconstruction stages.
Acoustics hit.
#define MAKE_STRING(A)
Make string.
Definition: JPrint.hh:127
Event selector.
Definition: JAAnet/JEvD.cc:67
void print(std::ostream &out=std::cout) const
Print hit.
Definition: Hit.hh:60
static const char WILDCARD
Definition: JDAQTags.hh:50
Basic data structure for L0 hit.
Axis object.
Definition: JAxis3D.hh:38
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
Scanning of objects from a single file according a format that follows from the extension of each fil...
double getW(const Trk &track, const size_t index, const double value)
Get track information.
Auxiliary class for defining the range of iterations of objects.
Definition: JLimit.hh:41
static const int JMUONPREFIT
I/O formatting auxiliaries.
JAxis3D & rotate(const JRotation3D &R)
Rotate axis.
Definition: JAxis3D.hh:225
char get()
Get single character.
Definition: JKeypress.hh:74
JDirection3D getDirection(const Vec &dir)
Get direction.
JFunction1D_t::result_type result_type
Definition: JPDF_t.hh:145
JDirection3D & rotate(const JRotation3D &R)
Rotate.
Keyboard settings for unbuffered input.
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:1989
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
static const int JMUONGANDALF
double getTheta() const
Get theta angle.
Definition: JVersor3D.hh:128
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.
std::vector< double > fitinf
place to store additional fit info, see km3net-dataformat/definitions/fitparameters.csv
Definition: Trk.hh:32
#define ERROR(A)
Definition: JMessage.hh:66
then awk string
std::istream & getline(std::istream &in, JString &object)
Read string from input stream until end of line.
Definition: JString.hh:478
double len
length, if applicable [m]
Definition: Trk.hh:22
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.
double getY() const
Get y position.
Definition: JVector3D.hh:104
Auxiliary data structure for muon PDF.
Vec dir
hit direction; i.e. direction of the PMT
Definition: Hit.hh:26
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
Streaming of input and output from Linux command.
Definition: JProcess.hh:29
p2
Definition: module-Z:fit.sh:74
then JCookie sh JDataQuality D $DETECTOR_ID R
Definition: JDataQuality.sh:41
Auxiliary class for a hit with background rate value.
Definition: JHitW0.hh:21
Definition: Hit.hh:8
const double getSpeedOfLight()
Get speed of light.
static const int JSTART_LENGTH_METRES
distance between first and last hits in metres from JStart.cc
bool has_muon(const Evt &evt)
Test whether given event has a muon.
double getT(const JVector3D &pos) const
Get arrival time of Cherenkov light at given position.
Definition: JLine1Z.hh:114
static const int JMUONSIMPLEX
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 JAcoustics sh $DETECTOR_ID source JAcousticsToolkit sh CHECK_EXIT_CODE typeset A EMITTERS get_tripods $WORKDIR tripod txt EMITTERS get_transmitters $WORKDIR transmitter txt EMITTERS for EMITTER in
Definition: JCanberra.sh:48
Utility class to parse command line options.
double t
hit time (from tdc+calibration or MC truth)
Definition: Hit.hh:23
Data structure for L0 hit.
Definition: JHitL0.hh:27
const double getInverseSpeedOfLight()
Get inverse speed of light.
int dom_id
module identifier from the data (unique in the detector).
Definition: Hit.hh:14
void setZ(const double z, const double velocity)
Set z-position of vertex.
Definition: JLine1Z.hh:75
then if[[!-f $DETECTOR]] then JDetector sh $DETECTOR fi cat $WORKDIR trigger_parameters txt<< EOFtrigger3DMuon.enabled=1;trigger3DMuon.numberOfHits=5;trigger3DMuon.gridAngle_deg=1;ctMin=0.0;TMaxLocal_ns=15.0;EOF set_variable TRIGGEREFFICIENCY_TRIGGERED_EVENTS_ONLY INPUT_FILES=() for((i=1;$i<=$NUMBER_OF_RUNS;++i));do JSirene.sh $DETECTOR $JPP_DATA/genhen.km3net_wpd_V2_0.evt.gz $WORKDIR/sirene_ ${i}.root JTriggerEfficiency.sh $DETECTOR $DETECTOR $WORKDIR/sirene_ ${i}.root $WORKDIR/trigger_efficiency_ ${i}.root $WORKDIR/trigger_parameters.txt $JPP_DATA/PMT_parameters.txt INPUT_FILES+=($WORKDIR/trigger_efficiency_ ${i}.root) done for ANGLE_DEG in $ANGLES_DEG[*];do set_variable SIGMA_NS 3.0 set_variable OUTLIERS 3 set_variable OUTPUT_FILE $WORKDIR/matrix\[${ANGLE_DEG}\deg\].root $JPP_DIR/examples/JReconstruction-f"$INPUT_FILES[*]"-o $OUTPUT_FILE-S ${SIGMA_NS}-A ${ANGLE_DEG}-O ${OUTLIERS}-d ${DEBUG}--!fiif[[$OPTION=="plot"]];then if((0));then for H1 in h0 h1;do JPlot1D-f"$WORKDIR/matrix["${^ANGLES_DEG}" deg].root:${H1}"-y"1 2e3"-Y-L TR-T""-\^"number of events [a.u.]"-> o chi2
Definition: JMatrixNZ.sh:106
Data structure for fit of straight line paralel to z-axis.
Definition: JLine1Z.hh:27
unsigned int tot
tot value as stored in raw data (int for pyroot)
Definition: Hit.hh:17
double getX() const
Get x position.
Definition: JVector3D.hh:94
no fit printf nominal n $STRING awk v X
double getTanThetaC()
Get average tangent of Cherenkov angle of water corresponding to group velocity.
unsigned int channel_id
PMT channel id {0,1, .., 30} local to moduke.
Definition: Hit.hh:15
Object reading from a list of files.
const JLimit & getLimit() const
Get limit.
Definition: JLimit.hh:84
void print(std::ostream &out=std::cout) const
Print track.
Definition: Trk.hh:182
KM3NeT DAQ constants, bit handling, etc.
static const int NUMBER_OF_PMTS
Total number of PMTs in module.
Definition: JDAQ.hh:26
Function adaptor.
The Trk class represents a Monte Carlo (MC) particle as well as a reconstructed track/shower.
Definition: Trk.hh:14
Auxiliary data structure for floating point format specification.
Definition: JManip.hh:486
Wrapper class around ROOT TStyle.
Definition: JStyle.hh:20
static bool select(const Trk &trk, const Evt &evt)
Default event selection.
Definition: JAAnet/JEvD.cc:77
static const int JMUONENERGY
int debug
debug level
The Evt class respresent a Monte Carlo (MC) event as well as an offline event.
Definition: Evt.hh:20
Data structure for size of TCanvas.
Definition: JCanvas.hh:26
#define DEBUG(A)
Message macros.
Definition: JMessage.hh:62
bool hasW(const Trk &trk, const int i)
Check availability of value.
Definition: JAAnet/JEvD.cc:101