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