CMS 3D CMS Logo

NuclearInteraction.cc
Go to the documentation of this file.
1 // system headers
2 #include <cmath>
3 #include <memory>
4 #include <vector>
5 #include <map>
6 #include <string>
7 #include <fstream>
8 #include <iostream>
9 #include <sys/stat.h>
10 #include <cmath>
11 #include <TFile.h>
12 #include <TTree.h>
13 #include <TROOT.h>
14 #include <Math/AxisAngle.h>
15 #include "Math/Boost.h"
17 
18 // Framework Headers
23 
25 
26 // Fast Sim headers
34 
35 // Math
37 
39 // Author: Patrick Janot
40 // Date: 25-Jan-2007
41 //
42 // Revision: Class structure modified to match SimplifiedGeometryPropagator
43 // Unified the way, the closest charged daughter is defined
44 // S. Kurz, 29 May 2017
46 
47 // TODO: save() function called in destructer. Does that actually make sense?
48 
51 
52 namespace fastsim {
54 
59  public:
62 
64  ~NuclearInteraction() override;
65 
67 
73  void interact(fastsim::Particle& particle,
75  std::vector<std::unique_ptr<fastsim::Particle> >& secondaries,
76  const RandomEngineAndDistribution& random) override;
77 
78  private:
80  unsigned index(int thePid);
81 
83  XYZVector orthogonal(const XYZVector& aVector) const;
84 
86  void save();
87 
90 
91  double theDistCut;
92  double theHadronEnergy;
94 
96  // Read/Save nuclear interactions from FullSim
98 
99  TFile* theFile = nullptr;
100  std::vector<std::vector<TTree*> > theTrees;
101  std::vector<std::vector<TBranch*> > theBranches;
102  std::vector<std::vector<NUEvent*> > theNUEvents;
103  std::vector<std::vector<unsigned> > theCurrentEntry;
104  std::vector<std::vector<unsigned> > theCurrentInteraction;
105  std::vector<std::vector<unsigned> > theNumberOfEntries;
106  std::vector<std::vector<unsigned> > theNumberOfInteractions;
107  std::vector<std::vector<std::string> > theFileNames;
108  std::vector<std::vector<double> > theHadronCM;
109 
110  unsigned ien4;
111 
112  std::ofstream myOutputFile;
113  unsigned myOutputBuffer;
114 
116 
118  // Properties of the Hadrons
120 
122  static std::vector<std::vector<double> > theRatiosMap;
124  static std::map<int, int> theIDMap;
125 
127  const std::vector<double> theHadronEN = {
128  1.0, 2.0, 3.0, 4.0, 5.0, 7.0, 9.0, 12.0, 15.0, 20.0, 30.0, 50.0, 100.0, 200.0, 300.0, 500.0, 700.0, 1000.0};
129 
130  // Use same order for all vectors below (Properties for "piplus", "piminus", "K0L", "Kplus", "Kminus", "p", "pbar", "n", "nbar")
131 
133  const std::vector<int> theHadronID = {211, -211, 130, 321, -321, 2212, -2212, 2112, -2112};
135  const std::vector<std::string> theHadronNA = {
136  "piplus", "piminus", "K0L", "Kplus", "Kminus", "p", "pbar", "n", "nbar"};
138  const std::vector<double> theHadronMA = {
139  0.13957, 0.13957, 0.497648, 0.493677, 0.493677, 0.93827, 0.93827, 0.939565, 0.939565};
141  const std::vector<double> theHadronPMin = {0.7, 0.0, 1.0, 1.0, 0.0, 1.1, 0.0, 1.1, 0.0};
143  const std::vector<double> theLengthRatio = {//pi+ pi- K0L K+ K- p pbar n nbar
144  0.2257,
145  0.2294,
146  0.3042,
147  0.2591,
148  0.2854,
149  0.3101,
150  0.5216,
151  0.3668,
152  0.4898};
154  const std::vector<double> theRatios = {// pi+ (211)
155  0.031390573,
156  0.531842852,
157  0.819614219,
158  0.951251711,
159  0.986382750,
160  1.000000000,
161  0.985087033,
162  0.982996773,
163  0.990832192,
164  0.992237923,
165  0.994841580,
166  0.973816742,
167  0.967264815,
168  0.971714258,
169  0.969122824,
170  0.978681792,
171  0.977312732,
172  0.984255819,
173  // pi- (-211)
174  0.035326512,
175  0.577356403,
176  0.857118809,
177  0.965683504,
178  0.989659360,
179  1.000000000,
180  0.989599240,
181  0.980665408,
182  0.988384816,
183  0.981038152,
184  0.975002104,
185  0.959996152,
186  0.953310808,
187  0.954705592,
188  0.957615400,
189  0.961150456,
190  0.965022184,
191  0.960573304,
192  // K0L (130)
193  0.000000000,
194  0.370261189,
195  0.649793096,
196  0.734342408,
197  0.749079499,
198  0.753360057,
199  0.755790543,
200  0.755872164,
201  0.751337674,
202  0.746685288,
203  0.747519634,
204  0.739357554,
205  0.735004444,
206  0.803039922,
207  0.832749896,
208  0.890900187,
209  0.936734805,
210  1.000000000,
211  // K+ (321)
212  0.000000000,
213  0.175571717,
214  0.391683394,
215  0.528946472,
216  0.572818635,
217  0.614210280,
218  0.644125538,
219  0.670304050,
220  0.685144573,
221  0.702870161,
222  0.714708513,
223  0.730805263,
224  0.777711536,
225  0.831090576,
226  0.869267129,
227  0.915747562,
228  0.953370523,
229  1.000000000,
230  // K- (-321)
231  0.000000000,
232  0.365353210,
233  0.611663677,
234  0.715315908,
235  0.733498956,
236  0.738361302,
237  0.745253654,
238  0.751459671,
239  0.750628335,
240  0.746442657,
241  0.750850669,
242  0.744895986,
243  0.735093960,
244  0.791663444,
245  0.828609543,
246  0.889993040,
247  0.940897842,
248  1.000000000,
249  // proton (2212)
250  0.000000000,
251  0.042849136,
252  0.459103223,
253  0.666165343,
254  0.787930873,
255  0.890397011,
256  0.920999533,
257  0.937832788,
258  0.950920131,
259  0.966595049,
260  0.979542270,
261  0.988061653,
262  0.983260159,
263  0.988958431,
264  0.991723494,
265  0.995273237,
266  1.000000000,
267  0.999962634,
268  // anti-proton (-2212)
269  1.000000000,
270  0.849956907,
271  0.775625988,
272  0.802018230,
273  0.816207485,
274  0.785899785,
275  0.754998487,
276  0.728977244,
277  0.710010673,
278  0.670890339,
279  0.665627872,
280  0.652682888,
281  0.613334247,
282  0.647534574,
283  0.667910938,
284  0.689919693,
285  0.709200185,
286  0.724199928,
287  // neutron (2112)
288  0.000000000,
289  0.059216484,
290  0.437844536,
291  0.610370629,
292  0.702090648,
293  0.780076890,
294  0.802143073,
295  0.819570432,
296  0.825829666,
297  0.840079750,
298  0.838435509,
299  0.837529986,
300  0.835687165,
301  0.885205014,
302  0.912450156,
303  0.951451221,
304  0.973215562,
305  1.000000000,
306  // anti-neutron
307  1.000000000,
308  0.849573257,
309  0.756479495,
310  0.787147094,
311  0.804572414,
312  0.791806302,
313  0.760234588,
314  0.741109531,
315  0.724118186,
316  0.692829761,
317  0.688465897,
318  0.671806061,
319  0.636461171,
320  0.675314029,
321  0.699134460,
322  0.724305037,
323  0.742556115,
324  0.758504713};
325 
327  // Used to build the ID map 'theIDMap' (i.e., what is to be considered as a proton, etc...)
329 
330  const std::vector<int> protonsID = {2212, 3222, -101, -102, -103, -104};
331  const std::vector<int> antiprotonsID = {-2212, -3222};
332  const std::vector<int> neutronsID = {2112, 3122, 3112, 3312, 3322, 3334, -3334};
333  const std::vector<int> antineutronsID = {-2112, -3122, -3112, -3312, -3322};
334  const std::vector<int> K0LsID = {130, 310};
335  const std::vector<int> KplussesID = {321};
336  const std::vector<int> KminussesID = {-321};
337  const std::vector<int> PiplussesID = {211};
338  const std::vector<int> PiminussesID = {-211};
339  };
340 
341  // TODO: Is this correct?
342  // Thread safety
343  static std::once_flag initializeOnce;
344  CMS_THREAD_GUARD(initializeOnce) std::vector<std::vector<double> > NuclearInteraction::theRatiosMap;
346 } // namespace fastsim
347 
349  : fastsim::InteractionModel(name), currentValuesWereSet(false) {
350  // Full path to FullSim root file
352 
353  // Read from config file
354  theDistCut = cfg.getParameter<double>("distCut");
355  theHadronEnergy = cfg.getParameter<double>("hadronEnergy");
356  inputFile = cfg.getUntrackedParameter<std::string>("inputFile", "");
357 
358  // The evolution of the interaction lengths with energy
359  // initialize once for all possible instances
360  std::call_once(initializeOnce, [this]() {
361  theRatiosMap.resize(theHadronID.size());
362  for (unsigned i = 0; i < theHadronID.size(); ++i) {
363  for (unsigned j = 0; j < theHadronEN.size(); ++j) {
364  theRatiosMap[i].push_back(theRatios[i * theHadronEN.size() + j]);
365  }
366  }
367 
368  // Build the ID map (i.e., what is to be considered as a proton, etc...)
369  // Protons
370  for (const auto& id : protonsID)
371  theIDMap[id] = 2212;
372  // Anti-Protons
373  for (const auto& id : antiprotonsID)
374  theIDMap[id] = -2212;
375  // Neutrons
376  for (const auto& id : neutronsID)
377  theIDMap[id] = 2112;
378  // Anti-Neutrons
379  for (const auto& id : antineutronsID)
380  theIDMap[id] = -2112;
381  // K0L's
382  for (const auto& id : K0LsID)
383  theIDMap[id] = 130;
384  // K+'s
385  for (const auto& id : KplussesID)
386  theIDMap[id] = 321;
387  // K-'s
388  for (const auto& id : KminussesID)
389  theIDMap[id] = -321;
390  // pi+'s
391  for (const auto& id : PiplussesID)
392  theIDMap[id] = 211;
393  // pi-'s
394  for (const auto& id : PiminussesID)
395  theIDMap[id] = -211;
396  });
397 
398  // Prepare the map of files
399  // Loop over the particle names
400  TFile* aVFile = nullptr;
401  std::vector<TTree*> aVTree(theHadronEN.size());
402  std::vector<TBranch*> aVBranch(theHadronEN.size());
403  std::vector<NUEvent*> aVNUEvents(theHadronEN.size());
404  std::vector<unsigned> aVCurrentEntry(theHadronEN.size());
405  std::vector<unsigned> aVCurrentInteraction(theHadronEN.size());
406  std::vector<unsigned> aVNumberOfEntries(theHadronEN.size());
407  std::vector<unsigned> aVNumberOfInteractions(theHadronEN.size());
408  std::vector<std::string> aVFileName(theHadronEN.size());
409  std::vector<double> aVHadronCM(theHadronEN.size());
410  theTrees.resize(theHadronNA.size());
411  theBranches.resize(theHadronNA.size());
412  theNUEvents.resize(theHadronNA.size());
413  theCurrentEntry.resize(theHadronNA.size());
414  theCurrentInteraction.resize(theHadronNA.size());
415  theNumberOfEntries.resize(theHadronNA.size());
416  theNumberOfInteractions.resize(theHadronNA.size());
417  theFileNames.resize(theHadronNA.size());
418  theHadronCM.resize(theHadronNA.size());
419  theFile = aVFile;
420  for (unsigned iname = 0; iname < theHadronNA.size(); ++iname) {
421  theTrees[iname] = aVTree;
422  theBranches[iname] = aVBranch;
423  theNUEvents[iname] = aVNUEvents;
424  theCurrentEntry[iname] = aVCurrentEntry;
425  theCurrentInteraction[iname] = aVCurrentInteraction;
426  theNumberOfEntries[iname] = aVNumberOfEntries;
427  theNumberOfInteractions[iname] = aVNumberOfInteractions;
428  theFileNames[iname] = aVFileName;
429  theHadronCM[iname] = aVHadronCM;
430  }
431 
432  // Read the information from a previous run (to keep reproducibility)
433  currentValuesWereSet = this->read(inputFile);
435  std::cout << "***WARNING*** You are reading nuclear-interaction information from the file " << inputFile
436  << " created in an earlier run." << std::endl;
437 
438  // Open the file for saving the information of the current run
439  myOutputFile.open("NuclearInteractionOutputFile.txt");
440  myOutputBuffer = 0;
441 
442  // Open the root file
443  edm::FileInPath myDataFile("FastSimulation/MaterialEffects/data/NuclearInteractions.root");
444  fullPath = myDataFile.fullPath();
445  theFile = TFile::Open(fullPath.c_str());
446 
447  // Open the trees
448  unsigned fileNb = 0;
449  for (unsigned iname = 0; iname < theHadronNA.size(); ++iname) {
450  for (unsigned iene = 0; iene < theHadronEN.size(); ++iene) {
451  std::ostringstream filename;
452  double theEne = theHadronEN[iene];
453  filename << "NuclearInteractionsVal_" << theHadronNA[iname] << "_E" << theEne << ".root";
454  theFileNames[iname][iene] = filename.str();
455 
456  ++fileNb;
457  std::string treeName = "NuclearInteractions_" + theHadronNA[iname] + "_E" + std::to_string(int(theEne));
458  theTrees[iname][iene] = (TTree*)theFile->Get(treeName.c_str());
459 
460  if (!theTrees[iname][iene])
461  throw cms::Exception("FastSimulation/MaterialEffects") << "Tree with name " << treeName << " not found ";
462  theBranches[iname][iene] = theTrees[iname][iene]->GetBranch("nuEvent");
463  if (!theBranches[iname][iene])
464  throw cms::Exception("FastSimulation/MaterialEffects")
465  << "Branch with name nuEvent not found in " << theFileNames[iname][iene];
466 
467  theNUEvents[iname][iene] = new NUEvent();
468  theBranches[iname][iene]->SetAddress(&theNUEvents[iname][iene]);
469  theNumberOfEntries[iname][iene] = theTrees[iname][iene]->GetEntries();
470 
471  if (currentValuesWereSet) {
472  theTrees[iname][iene]->GetEntry(theCurrentEntry[iname][iene]);
473  unsigned NInteractions = theNUEvents[iname][iene]->nInteractions();
474  theNumberOfInteractions[iname][iene] = NInteractions;
475  }
476 
477  // Compute the corresponding cm energies of the nuclear interactions
478  XYZTLorentzVector Proton(0., 0., 0., 0.986);
479  XYZTLorentzVector Reference(
480  0.,
481  0.,
482  std::sqrt(theHadronEN[iene] * theHadronEN[iene] - theHadronMA[iname] * theHadronMA[iname]),
483  theHadronEN[iene]);
484  theHadronCM[iname][iene] = (Reference + Proton).M();
485  }
486  }
487 
488  // Find the index for which EN = 4. (or thereabout)
489  ien4 = 0;
490  while (theHadronEN[ien4] < 4.0)
491  ++ien4;
492 
493  gROOT->cd();
494 }
495 
497  // Close all local files
498  // Among other things, this allows the TROOT destructor to end up
499  // without crashing, while trying to close these files from outside
500  theFile->Close();
501  delete theFile;
502 
503  for (auto& vEvents : theNUEvents) {
504  for (auto evtPtr : vEvents) {
505  delete evtPtr;
506  }
507  }
508 
509  // Save data
510  save();
511  // Close the output file
512  myOutputFile.close();
513 }
514 
516  const SimplifiedGeometry& layer,
517  std::vector<std::unique_ptr<fastsim::Particle> >& secondaries,
518  const RandomEngineAndDistribution& random) {
519  int pdgId = particle.pdgId();
520  //
521  // no valid PDGid
522  //
523  if (abs(pdgId) <= 100 || abs(pdgId) >= 1000000) {
524  return;
525  }
526 
527  double radLengths = layer.getThickness(particle.position(), particle.momentum());
528  // TEC layers have some fudge factor for nuclear interactions
529  radLengths *= layer.getNuclearInteractionThicknessFactor();
530  //
531  // no material
532  //
533  if (radLengths < 1E-10) {
534  return;
535  }
536 
537  // In case the events are not read from (old) saved file, then pick a random event from FullSim file
538  if (!currentValuesWereSet) {
539  currentValuesWereSet = true;
540  for (unsigned iname = 0; iname < theHadronNA.size(); ++iname) {
541  for (unsigned iene = 0; iene < theHadronEN.size(); ++iene) {
542  theCurrentEntry[iname][iene] = (unsigned)(theNumberOfEntries[iname][iene] * random.flatShoot());
543 
544  theTrees[iname][iene]->GetEntry(theCurrentEntry[iname][iene]);
545  unsigned NInteractions = theNUEvents[iname][iene]->nInteractions();
546  theNumberOfInteractions[iname][iene] = NInteractions;
547 
548  theCurrentInteraction[iname][iene] = (unsigned)(theNumberOfInteractions[iname][iene] * random.flatShoot());
549  }
550  }
551  }
552 
553  // Momentum of interacting hadron
554  double pHadron = std::sqrt(particle.momentum().Vect().Mag2());
555 
556  //
557  // The hadron has not enough momentum to create some relevant final state
558  //
559  if (pHadron < theHadronEnergy) {
560  return;
561  }
562 
563  // The particle type
564  std::map<int, int>::const_iterator thePit = theIDMap.find(pdgId);
565  // Use map for unique ID (e.g. proton = {2212, 3222, -101, -102, -103, -104})
566  int thePid = (thePit != theIDMap.end() ? thePit->second : pdgId);
567 
568  // Is this particle type foreseen?
569  unsigned fPid = abs(thePid);
570  if (fPid != 211 && fPid != 130 && fPid != 321 && fPid != 2112 && fPid != 2212) {
571  return;
572  }
573 
574  // The hashed ID
575  unsigned thePidIndex = index(thePid);
576  // The inelastic interaction length at p(Hadron) = 5 GeV/c
577  double theInelasticLength = radLengths * theLengthRatio[thePidIndex];
578 
579  // The elastic interaction length
580  // The baroque parameterization is a fit to Fig. 40.13 of the PDG
581  double ee = pHadron > 0.6 ? exp(-std::sqrt((pHadron - 0.6) / 1.122)) : exp(std::sqrt((0.6 - pHadron) / 1.122));
582  double theElasticLength = (0.8753 * ee + 0.15) * theInelasticLength;
583 
584  // The total interaction length
585  double theTotalInteractionLength = theInelasticLength + theElasticLength;
586 
587  //
588  // Probability to interact is dl/L0 (maximum for 4 GeV Hadron)
589  //
590  double aNuclInteraction = -std::log(random.flatShoot());
591  if (aNuclInteraction > theTotalInteractionLength) {
592  return;
593  }
594 
595  // The elastic part
596  double elastic = random.flatShoot();
597  if (elastic < theElasticLength / theTotalInteractionLength) {
598  // Characteristic scattering angle for the elastic part
599  // A of silicon
600  double A = 28.0855;
601  double theta0 = std::sqrt(3.) / std::pow(A, 1. / 3.) * particle.momentum().mass() / pHadron;
602 
603  // Draw an angle with theta/theta0*exp[(-theta/2theta0)**2] shape
604  double theta = theta0 * std::sqrt(-2. * std::log(random.flatShoot()));
605  double phi = 2. * M_PI * random.flatShoot();
606 
607  // Rotate the particle accordingly
608  ROOT::Math::AxisAngle rotation1(orthogonal(particle.momentum().Vect()), theta);
609  ROOT::Math::AxisAngle rotation2(particle.momentum().Vect(), phi);
610  // Rotate!
611  XYZVector rotated = rotation2((rotation1(particle.momentum().Vect())));
612 
613  // Create a daughter if the kink is large enough
614  if (std::sin(theta) > theDistCut) {
615  secondaries.emplace_back(
617  particle.position(),
618  XYZTLorentzVector(rotated.X(), rotated.Y(), rotated.Z(), particle.momentum().E())));
619 
620  // Set the ClosestCharged Daughter thing for tracking
621  if (particle.charge() != 0) {
622  secondaries.back()->setMotherDeltaR(particle.momentum());
623  secondaries.back()->setMotherPdgId(pdgId);
624  secondaries.back()->setMotherSimTrackIndex(particle.simTrackIndex());
625  }
626 
627  // The particle is destroyed
628  particle.momentum().SetXYZT(0., 0., 0., 0.);
629  } else {
630  // If kink is small enough just rotate particle
631  particle.momentum().SetXYZT(rotated.X(), rotated.Y(), rotated.Z(), particle.momentum().E());
632  }
633  // The inelastic part
634  } else {
635  // Avoid multiple map access
636  const std::vector<double>& aHadronCM = theHadronCM[thePidIndex];
637  const std::vector<double>& aRatios = theRatiosMap[thePidIndex];
638  // Find the file with the closest c.m energy
639  // The target nucleon
640  XYZTLorentzVector Proton(0., 0., 0., 0.939);
641  // The current particle
642  const XYZTLorentzVector& Hadron = (const XYZTLorentzVector&)particle.momentum();
643  // The smallest momentum for inelastic interactions
644  double pMin = theHadronPMin[thePidIndex];
645  // The corresponding smallest four vector
646  XYZTLorentzVector Hadron0(0., 0., pMin, std::sqrt(pMin * pMin + particle.momentum().M2()));
647 
648  // The current centre-of-mass energy
649  double ecm = (Proton + Hadron).M();
650 
651  // Get the files of interest (closest c.m. energies)
652  unsigned ene1 = 0;
653  unsigned ene2 = 0;
654  // The smallest centre-of-mass energy
655  double ecm1 = (Proton + Hadron0).M();
656 
657  double ecm2 = aHadronCM[0];
658  double ratio1 = 0.;
659  double ratio2 = aRatios[0];
660  if (ecm > aHadronCM[0] && ecm < aHadronCM[aHadronCM.size() - 1]) {
661  for (unsigned ene = 1; ene < aHadronCM.size() && ecm > aHadronCM[ene - 1]; ++ene) {
662  if (ecm < aHadronCM[ene]) {
663  ene2 = ene;
664  ene1 = ene2 - 1;
665  ecm1 = aHadronCM[ene1];
666  ecm2 = aHadronCM[ene2];
667  ratio1 = aRatios[ene1];
668  ratio2 = aRatios[ene2];
669  }
670  }
671  } else if (ecm > aHadronCM[aHadronCM.size() - 1]) {
672  ene1 = aHadronCM.size() - 1;
673  ene2 = aHadronCM.size() - 2;
674  ecm1 = aHadronCM[ene1];
675  ecm2 = aHadronCM[ene2];
676  ratio1 = aRatios[ene2];
677  ratio2 = aRatios[ene2];
678  }
679 
680  // The inelastic part of the cross section depends cm energy
681  double slope = (std::log10(ecm) - std::log10(ecm1)) / (std::log10(ecm2) - std::log10(ecm1));
682  double inelastic = ratio1 + (ratio2 - ratio1) * slope;
683  double inelastic4 = pHadron < 4. ? aRatios[ien4] : 1.;
684 
685  // Simulate an inelastic interaction
686  if (elastic > 1. - (inelastic * theInelasticLength) / theTotalInteractionLength) {
687  // Avoid mutliple map access
688  std::vector<unsigned>& aCurrentInteraction = theCurrentInteraction[thePidIndex];
689  std::vector<unsigned>& aNumberOfInteractions = theNumberOfInteractions[thePidIndex];
690  std::vector<NUEvent*>& aNUEvents = theNUEvents[thePidIndex];
691 
692  // Choice of the file to read according the the log10(ecm) distance
693  // and protection against low momentum proton and neutron that never interacts
694  // (i.e., empty files)
695  unsigned ene;
696  if (random.flatShoot() < slope || aNumberOfInteractions[ene1] == 0) {
697  ene = ene2;
698  } else {
699  ene = ene1;
700  }
701 
702  // The boost characteristics
703  XYZTLorentzVector theBoost = Proton + Hadron;
704  theBoost /= theBoost.e();
705 
706  // Check we are not either at the end of an interaction bunch
707  // or at the end of a file
708  if (aCurrentInteraction[ene] == aNumberOfInteractions[ene]) {
709  std::vector<unsigned>& aCurrentEntry = theCurrentEntry[thePidIndex];
710  std::vector<unsigned>& aNumberOfEntries = theNumberOfEntries[thePidIndex];
711  std::vector<TTree*>& aTrees = theTrees[thePidIndex];
712  ++aCurrentEntry[ene];
713 
714  aCurrentInteraction[ene] = 0;
715  if (aCurrentEntry[ene] == aNumberOfEntries[ene]) {
716  aCurrentEntry[ene] = 0;
717  }
718 
719  unsigned myEntry = aCurrentEntry[ene];
720  aTrees[ene]->GetEntry(myEntry);
721  aNumberOfInteractions[ene] = aNUEvents[ene]->nInteractions();
722  }
723 
724  // Read the interaction
725  NUEvent::NUInteraction anInteraction = aNUEvents[ene]->theNUInteractions()[aCurrentInteraction[ene]];
726 
727  unsigned firstTrack = anInteraction.first;
728  unsigned lastTrack = anInteraction.last;
729 
730  // Some rotation around the boost axis, for more randomness
731  XYZVector theAxis = theBoost.Vect().Unit();
732  double theAngle = random.flatShoot() * 2. * M_PI;
733  ROOT::Math::AxisAngle axisRotation(theAxis, theAngle);
734  ROOT::Math::Boost axisBoost(theBoost.x(), theBoost.y(), theBoost.z());
735 
736  // A rotation to bring the particles back to the Hadron direction
737  XYZVector zAxis(0., 0., 1.);
738  XYZVector orthAxis = (zAxis.Cross(theBoost.Vect())).Unit();
739  double orthAngle = acos(theBoost.Vect().Unit().Z());
740  ROOT::Math::AxisAngle orthRotation(orthAxis, orthAngle);
741 
742  // Loop on the nuclear interaction products
743  for (unsigned iTrack = firstTrack; iTrack <= lastTrack; ++iTrack) {
744  NUEvent::NUParticle aParticle = aNUEvents[ene]->theNUParticles()[iTrack];
745 
746  // Add a RawParticle with the proper energy in the c.m. frame of
747  // the nuclear interaction
748  double energy = std::sqrt(aParticle.px * aParticle.px + aParticle.py * aParticle.py +
749  aParticle.pz * aParticle.pz + aParticle.mass * aParticle.mass / (ecm * ecm));
750 
751  XYZTLorentzVector daugtherMomentum(aParticle.px * ecm, aParticle.py * ecm, aParticle.pz * ecm, energy * ecm);
752 
753  // Rotate to the collision axis
754  XYZVector rotated = orthRotation(daugtherMomentum.Vect());
755  // Rotate around the boost axis for more randomness
756  rotated = axisRotation(rotated);
757 
758  // Rotated the daughter
759  daugtherMomentum.SetXYZT(rotated.X(), rotated.Y(), rotated.Z(), daugtherMomentum.E());
760 
761  // Boost it in the lab frame
762  daugtherMomentum = axisBoost(daugtherMomentum);
763 
764  // Create secondary
765  secondaries.emplace_back(new fastsim::Particle(aParticle.id, particle.position(), daugtherMomentum));
766 
767  // The closestCharged Daughter thing for tracking
768  // BUT: 'aParticle' also has to be charged, only then the mother should be set
769  // Unfortunately, NUEvent::NUParticle does not contain any info about the charge
770  // Did some tests and effect is absolutely negligible!
771  if (particle.charge() != 0) {
772  secondaries.back()->setMotherDeltaR(particle.momentum());
773  secondaries.back()->setMotherPdgId(pdgId);
774  secondaries.back()->setMotherSimTrackIndex(particle.simTrackIndex());
775  }
776  }
777 
778  // The particle is destroyed
779  particle.momentum().SetXYZT(0., 0., 0., 0.);
780 
781  // This is a note from previous version of code but I don't understand it:
782  // ERROR The way this loops through the events breaks
783  // replay. Which events are retrieved depends on
784  // which previous events were processed.
785 
786  // Increment for next time
787  ++aCurrentInteraction[ene];
788 
789  // Simulate a stopping hadron (low momentum)
790  } else if (pHadron < 4. && elastic > 1. - (inelastic4 * theInelasticLength) / theTotalInteractionLength) {
791  // The particle is destroyed
792  particle.momentum().SetXYZT(0., 0., 0., 0.);
793  }
794  }
795 }
796 
798  // Size of buffer
799  ++myOutputBuffer;
800 
801  // Periodically close the current file and open a new one
802  if (myOutputBuffer / 1000 * 1000 == myOutputBuffer) {
803  myOutputFile.close();
804  myOutputFile.open("NuclearInteractionOutputFile.txt");
805  }
806  //
807  unsigned size1 = theCurrentEntry.size() * theCurrentEntry.begin()->size();
808  std::vector<unsigned> theCurrentEntries;
809  theCurrentEntries.resize(size1);
810  size1 *= sizeof(unsigned);
811  //
812  unsigned size2 = theCurrentInteraction.size() * theCurrentInteraction.begin()->size();
813  std::vector<unsigned> theCurrentInteractions;
814  theCurrentInteractions.resize(size2);
815  size2 *= sizeof(unsigned);
816 
817  // Save the current entries
818  std::vector<std::vector<unsigned> >::const_iterator aCurrentEntry = theCurrentEntry.begin();
819  std::vector<std::vector<unsigned> >::const_iterator lastCurrentEntry = theCurrentEntry.end();
820  unsigned allEntries = 0;
821  for (; aCurrentEntry != lastCurrentEntry; ++aCurrentEntry) {
822  unsigned size = aCurrentEntry->size();
823  for (unsigned iene = 0; iene < size; ++iene)
824  theCurrentEntries[allEntries++] = (*aCurrentEntry)[iene];
825  }
826 
827  // Save the current interactions
828  std::vector<std::vector<unsigned> >::const_iterator aCurrentInteraction = theCurrentInteraction.begin();
829  std::vector<std::vector<unsigned> >::const_iterator lastCurrentInteraction = theCurrentInteraction.end();
830  unsigned allInteractions = 0;
831  for (; aCurrentInteraction != lastCurrentInteraction; ++aCurrentInteraction) {
832  unsigned size = aCurrentInteraction->size();
833  for (unsigned iene = 0; iene < size; ++iene)
834  theCurrentInteractions[allInteractions++] = (*aCurrentInteraction)[iene];
835  }
836  // Write to file
837  myOutputFile.write((const char*)(&theCurrentEntries.front()), size1);
838  myOutputFile.write((const char*)(&theCurrentInteractions.front()), size2);
839  myOutputFile.flush();
840 }
841 
843  std::ifstream myInputFile;
844  struct stat results;
845  //
846  unsigned size1 = theCurrentEntry.size() * theCurrentEntry.begin()->size();
847  std::vector<unsigned> theCurrentEntries;
848  theCurrentEntries.resize(size1);
849  size1 *= sizeof(unsigned);
850  //
851  unsigned size2 = theCurrentInteraction.size() * theCurrentInteraction.begin()->size();
852  std::vector<unsigned> theCurrentInteractions;
853  theCurrentInteractions.resize(size2);
854  size2 *= sizeof(unsigned);
855  //
856  unsigned size = 0;
857 
858  // Open the file (if any), otherwise return false
859  myInputFile.open(inputFile.c_str());
860  if (myInputFile.is_open()) {
861  // Get the size of the file
862  if (stat(inputFile.c_str(), &results) == 0)
863  size = results.st_size;
864  else
865  return false; // Something is wrong with that file !
866 
867  // Position the pointer just before the last record
868  myInputFile.seekg(size - size1 - size2);
869  myInputFile.read((char*)(&theCurrentEntries.front()), size1);
870  myInputFile.read((char*)(&theCurrentInteractions.front()), size2);
871  myInputFile.close();
872 
873  // Read the current entries
874  std::vector<std::vector<unsigned> >::iterator aCurrentEntry = theCurrentEntry.begin();
875  std::vector<std::vector<unsigned> >::iterator lastCurrentEntry = theCurrentEntry.end();
876  unsigned allEntries = 0;
877  for (; aCurrentEntry != lastCurrentEntry; ++aCurrentEntry) {
878  unsigned size = aCurrentEntry->size();
879  for (unsigned iene = 0; iene < size; ++iene)
880  (*aCurrentEntry)[iene] = theCurrentEntries[allEntries++];
881  }
882 
883  // Read the current interactions
884  std::vector<std::vector<unsigned> >::iterator aCurrentInteraction = theCurrentInteraction.begin();
885  std::vector<std::vector<unsigned> >::iterator lastCurrentInteraction = theCurrentInteraction.end();
886  unsigned allInteractions = 0;
887  for (; aCurrentInteraction != lastCurrentInteraction; ++aCurrentInteraction) {
888  unsigned size = aCurrentInteraction->size();
889  for (unsigned iene = 0; iene < size; ++iene)
890  (*aCurrentInteraction)[iene] = theCurrentInteractions[allInteractions++];
891  }
892 
893  return true;
894  }
895 
896  return false;
897 }
898 
900  // Find hashed particle ID
901  unsigned myIndex = 0;
902  while (thePid != theHadronID[myIndex])
903  ++myIndex;
904  return myIndex;
905 }
906 
908  double x = fabs(aVector.X());
909  double y = fabs(aVector.Y());
910  double z = fabs(aVector.Z());
911 
912  if (x < y)
913  return x < z ? XYZVector(0., aVector.Z(), -aVector.Y()) : XYZVector(aVector.Y(), -aVector.X(), 0.);
914  else
915  return y < z ? XYZVector(-aVector.Z(), 0., aVector.X()) : XYZVector(aVector.Y(), -aVector.X(), 0.);
916 }
917 
DDAxes::y
NUEvent::NUParticle
Definition: NUEvent.h:14
fastsim::NuclearInteraction::theCurrentEntry
std::vector< std::vector< unsigned > > theCurrentEntry
Necessary to read the FullSim interactions.
Definition: NuclearInteraction.cc:103
fastsim::NuclearInteraction::theDistCut
double theDistCut
Cut on deltaR for the FastSim Tracking (ClosestChargedDaughter algorithm)
Definition: NuclearInteraction.cc:91
fastsim::Particle::momentum
const math::XYZTLorentzVector & momentum() const
Return momentum of the particle.
Definition: Particle.h:143
mps_fire.i
i
Definition: mps_fire.py:428
fastsim::NuclearInteraction::currentValuesWereSet
bool currentValuesWereSet
Read data from file that was created in a previous run.
Definition: NuclearInteraction.cc:115
MessageLogger.h
funct::false
false
Definition: Factorize.h:29
InteractionModel.h
fastsim::NuclearInteraction::theCurrentInteraction
std::vector< std::vector< unsigned > > theCurrentInteraction
Necessary to read the FullSim interactions.
Definition: NuclearInteraction.cc:104
contentValuesFiles.fullPath
fullPath
Definition: contentValuesFiles.py:64
gather_cfg.cout
cout
Definition: gather_cfg.py:144
NUEvent::NUParticle::py
float py
Definition: NUEvent.h:18
fastsim::SimplifiedGeometry
Implementation of a generic detector layer (base class for forward/barrel layers).
Definition: SimplifiedGeometry.h:35
XYZVector
math::XYZVector XYZVector
Definition: RawParticle.h:26
NUEvent::NUParticle::px
float px
Definition: NUEvent.h:17
fastsim::NuclearInteraction::theFile
TFile * theFile
Necessary to read the FullSim interactions.
Definition: NuclearInteraction.cc:99
fastsim::NuclearInteraction::theLengthRatio
const std::vector< double > theLengthRatio
Inelastic interaction length at p(Hadron) = 5 GeV/c (relativ to radionLength of material)
Definition: NuclearInteraction.cc:143
bookConverter.results
results
Definition: bookConverter.py:144
fastsim::NuclearInteraction::K0LsID
const std::vector< int > K0LsID
PdgID of K0.
Definition: NuclearInteraction.cc:334
DDAxes::x
fastsim::Particle::simTrackIndex
int simTrackIndex() const
Return index of the SimTrack.
Definition: Particle.h:153
XYZTLorentzVector
math::XYZTLorentzVector XYZTLorentzVector
Definition: RawParticle.h:25
fastsim::NuclearInteraction::theNUEvents
std::vector< std::vector< NUEvent * > > theNUEvents
Necessary to read the FullSim interactions.
Definition: NuclearInteraction.cc:102
fastsim::NuclearInteraction::theNumberOfEntries
std::vector< std::vector< unsigned > > theNumberOfEntries
Necessary to read the FullSim interactions.
Definition: NuclearInteraction.cc:105
NUEvent::NUParticle::mass
float mass
Definition: NUEvent.h:20
fastsim::NuclearInteraction::KminussesID
const std::vector< int > KminussesID
PdgID of K-.
Definition: NuclearInteraction.cc:336
RandomEngineAndDistribution.h
InteractionModelFactory.h
Particle.h
fastsim::NuclearInteraction::protonsID
const std::vector< int > protonsID
PdgID of protons.
Definition: NuclearInteraction.cc:330
funct::sin
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
fastsim::InteractionModel
Base class for any interaction model between a particle and a tracker layer.
Definition: InteractionModel.h:29
edm::FileInPath
Definition: FileInPath.h:61
fastsim::NuclearInteraction::theHadronCM
std::vector< std::vector< double > > theHadronCM
Necessary to read the FullSim interactions.
Definition: NuclearInteraction.cc:108
interactiveExample.theFile
theFile
Definition: interactiveExample.py:17
NUEvent::NUInteraction::last
unsigned last
Definition: NUEvent.h:28
fastsim::NuclearInteraction::KplussesID
const std::vector< int > KplussesID
PdgID of K+.
Definition: NuclearInteraction.cc:335
fastsim::NuclearInteraction::theFileNames
std::vector< std::vector< std::string > > theFileNames
Necessary to read the FullSim interactions.
Definition: NuclearInteraction.cc:107
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
NUEvent::NUParticle::id
int id
Definition: NUEvent.h:21
DDAxes::z
fastsim::NuclearInteraction::theRatiosMap
static std::vector< std::vector< double > > theRatiosMap
The evolution of the interaction lengths with energy.
Definition: NuclearInteraction.cc:122
fastsim::NuclearInteraction::myOutputFile
std::ofstream myOutputFile
Output file to save interactions.
Definition: NuclearInteraction.cc:112
HLT_FULL_cff.zAxis
zAxis
Definition: HLT_FULL_cff.py:46178
fastsim::NuclearInteraction::theRatios
const std::vector< double > theRatios
Filled into 'theRatiosMap' (evolution of the interaction lengths with energy)
Definition: NuclearInteraction.cc:154
corrVsCorr.filename
filename
Definition: corrVsCorr.py:123
theta
Geom::Theta< T > theta() const
Definition: Basic3DVectorLD.h:150
HCALHighEnergyHPDFilter_cfi.energy
energy
Definition: HCALHighEnergyHPDFilter_cfi.py:5
fastsim::NuclearInteraction::neutronsID
const std::vector< int > neutronsID
PdgID of neutrons.
Definition: NuclearInteraction.cc:332
fastsim::NuclearInteraction::save
void save()
Save the nuclear interactions to a file, so you can reproduce the events (e.g. in case of a crash).
Definition: NuclearInteraction.cc:797
DEFINE_EDM_PLUGIN
#define DEFINE_EDM_PLUGIN(factory, type, name)
Definition: PluginFactory.h:124
fastsim::Particle::position
const math::XYZTLorentzVector & position() const
Return position of the particle.
Definition: Particle.h:140
fastsim::NuclearInteraction::NuclearInteraction
NuclearInteraction(const std::string &name, const edm::ParameterSet &cfg)
Constructor.
Definition: NuclearInteraction.cc:348
phase1PixelTopology::layer
constexpr std::array< uint8_t, layerIndexSize > layer
Definition: phase1PixelTopology.h:99
FileInPath.h
fastsim::NuclearInteraction::theIDMap
static std::map< int, int > theIDMap
Build the ID map (i.e., what is to be considered as a proton, etc...)
Definition: NuclearInteraction.cc:124
edm::ParameterSet
Definition: ParameterSet.h:47
Constants.h
fastsim::NuclearInteraction::interact
void interact(fastsim::Particle &particle, const SimplifiedGeometry &layer, std::vector< std::unique_ptr< fastsim::Particle > > &secondaries, const RandomEngineAndDistribution &random) override
Perform the interaction.
Definition: NuclearInteraction.cc:515
ALCARECOTkAlMinBias_cff.pMin
pMin
GeV.
Definition: ALCARECOTkAlMinBias_cff.py:35
fastsim::NuclearInteraction::antiprotonsID
const std::vector< int > antiprotonsID
PdgID of anti-protons.
Definition: NuclearInteraction.cc:331
edmplugin::PluginFactory
Definition: PluginFactory.h:34
fastsim::NuclearInteraction::theHadronEnergy
double theHadronEnergy
Minimum energy for nuclear interaction.
Definition: NuclearInteraction.cc:92
XYZVector
math::XYZVector XYZVector
Definition: NuclearInteraction.cc:49
fastsim::Particle::pdgId
int pdgId() const
Return pdgId of the particle.
Definition: Particle.h:134
thread_safety_macros.h
NUEvent::NUInteraction::first
unsigned first
Definition: NUEvent.h:27
math::XYZVector
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
A
LorentzVector.h
dtResolutionTest_cfi.inputFile
inputFile
Definition: dtResolutionTest_cfi.py:14
fastsim::Particle::charge
double charge() const
Return charge of the particle.
Definition: Particle.h:137
fastsim::NuclearInteraction::~NuclearInteraction
~NuclearInteraction() override
Default destructor.
Definition: NuclearInteraction.cc:496
fastsim::Particle
Definition of a generic FastSim Particle which can be propagated through the detector (formerly Parti...
Definition: Particle.h:16
M_PI
#define M_PI
Definition: BXVectorInputProducer.cc:49
fastsim::NuclearInteraction
Implementation of nuclear interactions of hadrons in the tracker layers (based on fully simulated int...
Definition: NuclearInteraction.cc:58
trackerHitRTTI::vector
Definition: trackerHitRTTI.h:21
fastsim::NuclearInteraction::orthogonal
XYZVector orthogonal(const XYZVector &aVector) const
Return an orthogonal vector.
Definition: NuclearInteraction.cc:907
EgammaValidation_cff.pdgId
pdgId
Definition: EgammaValidation_cff.py:117
XYZTLorentzVector
math::XYZTLorentzVector XYZTLorentzVector
Definition: NuclearInteraction.cc:50
fastsim::NuclearInteraction::index
unsigned index(int thePid)
Return a hashed index for a given particle ID.
Definition: NuclearInteraction.cc:899
fastsim::initializeOnce
static std::once_flag initializeOnce
Definition: NuclearInteraction.cc:343
fastsim::NuclearInteraction::theHadronID
const std::vector< int > theHadronID
ID of the hadrons.
Definition: NuclearInteraction.cc:133
NUEvent.h
RandomEngineAndDistribution::flatShoot
double flatShoot(double xmin=0.0, double xmax=1.0) const
Definition: RandomEngineAndDistribution.h:27
SimplifiedGeometry.h
AlCaHLTBitMon_QueryRunRegistry.string
string string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
looper.cfg
cfg
Definition: looper.py:296
fastsim::NuclearInteraction::inputFile
std::string inputFile
Directory/Name of input file in case you want to read interactions from file.
Definition: NuclearInteraction.cc:93
fastsim::NuclearInteraction::theHadronMA
const std::vector< double > theHadronMA
Masses of the hadrons.
Definition: NuclearInteraction.cc:138
fastsim::NuclearInteraction::theHadronNA
const std::vector< std::string > theHadronNA
Names of the hadrons.
Definition: NuclearInteraction.cc:135
triggerObjects_cff.id
id
Definition: triggerObjects_cff.py:29
math::XYZTLorentzVector
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
Exception
Definition: hltDiff.cc:245
NUEvent::NUParticle::pz
float pz
Definition: NUEvent.h:19
fastsim::NuclearInteraction::ien4
unsigned ien4
Find the index for which EN = 4.
Definition: NuclearInteraction.cc:110
NUEvent::NUInteraction
Definition: NUEvent.h:24
fastsim::NuclearInteraction::PiminussesID
const std::vector< int > PiminussesID
PdgID of pi-.
Definition: NuclearInteraction.cc:338
Skims_PA_cff.name
name
Definition: Skims_PA_cff.py:17
fastsim::NuclearInteraction::antineutronsID
const std::vector< int > antineutronsID
PdgID of anti-neutrons.
Definition: NuclearInteraction.cc:333
Exception.h
fastsim::NuclearInteraction::read
bool read(std::string inputFile)
Read the nuclear interactions from a file, so you can reproduce the events (e.g. in case of a crash).
Definition: NuclearInteraction.cc:842
dqm-mbProfile.log
log
Definition: dqm-mbProfile.py:17
AlignmentPI::index
index
Definition: AlignmentPayloadInspectorHelper.h:46
fastsim::NuclearInteraction::theHadronPMin
const std::vector< double > theHadronPMin
Smallest momentum for inelastic interactions.
Definition: NuclearInteraction.cc:141
fastsim::NuclearInteraction::theNumberOfInteractions
std::vector< std::vector< unsigned > > theNumberOfInteractions
Necessary to read the FullSim interactions.
Definition: NuclearInteraction.cc:106
funct::pow
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
fastsim
Definition: BarrelSimplifiedGeometry.h:15
cuy.save
save
Definition: cuy.py:1164
ParameterSet.h
dqmiolumiharvest.j
j
Definition: dqmiolumiharvest.py:66
JetChargeProducer_cfi.exp
exp
Definition: JetChargeProducer_cfi.py:6
NUEvent
Definition: NUEvent.h:6
slope
static const double slope[3]
Definition: CastorTimeSlew.cc:6
edm_modernize_messagelogger.stat
stat
Definition: edm_modernize_messagelogger.py:27
Vector3D.h
fastsim::NuclearInteraction::theBranches
std::vector< std::vector< TBranch * > > theBranches
Necessary to read the FullSim interactions.
Definition: NuclearInteraction.cc:101
fastsim::NuclearInteraction::myOutputBuffer
unsigned myOutputBuffer
Needed to save interactions to file.
Definition: NuclearInteraction.cc:113
CMS_THREAD_GUARD
#define CMS_THREAD_GUARD(_var_)
Definition: thread_safety_macros.h:6
fastsim::NuclearInteraction::theTrees
std::vector< std::vector< TTree * > > theTrees
Necessary to read the FullSim interactions.
Definition: NuclearInteraction.cc:100
fastsim::NuclearInteraction::theHadronEN
const std::vector< double > theHadronEN
Filled into 'theRatiosMap' (evolution of the interaction lengths with energy)
Definition: NuclearInteraction.cc:127
fastsim::NuclearInteraction::PiplussesID
const std::vector< int > PiplussesID
PdgID of pt+.
Definition: NuclearInteraction.cc:337
edm::FileInPath::fullPath
std::string fullPath() const
Definition: FileInPath.cc:161
findQualityFiles.size
size
Write out results.
Definition: findQualityFiles.py:443
RandomEngineAndDistribution
Definition: RandomEngineAndDistribution.h:18