Jpp
JSirene.cc
Go to the documentation of this file.
1 #include <string>
2 #include <iostream>
3 #include <iomanip>
4 #include <vector>
5 #include <limits>
6 
7 #include "TROOT.h"
8 #include "TFile.h"
9 #include "TH1D.h"
10 #include "TProfile.h"
11 #include "TRandom3.h"
12 
13 #include "evt/Head.hh"
14 #include "evt/Evt.hh"
15 
16 #include "JLang/JSharedPointer.hh"
17 #include "JLang/Jpp.hh"
18 #include "JLang/JPredicate.hh"
19 #include "JSystem/JDate.hh"
20 
21 #include "JPhysics/JPDFToolkit.hh"
22 #include "JPhysics/JCDFTable.hh"
23 #include "JPhysics/JPDFTypes.hh"
25 #include "JTools/JFunction1D_t.hh"
27 
28 #include "JAAnet/JHit_t.hh"
29 #include "JAAnet/JTrk_t.hh"
30 #include "JAAnet/JAAnetToolkit.hh"
31 #include "JAAnet/JPDB.hh"
32 
33 #include "JSirene/km3-opa_efrac.hh"
34 //#include "JSirene/pythia.hh"
35 //#include "JSirene/JPythia.hh"
36 #include "JSirene/JSeaWater.hh"
37 #include "JSirene/JCDFTable1D.hh"
38 #include "JSirene/JCDFTable2D.hh"
40 
41 #include "JDetector/JDetector.hh"
45 
49 #include "JSupport/JSupport.hh"
50 #include "JSupport/JMeta.hh"
51 
52 #include "Jeep/JProperties.hh"
53 #include "Jeep/JPrint.hh"
54 #include "Jeep/JTimer.hh"
55 #include "Jeep/JParser.hh"
56 #include "Jeep/JMessage.hh"
57 
58 
59 int debug; //!< debug level
60 
61 namespace {
62 
63  using namespace JPP;
64 
65  typedef JHermiteSplineFunction1D_t JFunction1D_t;
68  JPolint1FunctionalGridMap>::maplist J3DMap_t;
72  JPolint1FunctionalGridMap>::maplist J4DMap_t;
73 
74  typedef JCDFTable<JFunction1D_t, J3DMap_t> JCDF4D_t; // muon
75  typedef JCDFTable<JFunction1D_t, J4DMap_t> JCDF5D_t; // shower
76 
79 
80 
81  const int NUMBER_OF_BINS = 200; //!< number of bins for CDF integral
82  const double SAFETY_FACTOR = 1.7; //!< safety factor for CDF integral
83 
84 
85  /**
86  * Auxiliary class for CDF tables.
87  */
88  template<class function_type, // CDF function
89  class integral_type> // CDF integral
90  struct JABC {
91  /**
92  * Constructor.
93  *
94  * \param file_descriptor file name descriptor
95  * \param type PDF type
96  */
97  JABC(const std::string& file_descriptor, const JPDFType_t type)
98  {
99  this->type = type;
100 
101  const std::string file_name = getFilename(file_descriptor, this->type);
102 
103  STATUS("loading input from file " << file_name << "... " << flush);
104 
105  try {
106  function.load(file_name.c_str());
107  }
108  catch(const JException& error) {
109  FATAL(error.what() << endl);
110  }
111 
112  integral = integral_type(function, NUMBER_OF_BINS, SAFETY_FACTOR);
113 
114  STATUS("OK" << endl);
115  }
116 
117  JPDFType_t type;
118  function_type function;
119  integral_type integral;
120  };
121 
122 
123  typedef JABC<JCDF4D_t, JCDF1D_t> JCDF_t; //!< muon light
124  typedef JABC<JCDF5D_t, JCDF2D_t> JCDG_t; //!< shower light
125 }
126 
127 
128 /**
129  * \file
130  *
131  * Main program to simulate detector response to muons and showers.
132  *
133  * \image html sirene.png "Picture by Claudine Colnard"
134  *
135  * Note that CDF file descriptor should contain the wild card character JPHYSICS::WILD_CARD;\n
136  * The file names are obtained by replacing JPHYSICS::WILD_CARD; with
137  * - JPHYSICS::DIRECT_LIGHT_FROM_MUON;
138  * - JPHYSICS::SCATTERED_LIGHT_FROM_MUON;
139  * - JPHYSICS::DIRECT_LIGHT_FROM_DELTARAYS;
140  * - JPHYSICS::SCATTERED_LIGHT_FROM_DELTARAYS;
141  * - JPHYSICS::DIRECT_LIGHT_FROM_EMSHOWER; and
142  * - JPHYSICS::SCATTERED_LIGHT_FROM_EMSHOWER,
143  * respectively.
144  *
145  * The CDF tables can be produced with the script <tt>JMakePDF.sh</tt>:
146  * <pre>
147  * JMakePDF.sh -P
148  * </pre>
149  * (option <tt>-h</tt> will print all available options).
150  * Note that the script will launch a large number of processes (<tt>JMakePDF</tt> and <tt>JMakePDG</tt>)\n
151  * which may take a considerable amount of time to completion.\n
152  * On a standard desktop, all jobs should be finished within 1/2 a day or so.
153  *
154  * The same script should then be run with option <tt>-M</tt> to merge the PDF files, i.e:
155  * <pre>
156  * JMakePDF.sh -M
157  * </pre>
158  *
159  * CDF tables are obtained by running the same script with option <tt>-C</tt>, i.e:
160  * <pre>
161  * JMakePDF.sh -C
162  * </pre>
163  *
164  * The various PDFs can be drawn using the <tt>JDrawPDF</tt> or <tt>JDrawPDG</tt> applications.\n
165  * The tabulated PDFs can be plotted using the <tt>JPlotPDF</tt> or <tt>JPlotPDG</tt> applications.\n
166  * The tabulated CDFs can be plotted using the <tt>JPlotCDF</tt> or <tt>JPlotCDG</tt> applications.
167  * \author mdejong
168  */
169 int main(int argc, char **argv)
170 {
171  using namespace std;
172  using namespace JPP;
173 
174  string fileDescriptor;
177  JLimit_t& numberOfEvents = inputFile.getLimit();
178  string detectorFile;
179  double Ecut_GeV = 0.1; // minimal energy for generation of light from shower [GeV]
180  double Emin_GeV = 1.0; // minimal energy of muon for shower generation [GeV]
181  double Dmin_m = 0.1; // minimal distance for positioning [m]
182  double Tmax_ns = 0.1; // maximal time between hits on same PMT to be merged
183  int Nmax_npe = numeric_limits<int>::max(); // maximal number of photo-electrons
184  bool writeEMShowers;
185  bool keep;
186  size_t numberOfHits;
187  double factor;
188  UInt_t seed;
189 
190  try {
191 
192  JProperties properties;
193 
194  properties.insert(gmake_property(Ecut_GeV));
195  properties.insert(gmake_property(Emin_GeV));
196  properties.insert(gmake_property(Dmin_m));
197  properties.insert(gmake_property(Tmax_ns));
198  properties.insert(gmake_property(Nmax_npe));
199 
200  JParser<> zap("Main program to simulate detector response to muons and showers.");
201 
202  zap['@'] = make_field(properties) = JPARSER::initialised();
203  zap['F'] = make_field(fileDescriptor, "file name descriptor for CDF tables");
204  zap['f'] = make_field(inputFile) = JPARSER::initialised();
205  zap['o'] = make_field(outputFile) = "sirene.root";
206  zap['n'] = make_field(numberOfEvents) = JLimit::max();
207  zap['a'] = make_field(detectorFile) = "";
208  zap['s'] = make_field(writeEMShowers, "store generated EM showers in event");
209  zap['k'] = make_field(keep, "keep position of tracks");
210  zap['N'] = make_field(numberOfHits, "minimum number of hits to output event") = 1;
211  zap['U'] = make_field(factor, "scaling factor applied to light yields") = 1.0;
212  zap['S'] = make_field(seed) = 0;
213  zap['d'] = make_field(debug) = 1;
214 
215  zap(argc, argv);
216  }
217  catch(const exception &error) {
218  FATAL(error.what() << endl);
219  }
220 
221 
222  gRandom->SetSeed(seed);
223 
224 
225  const JMeta meta(argc, argv);
226 
227  const double Zbed = 0.0; // level of seabed [m]
228 
229  vector<JCDF_t> CDF;
230  vector<JCDG_t> CDG;
231 
232  CDF.push_back(JCDF_t(fileDescriptor, DIRECT_LIGHT_FROM_MUON));
233  CDF.push_back(JCDF_t(fileDescriptor, SCATTERED_LIGHT_FROM_MUON));
234  CDF.push_back(JCDF_t(fileDescriptor, DIRECT_LIGHT_FROM_DELTARAYS));
235  CDF.push_back(JCDF_t(fileDescriptor, SCATTERED_LIGHT_FROM_DELTARAYS));
236 
237  CDG.push_back(JCDG_t(fileDescriptor, DIRECT_LIGHT_FROM_EMSHOWER));
238  CDG.push_back(JCDG_t(fileDescriptor, SCATTERED_LIGHT_FROM_EMSHOWER));
239 
240 
241  double maximal_road_width = 0.0; // road width [m]
242  double maximal_distance = 0.0; // road width [m]
243 
244  for (size_t i = 0; i != CDF.size(); ++i) {
245 
246  DEBUG("Range CDF["<< CDF[i].type << "] " << CDF[i].function.intensity.getXmax() << " m" << endl);
247 
248  maximal_road_width = max(maximal_road_width, CDF[i].function.intensity.getXmax());
249  }
250 
251  for (size_t i = 0; i != CDG.size(); ++i) {
252 
253  DEBUG("Range CDG["<< CDG[i].type << "] " << CDG[i].function.intensity.getXmax() << " m" << endl);
254 
255  if (!is_scattered(CDF[i].type)) {
256  maximal_road_width = max(maximal_road_width, CDG[i].function.intensity.getXmax());
257  }
258 
259  maximal_distance = max(maximal_distance, CDG[i].function.intensity.getXmax());
260  }
261 
262  NOTICE("Maximal road width [m] " << maximal_road_width << endl);
263  NOTICE("Maximal distance [m] " << maximal_distance << endl);
264 
265 
266  if (detectorFile == "" || inputFile.empty()) {
267  STATUS("Nothing to be done." << endl);
268  return 0;
269  }
270 
272 
273  try {
274 
275  STATUS("Load detector... " << flush);
276 
277  load(detectorFile, detector);
278 
279  STATUS("OK" << endl);
280  }
281  catch(const JException& error) {
282  FATAL(error);
283  }
284 
285 
287 
288  if (true) {
289 
290  STATUS("Setting up radiation tables... " << flush);
291 
292  const JRadiation hydrogen( 1.0, 1.0, 40, 0.01, 0.1, 0.1);
293  const JRadiation oxygen ( 8.0, 16.0, 40, 0.01, 0.1, 0.1);
294  const JRadiation chlorine(17.0, 35.0, 40, 0.01, 0.1, 0.1);
295 
296  JSharedPointer<JRadiation> Hydrogen(new JRadiationFunction(hydrogen, 300, 0.2, 1.0e11));
297  JSharedPointer<JRadiation> Oxygen (new JRadiationFunction(oxygen, 300, 0.2, 1.0e11));
298  JSharedPointer<JRadiation> Chlorine(new JRadiationFunction(chlorine, 300, 0.2, 1.0e11));
299 
300  JRadiationSource::source_type EErad(&JRadiation::TotalCrossSectionEErad, &JRadiation::EfromEErad);
301  JRadiationSource::source_type Brems(&JRadiation::TotalCrossSectionBrems, &JRadiation::EfromBrems);
302  JRadiationSource::source_type GNrad(&JRadiation::TotalCrossSectionGNrad, &JRadiation::EfromGNrad);
303 
304  radiation.push_back(new JRadiationSource(Oxygen, DENSITY_SEA_WATER * JSeaWater::O(), EErad));
305  radiation.push_back(new JRadiationSource(Chlorine, DENSITY_SEA_WATER * JSeaWater::Cl(), EErad));
306  radiation.push_back(new JRadiationSource(Hydrogen, DENSITY_SEA_WATER * JSeaWater::H(), EErad));
307 
308  radiation.push_back(new JRadiationSource(Oxygen, DENSITY_SEA_WATER * JSeaWater::O(), Brems));
309  radiation.push_back(new JRadiationSource(Chlorine, DENSITY_SEA_WATER * JSeaWater::Cl(), Brems));
310  radiation.push_back(new JRadiationSource(Hydrogen, DENSITY_SEA_WATER * JSeaWater::H(), Brems));
311 
312  radiation.push_back(new JRadiationSource(Oxygen, DENSITY_SEA_WATER * JSeaWater::O(), GNrad));
313  radiation.push_back(new JRadiationSource(Chlorine, DENSITY_SEA_WATER * JSeaWater::Cl(), GNrad));
314  radiation.push_back(new JRadiationSource(Hydrogen, DENSITY_SEA_WATER * JSeaWater::H(), GNrad));
315 
316  STATUS("OK" << endl);
317  }
318 
319 
320  Vec center(0,0,0);
321  Head header;
322 
323  try {
324 
325  header = inputFile.getHeader();
326 
327  JHead buffer(header);
328 
329  buffer.simul.push_back(JAANET::simul());
330 
331  buffer.simul.rbegin()->program = "JSirene";
332  buffer.simul.rbegin()->version = getGITVersion();
333  buffer.simul.rbegin()->date = getDate();
334  buffer.simul.rbegin()->time = getTime();
335 
336  buffer.push(&JHead::simul);
337 
338  center = get<Vec>(buffer);
339 
340  if (!keep) {
341 
342  buffer.coord_origin = coord_origin(0.0, 0.0, 0.0);
343  buffer.can.zmin += center.z;
344  buffer.can.zmax += center.z;
345 
346  buffer.push(&JHead::coord_origin);
347  buffer.push(&JHead::can);
348  }
349 
350  const JCircle2D circle(detector.begin(), detector.end());
351 
352  center += Vec(circle.getX(), circle.getY(), 0.0);
353 
354  copy(buffer, header);
355  }
356  catch(const JException& error) {
357  FATAL(error);
358  }
359 
360  if (!keep)
361  NOTICE("Offset applied to true tracks is: " << center << endl);
362  else
363  NOTICE("No offset applied to true tracks." << endl);
364 
365  JCylinder3D cylinder(detector.begin(), detector.end());
366 
367  cylinder.addMargin(maximal_distance);
368 
369  if (cylinder.getZmin() < Zbed) {
370  cylinder.setZmin(Zbed);
371  }
372 
373  NOTICE("Ligth generation volume: " << cylinder << endl);
374 
375  TProfile cpu("cpu", NULL, 14, 1.0, 8.0);
376  TH1D job("job", NULL, 400, 0.5, 400.5);
377 
378 
379  outputFile.open();
380 
381  if (!outputFile.is_open()) {
382  FATAL("Error opening file " << outputFile << endl);
383  }
384 
385  outputFile.put(meta);
386  outputFile.put(header);
387  outputFile.put(*gRandom);
388 
389 
390  JTimer timer;
391 
392  for (JMultipleFileScanner<Evt>& in = inputFile; in.hasNext(); ) {
393 
394  STATUS("event: " << setw(10) << in.getCounter() << '\r'); DEBUG(endl);
395 
396  job.Fill(1.0);
397 
398  Evt* evt = in.next();
399 
400  if (!keep) {
401  for (vector<Trk>::iterator track = evt->mc_trks.begin(); track != evt->mc_trks.end(); ++track) {
402  track->pos += center;
403  }
404  }
405 
406  Evt event(*evt); // output
407 
408  event.mc_hits.clear();
409 
410  JHits_t mc_hits; // temporary buffer
411 
412  timer.reset();
413  timer.start();
414 
415  for (vector<Trk>::const_iterator track = evt->mc_trks.begin(); track != evt->mc_trks.end(); ++track) {
416 
417  if (is_muon(*track)) {
418 
419  // -----------------------------------------------
420  // muon
421  // -----------------------------------------------
422 
423  job.Fill(2.0);
424 
426 
427  const pair<double, double> intersection = cylinder.getIntersection(getAxis(*track));
428 
429  double Zmin = intersection.first;
430  double Zmax = intersection.second;
431 
432  if (Zmax - Zmin <= Dmin_m) {
433  continue;
434  }
435 
436  JVertex vertex(0.0, track->t, track->E); // start of muon
437 
438  if (track->pos.z < Zbed) { // propagate muon through rock
439 
440  if (track->dir.z > 0.0)
441  vertex.step(gRock, (Zbed - track->pos.z) / track->dir.z);
442  else
443  continue;
444  }
445 
446  if (vertex.getZ() < Zmin) { // propagate muon through water
447  vertex.step(gWater, Zmin - vertex.getZ());
448  }
449 
450  if (vertex.getRange() <= Dmin_m) {
451  continue;
452  }
453 
454  job.Fill(3.0);
455 
456  const JDetectorSubset_t subdetector(detector, getAxis(*track), maximal_road_width);
457 
458  if (subdetector.empty()) {
459  continue;
460  }
461 
462  job.Fill(4.0);
463 
464  JTrack muon(vertex); // propagate muon trough detector
465 
466  while (vertex.getE() >= Emin_GeV && vertex.getZ() < Zmax) {
467 
468  const int N = radiation.size();
469 
470  double li[N]; // inverse interaction lengths
471  double ls = 1.0e-5; // minimal total inverse interaction length (100 km)^-1
472 
473  for (int i = 0; i != N; ++i) {
474  ls += li[i] = radiation[i]->getInverseInteractionLength(vertex.getE());
475  }
476 
477  const double ds = min(gRandom->Exp(1.0) / ls, vertex.getRange());
478 
479  vertex.step(ds);
480 
481  if (vertex.getE() >= Emin_GeV) {
482 
483  double Es = 0.0; // shower energy [GeV]
484  double y = gRandom->Uniform(ls);
485 
486  for (int i = 0; i != N; ++i) {
487 
488  y -= li[i];
489 
490  if (y < 0.0) {
491  Es = radiation[i]->getEnergyOfShower(vertex.getE()); // shower energy [GeV]
492  break;
493  }
494  }
495 
496  vertex.applyEloss(Es);
497 
498  if (Es >= Ecut_GeV) {
499  muon.push_back(vertex);
500  }
501  }
502  }
503 
504  // add muon end point
505 
506  if (vertex.getE() < Emin_GeV && vertex.getRange() > 0.0) {
507 
508  vertex.step(vertex.getRange());
509 
510  muon.push_back(vertex);
511  }
512 
513  // add information to output muon
514 
515  vector<Trk>::iterator trk = find_if(event.mc_trks.begin(),
516  event.mc_trks.end(),
517  make_predicate(&Trk::id, track->id));
518 
519  if (trk != event.mc_trks.end()) {
520  trk->len = (muon.rbegin()->getZ() < Zmax ? +1 : -1) * (muon.rbegin()->getZ() - muon.begin()->getZ());
521  trk->setusr(mc_usr_keys::energy_lost_in_can, muon.begin()->getE() - muon.rbegin()->getE());
522  }
523 
524  for (JDetector::const_iterator module = subdetector.begin(); module != subdetector.end(); ++module) {
525 
526  const double z0 = muon.begin()->getZ();
527  const double t0 = muon.begin()->getT();
528  const double R = module->getX();
529  const double Z = module->getZ() - R / getTanThetaC();
530 
531  if (Z >= muon.begin()->getZ() && Z <= muon.rbegin()->getZ()) {
532 
533  for (size_t i = 0; i != CDF.size(); ++i) {
534 
535  if (R < CDF[i].integral.getXmax()) {
536 
537  try {
538 
539  double W = 1.0; // mip
540 
541  if (is_deltarays(CDF[i].type)) {
542  W = getDeltaRaysFromMuon(muon.getE(Z)); // delta-rays
543  }
544 
545  const double NPE = CDF[i].integral.getNPE(R) * module->size() * factor * W;
546  const int N = gRandom->Poisson(NPE);
547 
548  if (N != 0) {
549 
550  double ns = 0.0;
551 
552  for (JModule::const_iterator pmt = module->begin(); pmt != module->end(); ++pmt) {
553 
554  const double R = pmt->getX();
555  const double Z = pmt->getZ() - z0;
556  const double theta = pmt->getTheta();
557  const double phi = fabs(pmt->getPhi());
558 
559  const double npe = CDF[i].function.getNPE(R, theta, phi) * factor * W;
560 
561  ns += npe;
562 
563  int n0 = min(getNumberOfPhotoElectrons(NPE, N, npe), Nmax_npe);
564 
565  job.Fill((double) (100 + CDF[i].type), (double) n0);
566 
567  while (n0 != 0) {
568 
569  const double t1 = CDF[i].function.getTime(R, theta, phi, gRandom->Rndm());
570  const int n1 = getNumberOfPhotoElectrons(n0);
571 
572  mc_hits.push_back(JHit_t(mc_hits.size() + 1,
573  pmt->getID(),
574  getHitType(CDF[i].type),
575  track->id,
576  t0 + (R * getTanThetaC() + Z) / C + t1,
577  n1));
578 
579  n0 -= n1;
580  }
581  }
582 
583  if (ns > NPE) {
584  job.Fill((double) (300 + CDF[i].type));
585  }
586  }
587  }
588  catch(const exception& error) {
589  job.Fill((double) (200 + CDF[i].type));
590  }
591  }
592  }
593  }
594  }
595 
596  for (JTrack::const_iterator vertex = muon.begin(); vertex != muon.end(); ++vertex) {
597 
598  const double z0 = vertex->getZ();
599  const double t0 = vertex->getT();
600  const double Es = vertex->getEs();
601 
602  int origin = track->id;
603 
604  if (writeEMShowers) {
605  origin = event.mc_trks.size() + 1;
606  }
607 
608  int number_of_hits = 0;
609 
610  JDetectorSubset_t::range_type range = subdetector.getRange(z0 - maximal_distance,
611  z0 + maximal_distance);
612 
613  for (JDetector::const_iterator module = range.begin(); module != range.end(); ++module) {
614 
615  const double dz = geanz.getMaximum(Es);
616  const double R = module->getX();
617  const double Z = module->getZ() - z0 - dz;
618  const double D = sqrt(R*R + Z*Z);
619  const double cd = Z / D;
620 
621  for (size_t i = 0; i != CDG.size(); ++i) {
622 
623  if (D < CDG[i].integral.getXmax()) {
624 
625  try {
626 
627  const double NPE = CDG[i].integral.getNPE(D, cd) * Es * module->size() * factor;
628  const int N = gRandom->Poisson(NPE);
629 
630  if (N != 0) {
631 
632  double ns = 0.0;
633 
634  for (JModule::const_iterator pmt = module->begin(); pmt != module->end(); ++pmt) {
635 
636  const double dz = geanz.getMaximum(Es);
637  const double R = pmt->getX();
638  const double Z = pmt->getZ() - z0 - dz;
639  const double theta = pmt->getTheta();
640  const double phi = fabs(pmt->getPhi());
641  const double D = sqrt(R*R + Z*Z);
642  const double cd = Z / D;
643 
644  const double npe = CDG[i].function.getNPE(D, cd, theta, phi) * Es * factor;
645 
646  ns += npe;
647 
648  int n0 = min(getNumberOfPhotoElectrons(NPE, N, npe), Nmax_npe);
649 
650  job.Fill((double) (100 + CDG[i].type), (double) n0);
651 
652  while (n0 != 0) {
653 
654  const double dz = geanz.getLength(Es, gRandom->Rndm());
655  const double Z = pmt->getZ() - z0 - dz;
656  const double D = sqrt(R*R + Z*Z);
657  const double cd = Z / D;
658 
659  const double t1 = CDG[i].function.getTime(D, cd, theta, phi, gRandom->Rndm());
660  const int n1 = getNumberOfPhotoElectrons(n0);
661 
662  mc_hits.push_back(JHit_t(mc_hits.size() + 1,
663  pmt->getID(),
664  getHitType(CDG[i].type),
665  origin,
666  t0 + (dz + D * getIndexOfRefraction()) / C + t1,
667  n1));
668 
669  n0 -= n1;
670 
671  number_of_hits += n1;
672  }
673  }
674 
675  if (ns > NPE) {
676  job.Fill((double) (300 + CDG[i].type));
677  }
678  }
679  }
680  catch(const exception& error) {
681  job.Fill((double) (200 + CDG[i].type));
682  }
683  }
684  }
685  }
686 
687  if (writeEMShowers && number_of_hits != 0) {
688 
689  event.mc_trks.push_back(JTrk_t(origin,
691  track->id,
692  track->pos + track->dir * vertex->getZ(),
693  track->dir,
694  vertex->getT(),
695  Es));
696  }
697  }
698 
699  } else if (is_tau(*track)) {
700 
701  // -----------------------------------------------
702  // tau
703  // -----------------------------------------------
704 
705  job.Fill(6.0);
706 
707  const double z0 = 0.0;
708  const double z1 = z0 + track->len;
709  const double t0 = track->t;
710  const double E = track->E;
711 
712  const JTransformation3D transformation = getTransformation(*track);
713 
714  JModule buffer;
715 
716  for (JDetector::const_iterator module = detector.begin(); module != detector.end(); ++module) {
717 
718  const JPosition3D pos = transformation.transform(module->getPosition());
719 
720  const double R = pos.getX();
721  const double Z = pos.getZ() - R / getTanThetaC();
722 
723  if (Z < z0 ||
724  Z > z1 ||
725  R > maximal_road_width) {
726  continue;
727  }
728 
729  for (size_t i = 0; i != CDF.size(); ++i) {
730 
731  if (R < CDF[i].integral.getXmax()) {
732 
733  try {
734 
735  double W = 1.0; // mip
736 
737  if (is_deltarays(CDF[i].type)) {
738  W = getDeltaRaysFromTau(E); // delta-rays
739  }
740 
741  const double NPE = CDF[i].integral.getNPE(R) * module->size() * factor * W;
742  const int N = gRandom->Poisson(NPE);
743 
744  if (N != 0) {
745 
746  buffer = *module;
747 
748  buffer.transform(transformation);
749 
750  double ns = 0.0;
751 
752  for (JModule::const_iterator pmt = buffer.begin(); pmt != buffer.end(); ++pmt) {
753 
754  const double R = pmt->getX();
755  const double Z = pmt->getZ() - z0;
756  const double theta = pmt->getTheta();
757  const double phi = fabs(pmt->getPhi());
758 
759  const double npe = CDF[i].function.getNPE(R, theta, phi) * factor * W;
760 
761  ns += npe;
762 
763  int n0 = min(getNumberOfPhotoElectrons(NPE, N, npe), Nmax_npe);
764 
765  job.Fill((double) (120 + CDF[i].type), (double) n0);
766 
767  while (n0 != 0) {
768 
769  const double t1 = CDF[i].function.getTime(R, theta, phi, gRandom->Rndm());
770  const int n1 = getNumberOfPhotoElectrons(n0);
771 
772  mc_hits.push_back(JHit_t(mc_hits.size() + 1,
773  pmt->getID(),
774  getHitType(CDF[i].type),
775  track->id,
776  t0 + (R * getTanThetaC() + Z) / C + t1,
777  n1));
778 
779  n0 -= n1;
780  }
781  }
782 
783  if (ns > NPE) {
784  job.Fill((double) (320 + CDF[i].type));
785  }
786  }
787  }
788  catch(const exception& error) {
789  job.Fill((double) (220 + CDF[i].type));
790  }
791 
792  }
793  }
794  }
795 
796  if (!buffer.empty()) {
797  job.Fill(7.0);
798  }
799 
800  } else if (!is_neutrino(*track)) {
801 
802  if (JPDB::getInstance().hasPDG(track->type)) {
803 
804  // -----------------------------------------------
805  // electron or hadron
806  // -----------------------------------------------
807 
808  job.Fill(8.0);
809 
810  double E = track->E;
811 
812  try {
813  E = getKineticEnergy(E, JPDB::getInstance().getPDG(track->type).mass);
814  }
815  catch(const exception& error) {
816  ERROR(error.what() << endl);
817  }
818 
819  E = pythia(track->type, E);
820 
821  if (E < Ecut_GeV || cylinder.getDistance(getPosition(*track)) > Dmin_m) {
822  continue;
823  }
824 
825  const double z0 = 0.0;
826  const double t0 = track->t;
827 
828  const JTransformation3D transformation = getTransformation(*track);
829 
830  JModule buffer;
831 
832  for (JDetector::const_iterator module = detector.begin(); module != detector.end(); ++module) {
833 
834  const JPosition3D pos = transformation.transform(module->getPosition());
835 
836  const double dz = geanz.getMaximum(E);
837  const double R = pos.getX();
838  const double Z = pos.getZ() - z0 - dz;
839  const double D = sqrt(R*R + Z*Z);
840  const double cd = Z / D;
841 
842  if (D > maximal_distance) {
843  continue;
844  }
845 
846  for (size_t i = 0; i != CDG.size(); ++i) {
847 
848  if (D < CDG[i].integral.getXmax()) {
849 
850  try {
851 
852  const double NPE = CDG[i].integral.getNPE(D, cd) * E * module->size() * factor;
853  const int N = gRandom->Poisson(NPE);
854 
855  if (N != 0) {
856 
857  buffer = *module;
858 
859  buffer.transform(transformation);
860 
861  double ns = 0.0;
862 
863  for (JModule::const_iterator pmt = buffer.begin(); pmt != buffer.end(); ++pmt) {
864 
865  const double dz = geanz.getMaximum(E);
866  const double R = pmt->getX();
867  const double Z = pmt->getZ() - z0 - dz;
868  const double theta = pmt->getTheta();
869  const double phi = fabs(pmt->getPhi());
870  const double D = sqrt(R*R + Z*Z);
871  const double cd = Z / D;
872 
873  const double npe = CDG[i].function.getNPE(D, cd, theta, phi) * E * factor;
874 
875  ns += npe;
876 
877  int n0 = min(getNumberOfPhotoElectrons(NPE, N, npe), Nmax_npe);
878 
879  job.Fill((double) (140 + CDG[i].type), (double) n0);
880 
881  while (n0 != 0) {
882 
883  const double dz = geanz.getLength(E, gRandom->Rndm());
884  const double Z = pmt->getZ() - z0 - dz;
885  const double D = sqrt(R*R + Z*Z);
886  const double cd = Z / D;
887 
888  const double t1 = CDG[i].function.getTime(D, cd, theta, phi, gRandom->Rndm());
889  const int n1 = getNumberOfPhotoElectrons(n0);
890 
891  mc_hits.push_back(JHit_t(mc_hits.size() + 1,
892  pmt->getID(),
893  getHitType(CDG[i].type, true),
894  track->id,
895  t0 + (dz + D * getIndexOfRefraction()) / C + t1,
896  n1));
897 
898  n0 -= n1;
899  }
900  }
901 
902  if (ns > NPE) {
903  job.Fill((double) (340 + CDG[i].type));
904  }
905  }
906  }
907  catch(const exception& error) {
908  job.Fill((double) (240 + CDG[i].type));
909  }
910  }
911  }
912  }
913 
914  if (!buffer.empty()) {
915  job.Fill(9.0);
916  }
917 
918  } else {
919 
920  job.Fill(21.0);
921  }
922  }
923  }
924 
925  if (!mc_hits.empty()) {
926 
927  mc_hits.merge(Tmax_ns);
928 
929  event.mc_hits.resize(mc_hits.size());
930 
931  copy(mc_hits.begin(), mc_hits.end(), event.mc_hits.begin());
932  }
933 
934  timer.stop();
935 
936  if (has_neutrino(event)) {
937  cpu.Fill(log10(get_neutrino(event).E), (double) timer.usec_ucpu * 1.0e-3);
938  }
939 
940  if (event.mc_hits.size() >= numberOfHits) {
941 
942  outputFile.put(event);
943 
944  job.Fill(10.0);
945  }
946  }
947  STATUS(endl);
948 
949  outputFile.put(cpu);
950  outputFile.put(job);
951  outputFile.put(*gRandom);
952  outputFile.close();
953 }
JAANET::JHead::simul
std::vector< JAANET::simul > simul
Definition: JHead.hh:1077
JPHYSICS::DIRECT_LIGHT_FROM_DELTARAYS
direct light from delta-rays
Definition: JPDFTypes.hh:35
JMeta.hh
JTOOLS::JPolint1FunctionalGridMap
Type definition of a 1st degree polynomial interpolation based on a JGridMap implementation.
Definition: JFunctionalMap_t.hh:91
JCDFTable1D.hh
JPredicate.hh
JSIRENE::JVertex::step
JVertex & step(const double ds)
Step.
Definition: JSireneToolkit.hh:266
JAANET::JHit_t
Auxiliary class to set-up Hit.
Definition: JHit_t.hh:25
JDETECTOR::JModule::transform
void transform(const JRotation3D &R, const JVector3D &pos)
Transformation of geometry (see method JGEOMETRY3D::JPosition3D::transform(const JRotation3D&,...
Definition: JModule.hh:261
JEEP::JTimer::reset
void reset()
Reset timer.
Definition: JTimer.hh:76
JPHYSICS::JRadiationFunction
Fast implementation of class JRadiation.
Definition: JRadiationSource.hh:35
JFileRecorder.hh
JEEP::JTimer
Auxiliary class for CPU timing and usage.
Definition: JTimer.hh:32
JSUPPORT::JLimit
Auxiliary class for defining the range of iterations of objects.
Definition: JLimit.hh:41
JGEOMETRY3D::JTransformation3D
Transformation.
Definition: JTransformation3D.hh:24
JPHYSICS::JCDFTable
Multi-dimensional CDF table for arrival time of Cherenkov light.
Definition: JCDFTable.hh:51
JSIRENE::JVertex
Vertex of energy loss of muon.
Definition: JSireneToolkit.hh:146
JMessage.hh
JSIRENE::JVertex::getEs
double getEs() const
Get shower energy.
Definition: JSireneToolkit.hh:214
JPHYSICS::gWater
static const JGeaneWater gWater
Function object for energy loss of muon in sea water.
Definition: JGeane.hh:328
JPrint.hh
JAANET::getAxis
JAxis3D getAxis(const Trk &track)
Get axis.
Definition: JAAnetToolkit.hh:244
JLANG::getInstance
T & getInstance(const T &object)
Get static instance from temporary object.
Definition: JObject.hh:73
JTOOLS::DENSITY_SEA_WATER
static const double DENSITY_SEA_WATER
Fixed environment values.
Definition: JConstants.hh:34
JAANET::JHead::coord_origin
JAANET::coord_origin coord_origin
Definition: JHead.hh:1086
JSIRENE::getHitType
JHitType_t getHitType(const JPDFType_t pdf, const bool shower=false)
Get hit type corresponding to given PDF type.
Definition: JSireneToolkit.hh:112
JPHYSICS::is_deltarays
bool is_deltarays(const int pdf)
Test if given PDF type corresponds to Cherenkov light from delta-rays.
Definition: JPDFTypes.hh:172
JGEOMETRY2D::JCircle2D
Data structure for circle in two dimensions.
Definition: JCircle2D.hh:29
JGEOMETRY3D::JVector3D::getZ
double getZ() const
Get z position.
Definition: JVector3D.hh:114
JEEP::JTimer::start
void start()
Start timer.
Definition: JTimer.hh:89
JPARSER::initialised
Empty structure for specification of parser element that is initialised (i.e.
Definition: JParser.hh:63
JDETECTOR::load
void load(const JString &file_name, JDetector &detector)
Load detector from input file.
Definition: JDetectorToolkit.hh:456
JCDFTable.hh
JAANET::is_tau
bool is_tau(const Trk &track)
Test whether given track is a (anti-)tau.
Definition: JAAnetToolkit.hh:375
JSharedPointer.hh
JPHYSICS::DIRECT_LIGHT_FROM_EMSHOWER
direct light from EM shower
Definition: JPDFTypes.hh:40
JGEOMETRY3D::JTransformation3D::transform
JPosition3D transform(const JPosition3D &pos) const
Transform position.
Definition: JTransformation3D.hh:49
JPDFToolkit.hh
std::vector
Definition: JSTDTypes.hh:12
debug
int debug
debug level
Definition: JSirene.cc:59
JSIRENE::pythia
static const JPythia pythia
Function object for relative light yield as a function of GEANT particle code.
Definition: JPythia.hh:96
main
int main(int argc, char **argv)
Definition: JSirene.cc:169
JTOOLS::H
static const double H
Planck constant [eV s].
Definition: JConstants.hh:25
JPARSER::JParser
Utility class to parse command line options.
Definition: JParser.hh:1493
JTOOLS::JMAPLIST
Auxiliary class for recursive map list generation.
Definition: JMapList.hh:108
JPHYSICS::getDeltaRaysFromMuon
double getDeltaRaysFromMuon(const double E)
Equivalent EM-shower energy due to delta-rays per unit muon track length.
Definition: JPDFToolkit.hh:81
JAANET::getTransformation
JTransformation3D getTransformation(const Trk &track)
Get transformation.
Definition: JAAnetToolkit.hh:280
JDate.hh
JTOOLS::JHermiteSplineFunction1D_t
Type definition of a spline interpolation method based on a JCollection with double result type.
Definition: JFunction1D_t.hh:92
NOTICE
#define NOTICE(A)
Definition: JMessage.hh:64
JAANET::coord_origin
Coordinate origin.
Definition: JHead.hh:489
JPHYSICS::JCDFTable2D
Custom class for CDF table in 2 dimensions.
Definition: JCDFTable2D.hh:39
JTOOLS::C
static const double C
Speed of light in vacuum [m/ns].
Definition: JConstants.hh:22
JHit_t.hh
JGEOMETRY3D::JCylinder3D::addMargin
void addMargin(const double D)
Add (safety) margin.
Definition: JCylinder3D.hh:163
JAANET::JHits_t::merge
void merge(const double Tmax_ns)
Merge hits on same PMT that are within given time window.
Definition: JHit_t.hh:118
JAANET::copy
void copy(const Head &from, JHead &to)
Copy header from from to to.
Definition: JHead.cc:152
JAANET::is_muon
bool is_muon(const Trk &track)
Test whether given track is a (anti-)muon.
Definition: JAAnetToolkit.hh:367
JAANET::JHead
Monte Carlo run header.
Definition: JHead.hh:839
JCDFTable2D.hh
JAANET::can::zmin
double zmin
Bottom [m].
Definition: JHead.hh:369
JSupport.hh
JPP
This name space includes all other name spaces (except KM3NETDAQ, KM3NET and ANTARES).
Definition: JAAnetToolkit.hh:37
ERROR
#define ERROR(A)
Definition: JMessage.hh:66
JFunction1D_t.hh
JPHYSICS::JGeanz::getMaximum
double getMaximum(const double E) const
Get depth of shower maximum.
Definition: JGeanz.hh:162
JGEOMETRY3D::JCylinder3D
Cylinder object.
Definition: JCylinder3D.hh:37
JPHYSICS::gRock
static const JGeane_t gRock(2.67e-1 *0.9 *DENSITY_ROCK, 3.40e-4 *1.2 *DENSITY_ROCK)
Function object for energy loss of muon in rock.
JSIRENE::getKineticEnergy
double getKineticEnergy(const double E, const double m)
Get kinetic energy of particle with given mass.
Definition: JSireneToolkit.hh:348
JDETECTOR::JDetectorSubset
Detector subset with binary search functionality.
Definition: JDetectorSubset.hh:118
JTOOLS::JPolint1FunctionalMap
Type definition of a 1st degree polynomial interpolation based on a JMap implementation.
Definition: JFunctionalMap_t.hh:55
JGEOMETRY3D::JPosition3D
Data structure for position in three dimensions.
Definition: JPosition3D.hh:35
JTOOLS::getIndexOfRefraction
double getIndexOfRefraction()
Get average index of refraction of water.
Definition: JConstants.hh:111
JFunctionalMap_t.hh
JAANET::JTrk_t
Auxiliary class to set-up Trk.
Definition: JTrk_t.hh:19
JAANET::get_neutrino
const Trk & get_neutrino(const Evt &evt)
Get incoming neutrino.
Definition: JAAnetToolkit.hh:438
JAAnetToolkit.hh
JPHYSICS::SCATTERED_LIGHT_FROM_EMSHOWER
scattered light from EM shower
Definition: JPDFTypes.hh:41
JPHYSICS::getDeltaRaysFromTau
double getDeltaRaysFromTau(const double E)
Equivalent EM-shower energy due to delta-rays per unit tau track length.
Definition: JPDFToolkit.hh:109
JAANET::has_neutrino
bool has_neutrino(const Evt &evt)
Test whether given event has an incoming neutrino.
Definition: JAAnetToolkit.hh:427
JSIRENE::JVertex::applyEloss
void applyEloss(const double Es)
Apply energy loss energy.
Definition: JSireneToolkit.hh:252
JAANET::getTime
double getTime(const Hit &hit)
Get true time of hit.
Definition: JAAnetToolkit.hh:87
JSYSTEM::ls
Auxiliary class to list files.
Definition: JSystemToolkit.hh:300
JDETECTOR::JModule
Data structure for a composite optical module.
Definition: JModule.hh:49
JEEP::getFilename
std::string getFilename(const std::string &file_name)
Get file name part, i.e.
Definition: JeepToolkit.hh:86
JLANG::getGITVersion
const char * getGITVersion()
Get GIT version.
Definition: Jpp.hh:58
JMultipleFileScanner.hh
JEEP::JTimer::usec_ucpu
unsigned long long usec_ucpu
Definition: JTimer.hh:225
std::pair
Definition: JSTDTypes.hh:15
STATUS
#define STATUS(A)
Definition: JMessage.hh:63
JSUPPORT::JMultipleFileScanner::hasNext
virtual bool hasNext()
Check availability of next element.
Definition: JMultipleFileScanner.hh:350
JAANET::can::zmax
double zmax
Top [m].
Definition: JHead.hh:370
Jpp.hh
km3-opa_efrac.hh
JParser.hh
JDetectorToolkit.hh
JSIRENE::JVertex::getZ
double getZ() const
Get position.
Definition: JSireneToolkit.hh:181
JPHYSICS::JPDFType_t
JPDFType_t
PDF types.
Definition: JPDFTypes.hh:27
JTOOLS::getTanThetaC
double getTanThetaC()
Get average tangent of Cherenkov angle of water.
Definition: JConstants.hh:133
JTrk_t.hh
JPHYSICS::is_scattered
bool is_scattered(const int pdf)
Test if given PDF type corresponds to scattered light.
Definition: JPDFTypes.hh:186
JAANET::JHead::can
JAANET::can can
Definition: JHead.hh:1083
JAANET::JHead::push
void push(T JHead::*pd)
Push given data member to Head.
Definition: JHead.hh:944
JSireneToolkit.hh
make_field
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:1954
gmake_property
#define gmake_property(A)
macro to convert (template) parameter to JPropertiesElement object
Definition: JProperties.hh:1319
JSUPPORT::JMultipleFileScanner
General purpose class for object reading from a list of file names.
Definition: JMultipleFileScanner.hh:167
JDETECTOR::JDetector
Detector data structure.
Definition: JDetector.hh:80
JAANET::is_neutrino
bool is_neutrino(const Trk &track)
Test whether given track is a neutrino.
Definition: JAAnetToolkit.hh:349
JEEP::JProperties
Utility class to parse parameter values.
Definition: JProperties.hh:496
JAANET::getPosition
JPosition3D getPosition(const Vec &v)
Get position.
Definition: JAAnetToolkit.hh:197
JAANET::detector
Detector file.
Definition: JHead.hh:130
JPHYSICS::JGeanz::getLength
double getLength(const double E, const double P, const double eps=1.0e-3) const
Get shower length for a given integrated probability.
Definition: JGeanz.hh:122
JSUPPORT::JMeta
Auxiliary class for ROOT I/O of application specific meta data.
Definition: JMeta.hh:71
DEBUG
#define DEBUG(A)
Message macros.
Definition: JMessage.hh:62
JSIRENE::getNumberOfPhotoElectrons
int getNumberOfPhotoElectrons(const double NPE, const int N, const double npe)
Get number of photo-electrons of a hit given the expectation values of the number of photo-electrons ...
Definition: JSireneToolkit.hh:54
std
Definition: jaanetDictionary.h:36
JCylinder3D.hh
JLANG::make_predicate
JPredicate< JResult_t T::*, JComparison::eq > make_predicate(JResult_t T::*member, const JResult_t value)
Helper method to create predicate for data member.
Definition: JPredicate.hh:128
JPDB.hh
JSYSTEM::getDate
const char * getDate()
Get ASCII formatted date.
Definition: JDate.hh:23
JSIRENE::JVertex::getE
double getE() const
Get muon energy.
Definition: JSireneToolkit.hh:203
JEEP::JTimer::stop
void stop()
Stop timer.
Definition: JTimer.hh:113
JTimer.hh
JSIRENE::JTrack::getE
double getE(const double z) const
Get muon energy at given position.
Definition: JSireneToolkit.hh:322
JProperties.hh
JMonteCarloFileSupportkit.hh
JRadiationSource.hh
JAANET::TRACK_TYPE_PHOTON
Definition: JParticleTypes.hh:70
JLANG::JException::what
virtual const char * what() const
Get error message.
Definition: JException.hh:65
JPHYSICS::DIRECT_LIGHT_FROM_MUON
direct light from muon
Definition: JPDFTypes.hh:29
JSIRENE::JVertex::getT
double getT() const
Get time.
Definition: JSireneToolkit.hh:192
JDetector.hh
JGEOMETRY3D::JVector3D::getX
double getX() const
Get x position.
Definition: JVector3D.hh:93
JDetectorSubset.hh
JSIRENE::JVertex::getRange
double getRange() const
Get range of muon.
Definition: JSireneToolkit.hh:226
FATAL
#define FATAL(A)
Definition: JMessage.hh:67
outputFile
string outputFile
Definition: JDAQTimesliceSelector.cc:37
JSIRENE::JTrack
Muon trajectory.
Definition: JSireneToolkit.hh:303
JPHYSICS::SCATTERED_LIGHT_FROM_DELTARAYS
scattered light from delta-rays
Definition: JPDFTypes.hh:36
JSUPPORT::JFileRecorder
Object writing to file.
Definition: JFileRecorder.hh:41
JPHYSICS::JRadiationSource
Source of radiation.
Definition: JRadiationSource.hh:231
JSeaWater.hh
JPDFTypes.hh
JAANET::JHits_t
Auxiliary data structure for list of hits with hit merging capability.
Definition: JHit_t.hh:107
JAANET::simul
Generator for simulation.
Definition: JHead.hh:308
JLANG::JException
General exception.
Definition: JException.hh:40
JDETECTOR::JDetectorSubset_t
Detector subset without binary search functionality.
Definition: JDetectorSubset.hh:37
JPHYSICS::JRadiation
Auxiliary class for the calculation of the muon radiative cross sections.
Definition: JRadiation.hh:39
JLANG::JSharedPointer
The template JSharedPointer class can be used to share a pointer to an object.
Definition: JSharedPointer.hh:28
JPHYSICS::JCDFTable1D
Custom class for CDF table in 1 dimension.
Definition: JCDFTable1D.hh:38
JPHYSICS::geanz
static const JGeanz geanz(1.85, 0.62, 0.54)
Function object for longitudinal EM-shower profile.
JPHYSICS::SCATTERED_LIGHT_FROM_MUON
scattered light from muon
Definition: JPDFTypes.hh:30