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