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