CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
NuclearInteractionSimulator.cc
Go to the documentation of this file.
1 //Framework Headers
4 
7 
9 
10 //#include "DQMServices/Core/interface/DaqMonitorBEInterface.h"
11 //#include "FWCore/ServiceRegistry/interface/Service.h"
12 
13 #include <iostream>
14 #include <sys/stat.h>
15 #include <cmath>
16 #include "TFile.h"
17 #include "TTree.h"
18 #include "TROOT.h"
19 
20 // Internal variable and vectors with name started frm "thePion" means
21 // vectors/variable not only for pions but for all type of hadrons
22 // treated inside this code
23 
25  std::vector<double>& hadronEnergies,
26  std::vector<int>& hadronTypes,
27  std::vector<std::string>& hadronNames,
28  std::vector<double>& hadronMasses,
29  std::vector<double>& hadronPMin,
30  double pionEnergy,
31  std::vector<double>& lengthRatio,
32  std::vector< std::vector<double> >& ratios,
33  std::map<int,int >& idMap,
35  unsigned int distAlgo,
36  double distCut)
37  :
39  thePionEN(hadronEnergies),
40  thePionID(hadronTypes),
41  thePionNA(hadronNames),
42  thePionMA(hadronMasses),
43  thePionPMin(hadronPMin),
44  thePionEnergy(pionEnergy),
45  theLengthRatio(lengthRatio),
46  theRatios(ratios),
47  theIDMap(idMap),
48  theDistAlgo(distAlgo),
49  theDistCut(distCut),
50  currentValuesWereSet(false)
51 {
53 
54  // Prepare the map of files
55  // Loop over the particle names
56  std::vector<TFile*> aVFile(thePionEN.size(),static_cast<TFile*>(0));
57  std::vector<TTree*> aVTree(thePionEN.size(),static_cast<TTree*>(0));
58  std::vector<TBranch*> aVBranch(thePionEN.size(),static_cast<TBranch*>(0));
59  std::vector<NUEvent*> aVNUEvents(thePionEN.size(),static_cast<NUEvent*>(0));
60  std::vector<unsigned> aVCurrentEntry(thePionEN.size(),static_cast<unsigned>(0));
61  std::vector<unsigned> aVCurrentInteraction(thePionEN.size(),static_cast<unsigned>(0));
62  std::vector<unsigned> aVNumberOfEntries(thePionEN.size(),static_cast<unsigned>(0));
63  std::vector<unsigned> aVNumberOfInteractions(thePionEN.size(),static_cast<unsigned>(0));
64  std::vector<std::string> aVFileName(thePionEN.size(),static_cast<std::string>(""));
65  std::vector<double> aVPionCM(thePionEN.size(),static_cast<double>(0));
66  theFiles.resize(thePionNA.size());
67  theTrees.resize(thePionNA.size());
68  theBranches.resize(thePionNA.size());
69  theNUEvents.resize(thePionNA.size());
70  theCurrentEntry.resize(thePionNA.size());
71  theCurrentInteraction.resize(thePionNA.size());
72  theNumberOfEntries.resize(thePionNA.size());
73  theNumberOfInteractions.resize(thePionNA.size());
74  theFileNames.resize(thePionNA.size());
75  thePionCM.resize(thePionNA.size());
76  for ( unsigned iname=0; iname<thePionNA.size(); ++iname ) {
77  theFiles[iname] = aVFile;
78  theTrees[iname] = aVTree;
79  theBranches[iname] = aVBranch;
80  theNUEvents[iname] = aVNUEvents;
81  theCurrentEntry[iname] = aVCurrentEntry;
82  theCurrentInteraction[iname] = aVCurrentInteraction;
83  theNumberOfEntries[iname] = aVNumberOfEntries;
84  theNumberOfInteractions[iname] = aVNumberOfInteractions;
85  theFileNames[iname] = aVFileName;
86  thePionCM[iname] = aVPionCM;
87  }
88 
89  // Read the information from a previous run (to keep reproducibility)
90  currentValuesWereSet = this->read(inputFile);
92  std::cout << "***WARNING*** You are reading nuclear-interaction information from the file "
93  << inputFile << " created in an earlier run."
94  << std::endl;
95 
96  // Open the file for saving the information of the current run
97  myOutputFile.open ("NuclearInteractionOutputFile.txt");
98  myOutputBuffer = 0;
99 
100 
101  // Open the root files
102  // for ( unsigned file=0; file<theFileNames.size(); ++file ) {
103  unsigned fileNb = 0;
104  for ( unsigned iname=0; iname<thePionNA.size(); ++iname ) {
105  for ( unsigned iene=0; iene<thePionEN.size(); ++iene ) {
106  //std::cout << "iname/iene " << iname << " " << iene << std::endl;
107  std::ostringstream filename;
108  double theEne = thePionEN[iene];
109  filename << "NuclearInteractionsVal_" << thePionNA[iname] << "_E"<< theEne << ".root";
110  theFileNames[iname][iene] = filename.str();
111  //std::cout << "thePid/theEne " << thePionID[iname] << " " << theEne << std::endl;
112 
113  edm::FileInPath myDataFile("FastSimulation/MaterialEffects/data/"+theFileNames[iname][iene]);
114  fullPath = myDataFile.fullPath();
115  // theFiles[file] = TFile::Open(theFileNames[file].c_str());
116  theFiles[iname][iene] = TFile::Open(fullPath.c_str());
117  if ( !theFiles[iname][iene] ) throw cms::Exception("FastSimulation/MaterialEffects")
118  << "File " << theFileNames[iname][iene] << " " << fullPath << " not found ";
119  ++fileNb;
120  //
121  theTrees[iname][iene] = (TTree*) theFiles[iname][iene]->Get("NuclearInteractions");
122  if ( !theTrees[iname][iene] ) throw cms::Exception("FastSimulation/MaterialEffects")
123  << "Tree with name NuclearInteractions not found in " << theFileNames[iname][iene];
124  //
125  theBranches[iname][iene] = theTrees[iname][iene]->GetBranch("nuEvent");
126  //std::cout << "The branch = " << theBranches[iname][iene] << std::endl;
127  if ( !theBranches[iname][iene] ) throw cms::Exception("FastSimulation/MaterialEffects")
128  << "Branch with name nuEvent not found in " << theFileNames[iname][iene];
129  //
130  theNUEvents[iname][iene] = new NUEvent();
131  //std::cout << "The branch = " << theBranches[iname][iene] << std::endl;
132  theBranches[iname][iene]->SetAddress(&theNUEvents[iname][iene]);
133  //
134  theNumberOfEntries[iname][iene] = theTrees[iname][iene]->GetEntries();
135 
137  theTrees[iname][iene]->GetEntry(theCurrentEntry[iname][iene]);
138  unsigned NInteractions = theNUEvents[iname][iene]->nInteractions();
139  theNumberOfInteractions[iname][iene] = NInteractions;
140  }
141 
142  //
143  // Compute the corresponding cm energies of the nuclear interactions
144  XYZTLorentzVector Proton(0.,0.,0.,0.986);
146  Reference(0.,
147  0.,
148  std::sqrt(thePionEN[iene]*thePionEN[iene]
149  -thePionMA[iname]*thePionMA[iname]),
150  thePionEN[iene]);
151  thePionCM[iname][iene] = (Reference+Proton).M();
152 
153  }
154 
155  }
156 
157  // Find the index for which EN = 4. (or thereabout)
158  ien4 = 0;
159  while ( thePionEN[ien4] < 4.0 ) ++ien4;
160 
161  gROOT->cd();
162 
163  // Information (Should be on LogInfo)
164 // std::cout << " ---> A total of " << fileNb
165 // << " nuclear-interaction files was sucessfully open" << std::endl;
166 
167  // dbe = edm::Service<DaqMonitorBEInterface>().operator->();
168  // htot = dbe->book1D("Total", "All particles",150,0.,150.);
169  // helas = dbe->book1D("Elastic", "Elastic interactions",150,0.,150.);
170  // hinel = dbe->book1D("Inelastic", "Inelastic interactions",150,0.,150.);
171  // hscatter = dbe->book1D("Scattering","Elastic Scattering angle",200,0.,2.);
172  // hscatter2 = dbe->book2D("Scattering2","Elastic Scattering angle vs p",100,0.,10.,200,0.,2.);
173  // hAfter = dbe->book1D("eAfter","Energy after collision",200,0.,4.);
174  // hAfter2 = dbe->book2D("eAfter2","Energy after collision",100,-2.5,2.5,100,0.,4);
175  // hAfter3 = dbe->book2D("eAfter3","Energy after collision",100,0.,1000.,100,0.,4);
176 
177 }
178 
180 
181  // Close all local files
182  // Among other things, this allows the TROOT destructor to end up
183  // without crashing, while trying to close these files from outside
184  for ( unsigned ifile=0; ifile<theFiles.size(); ++ifile ) {
185  for ( unsigned iene=0; iene<theFiles[ifile].size(); ++iene ) {
186  auto file = theFiles[ifile][iene];
187  if(file) {
188  // std::cout << "Closing file " << iene << " with name " << theFileNames[ifile][iene] << std::endl;
189  file->Close();
190  delete file;
191  }
192  }
193  }
194 
195  for(auto& vEvents: theNUEvents) {
196  for(auto evtPtr: vEvents) {
197  delete evtPtr;
198  }
199  }
200 
201  // Close the output file
202  myOutputFile.close();
203 
204  // dbe->save("test.root");
205 
206 }
207 
209 {
210  if(!currentValuesWereSet) {
211  currentValuesWereSet = true;
212  for ( unsigned iname=0; iname<thePionNA.size(); ++iname ) {
213  for ( unsigned iene=0; iene<thePionEN.size(); ++iene ) {
214  theCurrentEntry[iname][iene] = (unsigned) (theNumberOfEntries[iname][iene] * random->flatShoot());
215 
216  theTrees[iname][iene]->GetEntry(theCurrentEntry[iname][iene]);
217  unsigned NInteractions = theNUEvents[iname][iene]->nInteractions();
218  theNumberOfInteractions[iname][iene] = NInteractions;
219 
220  theCurrentInteraction[iname][iene] = (unsigned) (theNumberOfInteractions[iname][iene] * random->flatShoot());
221  }
222  }
223  }
224 
225  // Read a Nuclear Interaction in a random manner
226 
227  double pHadron = std::sqrt(Particle.Vect().Mag2());
228  // htot->Fill(pHadron);
229 
230  // The hadron has enough momentum to create some relevant final state
231  if ( pHadron > thePionEnergy ) {
232 
233  // The particle type
234  std::map<int,int>::const_iterator thePit = theIDMap.find(Particle.pid());
235 
236  int thePid = thePit != theIDMap.end() ? thePit->second : Particle.pid();
237 
238  // Is this particle type foreseen?
239  unsigned fPid = abs(thePid);
240  if ( fPid != 211 && fPid != 130 && fPid != 321 && fPid != 2112 && fPid != 2212 ) {
241  return;
242  //std::cout << "Unknown particle type = " << thePid << std::endl;
243  //thePid = 211;
244  }
245 
246  // The inelastic interaction length at p(pion) = 5 GeV/c
247  unsigned thePidIndex = index(thePid);
248  double theInelasticLength = radLengths * theLengthRatio[thePidIndex];
249 
250  // The elastic interaction length
251  // The baroque parameterization is a fit to Fig. 40.13 of the PDG
252  double ee = pHadron > 0.6 ?
253  exp(-std::sqrt((pHadron-0.6)/1.122)) : exp(std::sqrt((0.6-pHadron)/1.122));
254  double theElasticLength = ( 0.8753 * ee + 0.15 )
255  // double theElasticLength = ( 0.15 + 0.195 / log(pHadron/0.4) )
256  // double theElasticLength = ( 0.15 + 0.305 / log(pHadron/0.35) )
257  * theInelasticLength;
258 
259  // The total interaction length
260  double theTotalInteractionLength = theInelasticLength + theElasticLength;
261 
262  // Probability to interact is dl/L0 (maximum for 4 GeV pion)
263  double aNuclInteraction = -std::log(random->flatShoot());
264  if ( aNuclInteraction < theTotalInteractionLength ) {
265 
266  // The elastic part
267  double elastic = random->flatShoot();
268  if ( elastic < theElasticLength/theTotalInteractionLength ) {
269 
270  // helas->Fill(pHadron);
271 
272  // Characteristic scattering angle for the elastic part
273  double theta0 = std::sqrt(3.)/ std::pow(theA(),1./3.) * Particle.mass()/pHadron;
274 
275  // Draw an angle with theta/theta0*exp[(-theta/2theta0)**2] shape
276  double theta = theta0 * std::sqrt(-2.*std::log(random->flatShoot()));
277  double phi = 2. * 3.14159265358979323 * random->flatShoot();
278 
279  // Rotate the particle accordingly
280  RawParticle::Rotation rotation1(orthogonal(Particle.Vect()),theta);
281  RawParticle::Rotation rotation2(Particle.Vect(),phi);
282  Particle.rotate(rotation1);
283  Particle.rotate(rotation2);
284 
285  // Distance
286  double distance = std::sin(theta);
287 
288  // Create a daughter if the kink is large engough
289  if ( distance > theDistCut ) {
290  _theUpdatedState.resize(1);
291  _theUpdatedState[0].SetXYZT(Particle.Px(), Particle.Py(),
292  Particle.Pz(), Particle.E());
293  _theUpdatedState[0].setID(Particle.pid());
294  }
295 
296  // hscatter->Fill(myTheta);
297  // hscatter2->Fill(pHadron,myTheta);
298 
299  }
300 
301  // The inelastic part
302  else {
303 
304  // Avoid multiple map access
305  const std::vector<double>& aPionCM = thePionCM[thePidIndex];
306  const std::vector<double>& aRatios = theRatios[thePidIndex];
307  // Find the file with the closest c.m energy
308  // The target nucleon
309  XYZTLorentzVector Proton(0.,0.,0.,0.939);
310  // The current particle
311  const XYZTLorentzVector& Hadron = (const XYZTLorentzVector&)Particle;
312  // The smallest momentum for inelastic interactions
313  double pMin = thePionPMin[thePidIndex];
314  // The correspong smallest four vector
315  XYZTLorentzVector Hadron0(0.,0.,pMin,std::sqrt(pMin*pMin+Particle.M2()));
316 
317  // The current centre-of-mass energy
318  double ecm = (Proton+Hadron).M();
319  //std::cout << "Proton = " << Proton << std::endl;
320  //std::cout << "Hadron = " << Hadron << std::endl;
321  //std::cout << "ecm = " << ecm << std::endl;
322  // Get the files of interest (closest c.m. energies)
323  unsigned ene1=0;
324  unsigned ene2=0;
325  // The smallest centre-of-mass energy
326  // double ecm1=1.63;
327  double ecm1= (Proton+Hadron0).M();
328  //std::cout << "ecm1 = " << ecm1 << std::endl;
329  //std::cout << "ecm[0] = " << aPionCM[0] << std::endl;
330  //std::cout << "ecm[11] = " << aPionCM [ aPionCM.size()-1 ] << std::endl;
331  double ecm2=aPionCM[0];
332  double ratio1=0.;
333  double ratio2=aRatios[0];
334  if ( ecm > aPionCM[0] && ecm < aPionCM [ aPionCM.size()-1 ] ) {
335  for ( unsigned ene=1;
336  ene < aPionCM.size() && ecm > aPionCM[ene-1];
337  ++ene ) {
338  if ( ecm<aPionCM[ene] ) {
339  ene2 = ene;
340  ene1 = ene2-1;
341  ecm1 = aPionCM[ene1];
342  ecm2 = aPionCM[ene2];
343  ratio1 = aRatios[ene1];
344  ratio2 = aRatios[ene2];
345  }
346  }
347  } else if ( ecm > aPionCM[ aPionCM.size()-1 ] ) {
348  ene1 = aPionCM.size()-1;
349  ene2 = aPionCM.size()-2;
350  ecm1 = aPionCM[ene1];
351  ecm2 = aPionCM[ene2];
352  ratio1 = aRatios[ene2];
353  ratio2 = aRatios[ene2];
354  }
355 
356 
357  // The inelastic part of the cross section depends cm energy
358  double slope = (std::log10(ecm )-std::log10(ecm1))
359  / (std::log10(ecm2)-std::log10(ecm1));
360  double inelastic = ratio1 + (ratio2-ratio1) * slope;
361  double inelastic4 = pHadron < 4. ? aRatios[ien4] : 1.;
362 
363  //std::cout << "Inelastic = " << ratio1 << " " << ratio2 << " " << inelastic << std::endl;
364  // std::cout << "Energy/Inelastic : "
365  // << Hadron.e() << " " << inelastic << std::endl;
366 
367  // Simulate an inelastic interaction
368  if ( elastic > 1.- (inelastic*theInelasticLength)
369  /theTotalInteractionLength ) {
370 
371  // Avoid mutliple map access
372  std::vector<unsigned>& aCurrentInteraction = theCurrentInteraction[thePidIndex];
373  std::vector<unsigned>& aNumberOfInteractions = theNumberOfInteractions[thePidIndex];
374  std::vector<NUEvent*>& aNUEvents = theNUEvents[thePidIndex];
375  // hinel->Fill(pHadron);
376  // std::cout << "INELASTIC INTERACTION ! "
377  // << pHadron << " " << theInelasticLength << " "
378  // << inelastic * theInelasticLength << std::endl;
379  // Choice of the file to read according the the log10(ecm) distance
380  // and protection against low momentum proton and neutron that never interacts
381  // (i.e., empty files)
382  unsigned ene;
383  if ( random->flatShoot() < slope || aNumberOfInteractions[ene1] == 0 )
384  ene = ene2;
385  else
386  ene = ene1;
387 
388  //std::cout << "Ecm1/2 = " << ecm1 << " " << ecm2 << std::endl;
389  //std::cout << "Ratio1/2 = " << ratio1 << " " << ratio2 << std::endl;
390  //std::cout << "Ene = " << ene << " slope = " << slope << std::endl;
391 
392  //std::cout << "Pion energy = " << Hadron.E()
393  // << "File chosen " << theFileNames[thePidIndex][ene]
394  // << std::endl;
395 
396  // The boost characteristics
397  XYZTLorentzVector theBoost = Proton + Hadron;
398  theBoost /= theBoost.e();
399 
400  // std::cout << "File chosen : " << thePid << "/" << ene
401  // << " Current interaction = " << aCurrentInteraction[ene]
402  // << " Total interactions = " << aNumberOfInteractions[ene]
403  // << std::endl;
404  // theFiles[thePidIndex][ene]->cd();
405  // gDirectory->ls();
406 
407  // Check we are not either at the end of an interaction bunch
408  // or at the end of a file
409  if ( aCurrentInteraction[ene] == aNumberOfInteractions[ene] ) {
410  // std::cout << "End of interaction bunch ! ";
411  std::vector<unsigned>& aCurrentEntry = theCurrentEntry[thePidIndex];
412  std::vector<unsigned>& aNumberOfEntries = theNumberOfEntries[thePidIndex];
413  std::vector<TTree*>& aTrees = theTrees[thePidIndex];
414  ++aCurrentEntry[ene];
415  // std::cerr << "Read the next entry "
416  // << aCurrentEntry[ene] << std::endl;
417  aCurrentInteraction[ene] = 0;
418  if ( aCurrentEntry[ene] == aNumberOfEntries[ene] ) {
419  aCurrentEntry[ene] = 0;
420  // std::cout << "End of file - Rewind! " << std::endl;
421  }
422  unsigned myEntry = aCurrentEntry[ene];
423  // std::cout << "The new entry " << myEntry
424  // << " is read ... in TTree " << aTrees[ene] << " ";
425  aTrees[ene]->GetEntry(myEntry);
426  // std::cout << "The number of interactions in the new entry is ... ";
427  aNumberOfInteractions[ene] = aNUEvents[ene]->nInteractions();
428  // std::cout << aNumberOfInteractions[ene] << std::endl;
429  }
430 
431  // Read the interaction
432  NUEvent::NUInteraction anInteraction
433  = aNUEvents[ene]->theNUInteractions()[aCurrentInteraction[ene]];
434 
435  unsigned firstTrack = anInteraction.first;
436  unsigned lastTrack = anInteraction.last;
437  // std::cout << "First and last tracks are " << firstTrack << " " << lastTrack << std::endl;
438 
439  _theUpdatedState.resize(lastTrack-firstTrack+1);
440 
441  double distMin = 1E99;
442 
443  // Some rotation around the boost axis, for more randomness
444  XYZVector theAxis = theBoost.Vect().Unit();
445  double theAngle = random->flatShoot() * 2. * 3.14159265358979323;
446  RawParticle::Rotation axisRotation(theAxis,theAngle);
447  RawParticle::Boost axisBoost(theBoost.x(),theBoost.y(),theBoost.z());
448 
449  // A rotation to bring the particles back to the pion direction
450  XYZVector zAxis(0.,0.,1.);
451  XYZVector orthAxis = (zAxis.Cross(theBoost.Vect())).Unit();
452  double orthAngle = acos(theBoost.Vect().Unit().Z());
453  RawParticle::Rotation orthRotation(orthAxis,orthAngle);
454 
455  // A few checks
456  // double eAfter = 0.;
457 
458  // Loop on the nuclear interaction products
459  for ( unsigned iTrack=firstTrack; iTrack<=lastTrack; ++iTrack ) {
460 
461  unsigned idaugh = iTrack - firstTrack;
462  NUEvent::NUParticle aParticle = aNUEvents[ene]->theNUParticles()[iTrack];
463  // std::cout << "Track " << iTrack
464  // << " id/px/py/pz/mass "
465  // << aParticle.id << " "
466  // << aParticle.px << " "
467  // << aParticle.py << " "
468  // << aParticle.pz << " "
469  // << aParticle.mass << " " << endl;
470 
471  // Add a RawParticle with the proper energy in the c.m frame of
472  // the nuclear interaction
473  double energy = std::sqrt( aParticle.px*aParticle.px
474  + aParticle.py*aParticle.py
475  + aParticle.pz*aParticle.pz
476  + aParticle.mass*aParticle.mass/(ecm*ecm) );
477 
478  RawParticle& aDaughter = _theUpdatedState[idaugh];
479  aDaughter.SetXYZT(aParticle.px*ecm,aParticle.py*ecm,
480  aParticle.pz*ecm,energy*ecm);
481  aDaughter.setID(aParticle.id);
482 
483  // Rotate to the collision axis
484  aDaughter.rotate(orthRotation);
485 
486  // Rotate around the boost axis for more randomness
487  aDaughter.rotate(axisRotation);
488 
489  // Boost it in the lab frame
490  aDaughter.boost(axisBoost);
491 
492  // Store the closest daughter index (for later tracking purposes, so charged particles only)
493  double distance = distanceToPrimary(Particle,aDaughter);
494  // Find the closest daughter, if closer than a given upper limit.
495  if ( distance < distMin && distance < theDistCut ) {
496  distMin = distance;
498  }
499 
500  // eAfter += aDaughter.E();
501 
502  }
503 
504  /*
505  double eBefore = Particle.E();
506  double rapid = Particle.momentum().Eta();
507  if ( eBefore > 0. ) {
508  hAfter->Fill(eAfter/eBefore);
509  hAfter2->Fill(rapid,eAfter/eBefore);
510  hAfter3->Fill(eBefore,eAfter/eBefore);
511  }
512  */
513 
514  // ERROR The way this loops through the events breaks
515  // replay. Which events are retrieved depends on
516  // which previous events were processed.
517 
518  // Increment for next time
519  ++aCurrentInteraction[ene];
520 
521  // Simulate a stopping hadron (low momentum)
522  } else if ( pHadron < 4. &&
523  elastic > 1.- (inelastic4*theInelasticLength)
524  /theTotalInteractionLength ) {
525  // A fake particle with 0 momentum as a daughter!
526  _theUpdatedState.resize(1);
527  _theUpdatedState[0].SetXYZT(0.,0.,0.,0.);
528  _theUpdatedState[0].setID(22);
529  }
530 
531  }
532 
533  }
534 
535  }
536 
537 }
538 
539 double
541  const RawParticle& aDaughter) const {
542 
543  double distance = 2E99;
544 
545  // Compute the distance only for charged primaries
546  if ( fabs(Particle.charge()) > 1E-12 ) {
547 
548  // The secondary must have the same charge
549  double chargeDiff = fabs(aDaughter.charge()-Particle.charge());
550  if ( fabs(chargeDiff) < 1E-12 ) {
551 
552  // Here are two distance definitions * to be tuned *
553  switch ( theDistAlgo ) {
554 
555  case 1:
556  // sin(theta12)
557  distance = (aDaughter.Vect().Unit().Cross(Particle.Vect().Unit())).R();
558  break;
559 
560  case 2:
561  // sin(theta12) * p1/p2
562  distance = (aDaughter.Vect().Cross(Particle.Vect())).R()
563  /aDaughter.Vect().Mag2();
564  break;
565 
566  default:
567  // Should not happen
568  distance = 2E99;
569  break;
570 
571  }
572 
573  }
574 
575  }
576 
577  return distance;
578 
579 }
580 
581 void
583 
584  // Size of buffer
585  ++myOutputBuffer;
586 
587  // Periodically close the current file and open a new one
588  if ( myOutputBuffer/1000*1000 == myOutputBuffer ) {
589  myOutputFile.close();
590  myOutputFile.open ("NuclearInteractionOutputFile.txt");
591  // myOutputFile.seekp(0); // No need to rewind in that case
592  }
593  //
594  unsigned size1 =
595  theCurrentEntry.size()*
596  theCurrentEntry.begin()->size();
597  std::vector<unsigned> theCurrentEntries;
598  theCurrentEntries.resize(size1);
599  size1*=sizeof(unsigned);
600  //
601  unsigned size2 =
602  theCurrentInteraction.size()*
603  theCurrentInteraction.begin()->size();
604  std::vector<unsigned> theCurrentInteractions;
605  theCurrentInteractions.resize(size2);
606  size2 *= sizeof(unsigned);
607 
608  // Save the current entries
609  std::vector< std::vector<unsigned> >::const_iterator aCurrentEntry = theCurrentEntry.begin();
610  std::vector< std::vector<unsigned> >::const_iterator lastCurrentEntry = theCurrentEntry.end();
611  unsigned allEntries=0;
612  for ( ; aCurrentEntry!=lastCurrentEntry; ++aCurrentEntry ) {
613  unsigned size = aCurrentEntry->size();
614  for ( unsigned iene=0; iene<size; ++iene )
615  theCurrentEntries[allEntries++] = (*aCurrentEntry)[iene];
616  }
617 
618  // Save the current interactions
619  std::vector< std::vector<unsigned> >::const_iterator aCurrentInteraction = theCurrentInteraction.begin();
620  std::vector< std::vector<unsigned> >::const_iterator lastCurrentInteraction = theCurrentInteraction.end();
621  unsigned allInteractions=0;
622  for ( ; aCurrentInteraction!=lastCurrentInteraction; ++aCurrentInteraction ) {
623  unsigned size = aCurrentInteraction->size();
624  for ( unsigned iene=0; iene<size; ++iene )
625  theCurrentInteractions[allInteractions++] = (*aCurrentInteraction)[iene];
626  }
627  //
628  myOutputFile.write((const char*)(&theCurrentEntries.front()), size1);
629  myOutputFile.write((const char*)(&theCurrentInteractions.front()), size2);
630  myOutputFile.flush();
631 
632 }
633 
634 bool
636 
637  std::ifstream myInputFile;
638  struct stat results;
639  //
640  unsigned size1 =
641  theCurrentEntry.size()*
642  theCurrentEntry.begin()->size();
643  std::vector<unsigned> theCurrentEntries;
644  theCurrentEntries.resize(size1);
645  size1*=sizeof(unsigned);
646  //
647  unsigned size2 =
648  theCurrentInteraction.size()*
649  theCurrentInteraction.begin()->size();
650  std::vector<unsigned> theCurrentInteractions;
651  theCurrentInteractions.resize(size2);
652  size2 *= sizeof(unsigned);
653  //
654  unsigned size = 0;
655 
656 
657  // Open the file (if any)
658  myInputFile.open (inputFile.c_str());
659  if ( myInputFile.is_open() ) {
660 
661  // Get the size of the file
662  if ( stat(inputFile.c_str(), &results) == 0 ) size = results.st_size;
663  else return false; // Something is wrong with that file !
664 
665  // Position the pointer just before the last record
666  myInputFile.seekg(size-size1-size2);
667  myInputFile.read((char*)(&theCurrentEntries.front()),size1);
668  myInputFile.read((char*)(&theCurrentInteractions.front()),size2);
669  myInputFile.close();
670 
671  // Read the current entries
672  std::vector< std::vector<unsigned> >::iterator aCurrentEntry = theCurrentEntry.begin();
673  std::vector< std::vector<unsigned> >::iterator lastCurrentEntry = theCurrentEntry.end();
674  unsigned allEntries=0;
675  for ( ; aCurrentEntry!=lastCurrentEntry; ++aCurrentEntry ) {
676  unsigned size = aCurrentEntry->size();
677  for ( unsigned iene=0; iene<size; ++iene )
678  (*aCurrentEntry)[iene] = theCurrentEntries[allEntries++];
679  }
680 
681  // Read the current interactions
682  std::vector< std::vector<unsigned> >::iterator aCurrentInteraction = theCurrentInteraction.begin();
683  std::vector< std::vector<unsigned> >::iterator lastCurrentInteraction = theCurrentInteraction.end();
684  unsigned allInteractions=0;
685  for ( ; aCurrentInteraction!=lastCurrentInteraction; ++aCurrentInteraction ) {
686  unsigned size = aCurrentInteraction->size();
687  for ( unsigned iene=0; iene<size; ++iene )
688  (*aCurrentInteraction)[iene] = theCurrentInteractions[allInteractions++];
689  }
690 
691  return true;
692  }
693 
694  return false;
695 
696 }
697 
698 unsigned
700 
701  unsigned myIndex=0;
702  while ( thePid != thePionID[myIndex] ) ++myIndex;
703  // std::cout << "pid/index = " << thePid << " " << myIndex << std::endl;
704  return myIndex;
705 
706 }
void boost(double bx, double by, double bz)
Definition: RawParticle.cc:183
double flatShoot(double xmin=0.0, double xmax=1.0) const
dictionary ratios
Definition: plotFactory.py:68
static const double slope[3]
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
std::vector< std::vector< TFile * > > theFiles
Geom::Theta< T > theta() const
XYZVector orthogonal(const XYZVector &) const
A vector orthogonal to another one (because it&#39;s not in XYZTLorentzVector)
std::vector< std::vector< unsigned > > theNumberOfInteractions
std::vector< std::vector< TTree * > > theTrees
int pid() const
get the HEP particle ID number
Definition: RawParticle.h:265
TRandom random
Definition: MVATrainer.cc:138
double mass() const
get the MEASURED mass
Definition: RawParticle.h:283
ROOT::Math::Boost Boost
Definition: RawParticle.h:40
std::vector< std::vector< unsigned > > theCurrentEntry
std::vector< std::vector< std::string > > theFileNames
Definition: NUEvent.h:6
void compute(ParticlePropagator &Particle, RandomEngineAndDistribution const *)
Generate a nuclear interaction according to the probability that it happens.
math::XYZVector XYZVector
std::vector< std::vector< unsigned > > theNumberOfEntries
T sqrt(T t)
Definition: SSEVec.h:18
void setID(const int id)
Definition: RawParticle.cc:102
tuple results
Definition: mps_update.py:44
std::vector< std::vector< double > > thePionCM
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
void rotate(double rphi, const XYZVector &raxis)
Definition: RawParticle.cc:155
bool read(std::string inputFile)
Read former nuclear interaction (from previous run)
double charge() const
get the MEASURED charge
Definition: RawParticle.h:282
unsigned index(int thePid)
Return a hashed index for a given pid.
std::vector< std::string > thePionNA
std::vector< std::vector< TBranch * > > theBranches
std::vector< std::vector< unsigned > > theCurrentInteraction
std::vector< std::vector< double > > theRatios
~NuclearInteractionSimulator()
Default Destructor.
NuclearInteractionSimulator(std::vector< double > &hadronEnergies, std::vector< int > &hadronTypes, std::vector< std::string > &hadronNames, std::vector< double > &hadronMasses, std::vector< double > &hadronPMin, double pionEnergy, std::vector< double > &lengthRatio, std::vector< std::vector< double > > &ratios, std::map< int, int > &idMap, std::string inputFile, unsigned int distAlgo, double distCut)
Constructor.
std::vector< RawParticle > _theUpdatedState
tuple filename
Definition: lut2db_cfg.py:20
ROOT::Math::AxisAngle Rotation
Definition: RawParticle.h:35
tuple zAxis
Definition: MetAnalyzer.py:56
tuple cout
Definition: gather_cfg.py:145
std::vector< std::vector< NUEvent * > > theNUEvents
volatile std::atomic< bool > shutdown_flag false
std::string fullPath() const
Definition: FileInPath.cc:184
tuple size
Write out results.
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
void save()
Save current nuclear interaction (for later use)
double distanceToPrimary(const RawParticle &Particle, const RawParticle &aDaughter) const
Compute distance between secondary and primary.
math::XYZTLorentzVector XYZTLorentzVector
Definition: RawParticle.h:15