Jpp  16.0.1
the software that should make you happy
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Functions
JEvD.cc File Reference

Program to display hit probabilities. More...

#include <string>
#include <iostream>
#include <iomanip>
#include <vector>
#include <algorithm>
#include "TROOT.h"
#include "TApplication.h"
#include "TCanvas.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 "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 "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/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 JEvD.cc.

Function Documentation

int main ( int  argc,
char **  argv 
)

Definition at line 188 of file JEvD.cc.

189 {
190  using namespace std;
191  using namespace JPP;
192  using namespace KM3NETDAQ;
193 
194  typedef JParallelFileScanner< JTypeList<JDAQEvent, JEvt> > JParallelFileScanner_t;
195  typedef JParallelFileScanner_t::multi_pointer_type multi_pointer_type;
196 
197  JParallelFileScanner_t inputFile;
198  JLimit_t& numberOfEvents = inputFile.getLimit();
199  string detectorFile;
200  string pdfFile;
201  string outputFile;
203  int application;
204  JEventSelector event_selector;
205  JCanvas canvas;
206  bool batch;
207  double arrowSize;
208  string arrowType;
209  double arrowScale;
210  int debug;
211 
212 
213  try {
214 
215  parameters.numberOfPrefits = 1;
216 
217  JParser<> zap("Program to display hit probabilities.");
218 
219  zap['w'] = make_field(canvas, "size of canvas <nx>x<ny> [pixels]") = JCanvas(1200, 600);
220  zap['f'] = make_field(inputFile, "input file (output of JXXXMuonReconstruction.sh)");
221  zap['a'] = make_field(detectorFile);
222  zap['n'] = make_field(numberOfEvents) = JLimit::max();
223  zap['P'] = make_field(pdfFile);
224  zap['o'] = make_field(outputFile, "graphics output file name") = MAKE_STRING("display_" << WILDCARD << ".gif");
226  zap['A'] = make_field(application) = JMUONGANDALF, JMUONENERGY, JMUONSTART;
227  zap['L'] = make_field(event_selector) = JPARSER::initialised();
228  zap['S'] = make_field(arrowSize) = 0.003;
229  zap['T'] = make_field(arrowType) = "|->";
230  zap['F'] = make_field(arrowScale) = 250.0;
231  zap['B'] = make_field(batch, "batch processing");
232  zap['d'] = make_field(debug) = 1;
233 
234  zap(argc, argv);
235  }
236  catch(const exception& error) {
237  FATAL(error.what() << endl);
238  }
239 
240  if (batch && outputFile == "") {
241  FATAL("Missing output file name " << outputFile << " in batch mode." << endl);
242  }
243 
244  if (!batch && outputFile == "") {
245  outputFile = MAKE_STRING(WILDCARD << ".gif");
246  }
247 
248  if (outputFile.find(WILDCARD) == string::npos) {
249  FATAL("Output file name " << outputFile << " has no wild card '" << WILDCARD << "'" << endl);
250  }
251 
252 
254 
255  try {
256  load(detectorFile, detector);
257  }
258  catch(const JException& error) {
259  FATAL(error);
260  }
261 
262  const JModuleRouter router(detector);
263 
264  JSummaryFileRouter summary(inputFile, parameters.R_Hz);
265 
266  const JMuonPDF_t pdf(pdfFile, parameters.TTS_ns);
267 
268  const JTimeRange T_ns(parameters.TMin_ns, parameters.TMax_ns);
269 
270  typedef vector<JHitL0> JDataL0_t;
271  typedef vector<JHitW0> JDataW0_t;
272 
273  const JBuildL0<JHitL0> buildL0;
274 
275 
276  // Monte Carlo
277 
278  JVector3D center(0.0, 0.0, 0.0);
279 
280  try {
281  center = get<JPosition3D>(getHeader(inputFile));
282  } catch(const exception& error) {}
283 
284 
285  // ROOT
286 
287  gROOT->SetBatch(batch);
288 
289  TApplication* tp = new TApplication("user", NULL, NULL);
290  TCanvas* cv = new TCanvas("display", "", canvas.x, canvas.y);
291 
292  JSinglePointer<TStyle> gStyle(new JStyle("gplot", cv->GetWw(), cv->GetWh()));
293 
294  gROOT->SetStyle("gplot");
295  gROOT->ForceStyle();
296 
297  const size_t NUMBER_OF_PADS = 3;
298 
299  cv->SetFillStyle(4000);
300  cv->SetFillColor(kWhite);
301  cv->Divide(NUMBER_OF_PADS, 1);
302 
303  const double Dmax = getMaximalDistance(detector);
304  const double Rmin = 0.0;
305  const double Rmax = min(parameters.roadWidth_m, 0.4 * Dmax);
306  const double Tmin = min(parameters.TMin_ns, -10.0);
307  const double Tmax = max(parameters.TMax_ns, +100.0);
308  const double Amin = 0.002 * (Tmax - Tmin); // minimal arrow length [ns]
309  const double Amax = 0.8 * (Tmax - Tmin); // maximal arrow length [ns]
310  const double ymin = min(-Amax, Tmin - 0.3 * Amax);
311  const double ymax = max(+Amax, Tmax + 0.5 * Amax);
312 
313  const string Xlabel[NUMBER_OF_PADS] = { "R [m]", "#phi [rad]", "z [m]" };
314  const double Xmin [NUMBER_OF_PADS] = { Rmin, -PI, -0.5 * Dmax };
315  const double Xmax [NUMBER_OF_PADS] = { Rmax, +PI, +0.5 * Dmax };
316 
317  double Xs[NUMBER_OF_PADS];
318 
319  for (size_t i = 0; i != NUMBER_OF_PADS; ++i) {
320  Xs[i] = 0.003 * (Xmax[i] - Xmin[i]) * (0.5 * NUMBER_OF_PMTS); // x-offset arrow as function of PMT number
321  }
322 
323  TH2D H2[NUMBER_OF_PADS];
324  TGraph G2[NUMBER_OF_PADS];
325 
326  for (size_t i = 0; i != NUMBER_OF_PADS; ++i) {
327 
328  H2[i] = TH2D(MAKE_CSTRING("h" << i), NULL, 100, Xmin[i] - Xs[i], Xmax[i] + Xs[i], 100, ymin, ymax);
329 
330  H2[i].GetXaxis()->SetTitle(Xlabel[i].c_str());
331  H2[i].GetYaxis()->SetTitle("#Deltat [ns]");
332 
333  H2[i].GetXaxis()->CenterTitle(true);
334  H2[i].GetYaxis()->CenterTitle(true);
335 
336  H2[i].SetStats(kFALSE);
337 
338  G2[i].Set(2);
339 
340  G2[i].SetPoint(0, H2[i].GetXaxis()->GetXmin(), 0.0);
341  G2[i].SetPoint(1, H2[i].GetXaxis()->GetXmax(), 0.0);
342 
343  cv->cd(i+1);
344 
345  H2[i].Draw("AXIS");
346  G2[i].Draw("SAME");
347  }
348 
349 
350  for (JTreeScanner<Evt> mc(inputFile); inputFile.hasNext(); ) {
351 
352  cout << "\revent: " << setw(8) << inputFile.getCounter() << flush;
353 
354  multi_pointer_type ps = inputFile.next();
355 
356  JDAQEvent* tev = ps;
357  JEvt* in = ps;
358  Evt* event = NULL;
359 
360  if (mc.getEntries() != 0) {
361  event = mc.getEntry(tev->getCounter()); // Monte Carlo true information
362  }
363 
364 
365  in->select(JHistory::is_event(application));
366 
367  if (!in->empty()) {
368 
369  sort(in->begin(), in->end(), qualitySorter);
370 
371  if (!event_selector(*in, event)) {
372  continue;
373  }
374 
375 
376  JDataL0_t dataL0;
377 
378  buildL0(*tev, router, true, back_inserter(dataL0));
379 
380  summary.update(*tev);
381 
382 
383  JFit muon; // Monte Carlo true muon
384 
385  if (event != NULL) {
386 
387  const time_converter converter = time_converter(*event, *tev);
388 
389  for (const auto& t1 : event->mc_trks) {
390  if (is_muon(t1)) {
391  if (t1.E > muon.getE()) {
392 
393  JTrack3E ta = getTrack(t1);
394 
395  ta.add(center);
396  ta.add(converter.putTime());
397 
398  muon = getFit(0, ta, 0.0, 0, t1.E, 1);
399 
400  muon.setW(JSTART_LENGTH_METRES, fabs(t1.len));
401  }
402  }
403  }
404  }
405 
406 
407  bool next_event = false; // goto next event
408  bool monte_carlo = false; // show Monte Carlo true muon
409  size_t index = 0; // index of fit
410 
411  while (!next_event) {
412 
413  JFit fit;
414 
415  if (!monte_carlo)
416  fit = (*in)[index];
417  else
418  fit = muon;
419 
420 
421  JRotation3D R (getDirection(fit));
422  JLine1Z tz(getPosition (fit).rotate(R), fit.getT());
423  JRange<double> Z_m;
424  /*
425  if (fit.hasW(JSTART_LENGTH_METRES)) {
426  Z_m = JRange<double>(0.0, fit.getW(JSTART_LENGTH_METRES)) + JRange<double>(parameters.ZMin_m, parameters.ZMax_m);
427  }
428  */
429  const JModel<JLine1Z> match(tz, parameters.roadWidth_m, T_ns, Z_m);
430 
431  // hit selection based on fit result
432 
433  JDataW0_t data;
434 
435  for (JDataL0_t::const_iterator i = dataL0.begin(); i != dataL0.end(); ++i) {
436 
437  JHitW0 hit(*i, summary.getRate(i->getPMTIdentifier()));
438 
439  hit.rotate(R);
440 
441  if (match(hit)) {
442  data.push_back(hit);
443  }
444  }
445 
446  // select first hit in PMT
447 
448  sort(data.begin(), data.end(), JHitW0::compare);
449 
450  JDataW0_t::iterator __end = unique(data.begin(), data.end(), equal_to<JDAQPMTIdentifier>());
451 
452  double E_GeV = parameters.E_GeV;
453  /*
454  if (fit.getE() > 0.1) {
455  E_GeV = fit.getE();
456  }
457  */
458 
459  // move fit to geometrical center of hits
460 
461  JRange<double> zs(make_array(data.begin(), __end, &JHitW0::getZ));
462 
463  const double z0 = tz.getZ();
464  const double z1 = 0.5 * (zs.getLowerLimit() + zs.getUpperLimit());
465 
466  tz.setZ(z1, getSpeedOfLight());
467 
468 
469  // graphics
470 
471  TLatex latex [NUMBER_OF_PADS];
472  vector<TArrow> arrow [NUMBER_OF_PADS];
473  vector<TMarker> marker[NUMBER_OF_PADS];
474 
475  for (int i = 0; i != NUMBER_OF_PADS; ++i) {
476 
477  latex[i].SetTextAlign(12);
478  latex[i].SetTextFont(42);
479  latex[i].SetTextSize(0.06);
480 
481  latex[i].SetX(H2[i].GetXaxis()->GetXmin() + 0.05 * (H2[i].GetXaxis()->GetXmax() - H2[i].GetXaxis()->GetXmin()));
482  latex[i].SetY(H2[i].GetYaxis()->GetXmax() - 0.05 * (H2[i].GetYaxis()->GetXmax() - H2[i].GetYaxis()->GetXmin()));
483  }
484 
485  if (fit.hasW(JSTART_LENGTH_METRES) && fit.getW(JSTART_LENGTH_METRES) > 0.0) {
486 
487  marker[2].push_back(TMarker(z0 - tz.getZ(), 0.0, kFullCircle));
488  marker[2].push_back(TMarker(z0 - tz.getZ() + fit.getW(JSTART_LENGTH_METRES), 0.0, kFullCircle));
489 
490  static_cast<TAttMarker&>(marker[2][0]) = TAttMarker(kRed, kFullCircle, 0.7);
491  static_cast<TAttMarker&>(marker[2][1]) = TAttMarker(kRed, kFullCircle, 0.7);
492  }
493 
494  double chi2 = 0;
495 
496  for (JDataW0_t::const_iterator hit = data.begin(); hit != __end; ++hit) {
497 
498  const double x = hit->getX() - tz.getX();
499  const double y = hit->getY() - tz.getY();
500  const double z = hit->getZ() - tz.getZ();
501  const double R = sqrt(x*x + y*y);
502 
503  const double t1 = tz.getT() + (z + R * getTanThetaC()) * getInverseSpeedOfLight();
504 
505  JDirection3D dir(hit->getDX(), hit->getDY(), hit->getDZ()); // PMT orientation
506 
507  dir.rotate(JRotation3Z(-atan2(y,x))); // rotate PMT axis to x-z plane
508 
509  const double theta = dir.getTheta();
510  const double phi = fabs(dir.getPhi()); // rotational symmetry of Cherenkov cone
511 
512  //const double E = gWater.getE(E_GeV, z); // correct for energy loss
513  const double E = E_GeV;
514  const double dt = T_ns.constrain(hit->getT() - t1);
515 
516  JMuonPDF_t::result_type H1 = pdf.calculate(E, R, theta, phi, dt);
517  JMuonPDF_t::result_type H0(hit->getR() * 1e-9, 0.0, T_ns);
518 
519  if (H1.V >= parameters.VMax_npe) {
520  H1 *= parameters.VMax_npe / H1.V;
521  }
522 
523  H1 += H0; // signal + background
524 
525  chi2 += H1.getChi2() - H0.getChi2();
526 
527  const double derivative = H1.getDerivativeOfChi2() - H0.getDerivativeOfChi2();
528 
529  double size = derivative * arrowScale; // size of arrow
530 
531  if (fabs(size) < Amin) {
532  size = (size > 0.0 ? +Amin : -Amin);
533  } else if (fabs(size) > Amax) {
534  size = (size > 0.0 ? +Amax : -Amax);
535  }
536 
537  const double X[NUMBER_OF_PADS] = { R, atan2(y,x), z - R/getTanThetaC() };
538 
539  const double xs = (double) (NUMBER_OF_PMTS - 2 * hit->getPMTAddress()) / (double) NUMBER_OF_PMTS;
540 
541  for (size_t i = 0; i != NUMBER_OF_PADS; ++i) {
542  arrow[i].push_back(TArrow(X[i] + xs*Xs[i], dt, X[i] + xs*Xs[i], dt + size, arrowSize, arrowType.c_str()));
543  }
544  }
545 
546  latex[0].SetTitle(MAKE_CSTRING("" << FILL(6,'0') << tev->getRunNumber() << FILL() <<
547  ":" << tev->getFrameIndex() <<
548  " " << tev->getCounter()));
549 
550  if (!monte_carlo)
551  latex[1].SetTitle(MAKE_CSTRING("[" << index << "]" << " "
552  "Q = " << FIXED(4,0) << getQuality(chi2) << " "
553  "#beta = " << FIXED(6,2) << fit.getW(JGANDALF_BETA0_RAD, 0.0) * 180.0/PI << " [deg]"));
554  else
555  latex[1].SetTitle(MAKE_CSTRING("Q = " << FIXED(4,0) << getQuality(chi2)));
556 
557  if (muon.getStatus() >= 0) {
558 
559  if (!monte_carlo)
560  latex[2].SetTitle(MAKE_CSTRING("#Delta#alpha = " << FIXED(6,2) << getAngle(getDirection(muon), getDirection(fit)) << " [deg]"));
561  else
562  latex[2].SetTitle("Monte Carlo");
563  }
564 
565 
566  // draw
567 
568  for (int i = 0; i != NUMBER_OF_PADS; ++i) {
569 
570  cv->cd(i+1);
571 
572  latex[i].Draw();
573 
574  for (auto& a1 : arrow[i]) {
575  a1.Draw();
576  }
577 
578  for (auto& m1 : marker[i]) {
579  m1.Draw();
580  }
581  }
582 
583  cv->Update();
584 
585 
586  // action
587 
588  if (batch) {
589 
590  cv->SaveAs(replace(outputFile, WILDCARD, MAKE_STRING(inputFile.getCounter())).c_str());
591 
592  next_event = true;
593 
594  } else {
595 
596  for (int c = 0; c == 0; ) {
597 
598  if (!batch) {
599 
600  static bool first = true;
601 
602  if (first) {
603  cout << endl << "Type '?' for possible options." << endl;
604  }
605 
606  first = false;
607  }
608 
609  cout << "\n> " << flush;
610 
611  switch((c = JKeypress(true).get())) {
612 
613  case '?':
614  cout << endl
615  << "possible options: " << endl
616  << 'q' << " -> " << "exit application" << endl
617  << 'u' << " -> " << "update canvas" << endl
618  << 's' << " -> " << "save graphics to file" << endl
619  << '+' << " -> " << "next fit" << endl
620  << '-' << " -> " << "previous fit" << endl
621  << 'M' << " -> " << "Monte Carlo true muon information" << endl
622  << 'F' << " -> " << "fit information" << endl
623  << 'x' << " -> " << "remake standard event selector" << endl
624  << 'L' << " -> " << "reload event selector" << endl
625  << 'r' << " -> " << "rewind input file" << endl
626  << 'R' << " -> " << "switch to ROOT mode (quit ROOT to continue)" << endl
627  << ' ' << " -> " << "next event (as well as any other key)" << endl;
628  c = 0;
629  break;
630 
631  case 'q':
632  cout << endl;
633  return 0;
634 
635  case 'u':
636  cv->Update();
637  c = 0;
638  break;
639 
640  case 's':
641  cv->SaveAs(replace(outputFile, WILDCARD, MAKE_STRING(inputFile.getCounter())).c_str());
642  c = 0;
643  break;
644 
645  case '+':
646  monte_carlo = false;
647  index = (index != in->size() - 1 ? index + 1 : 0);
648  break;
649 
650  case '-':
651  monte_carlo = false;
652  index = (index != 0 ? index - 1 : in->size() - 1);
653  break;
654 
655  case 'M':
656  if (muon.getStatus() >= 0)
657  monte_carlo = true;
658  else
659  ERROR(endl << "No Monte Carlo muon available." << endl);
660  break;
661 
662  case 'F':
663  monte_carlo = false;
664  break;
665 
666  case 'r':
667  inputFile.rewind();
668  next_event = true;
669  break;
670 
671  case 'x':
672  execute(MAKE_STRING("make -C " << getPath(argv[0]) << "./" << " libs"), debug);
673  c = 0;
674  break;
675 
676  case 'L':
677  if (event_selector.is_valid())
678  event_selector.reload();
679  else
680  ERROR(endl << "No event selector (use command line option -L)." << endl);
681  break;
682 
683  case 'R':
684  tp->Run(kTRUE);
685  c = 0;
686  break;
687 
688  default:
689  next_event = true;
690  break;
691  }
692  }
693  }
694  }
695  }
696  }
697  cout << endl;
698 }
static const int JMUONSTART
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
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 E
Definition: JMuonPostfit.sh:35
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
std::string getPath(const std::string &file_name)
Get path, i.e. part before last JEEP::PATHNAME_SEPARATOR if any.
Definition: JeepToolkit.hh:108
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
*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
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.
Rotation around Z-axis.
Definition: JRotation3D.hh:85
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.
Detector file.
Definition: JHead.hh:224
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.
#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
JPosition3D getPosition(const Vec &pos)
Get position.
#define ERROR(A)
Definition: JMessage.hh:66
then break fi done getCenter read X Y Z let X
static const int JGANDALF_BETA0_RAD
KM3NeT Data Definitions v2.1.0-12-g9520e6e https://git.km3net.de/common/km3net-dataformat.
double putTime() const
Get Monte Carlo time minus DAQ/trigger time.
static const double PI
Mathematical constants.
File router for fast addressing of summary data.
int debug
debug level
Definition: JSirene.cc:63
then usage $script[distance] fi case set_variable R
Definition: JDrawLED.sh:43
Auxiliary data structure for sequence of same character.
Definition: JManip.hh:328
Auxiliary data structure for muon PDF.
Definition: JPDF_t.hh:135
#define FATAL(A)
Definition: JMessage.hh:67
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.
static const int JSTART_LENGTH_METRES
distance between first and last hits in metres from JStart.cc
void load(const std::string &file_name, JDetector &detector)
Load detector from input file.
const double getInverseSpeedOfLight()
Get inverse speed of light.
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
$WORKDIR ev_configure_domsimulator txt echo process $DOM_SIMULATOR $i $SOURCE_HOST[$index] csh c(setenv ROOTSYS $ROOTSYS &&source $JPP_DIR/setenv.csh $JPP_DIR &&($DOM_SIMULATOR\-u\$NAME\$\-H\$SERVER\$\-M\$LOGGER\$\-d $DEBUG</dev/null > &/dev/null &))'
do set_variable DETECTOR_TXT $WORKDIR detector
static const int NUMBER_OF_PMTS
Total number of PMTs in module.
Definition: JDAQ.hh:26
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:42
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