CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Member Functions | Private Attributes
NuclearInteractionSimulator Class Reference

#include <NuclearInteractionSimulator.h>

Inheritance diagram for NuclearInteractionSimulator:
MaterialEffectsSimulator

Public Member Functions

 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. More...
 
bool read (std::string inputFile)
 Read former nuclear interaction (from previous run) More...
 
void save ()
 Save current nuclear interaction (for later use) More...
 
 ~NuclearInteractionSimulator ()
 Default Destructor. More...
 
- Public Member Functions inherited from MaterialEffectsSimulator
RHEP_const_iter beginDaughters () const
 Returns const iterator to the beginning of the daughters list. More...
 
int closestDaughterId ()
 The id of the closest charged daughter (filled for nuclear interactions only) More...
 
double eMass () const
 Electron mass in GeV/c2. More...
 
RHEP_const_iter endDaughters () const
 Returns const iterator to the end of the daughters list. More...
 
double excitE () const
 Mean excitation energy (in GeV) More...
 
 MaterialEffectsSimulator (double A=28.0855, double Z=14.0000, double density=2.329, double radLen=9.360)
 
unsigned nDaughters () const
 Returns the number of daughters. More...
 
XYZVector orthogonal (const XYZVector &) const
 A vector orthogonal to another one (because it's not in XYZTLorentzVector) More...
 
double radLenIncm () const
 One radiation length in cm. More...
 
double rho () const
 Density in g/cm3. More...
 
void setNormalVector (const GlobalVector &normal)
 Sets the vector normal to the surface traversed. More...
 
double theA () const
 A. More...
 
double theZ () const
 Z. More...
 
void updateState (ParticlePropagator &myTrack, double radlen, RandomEngineAndDistribution const *)
 Compute the material effect (calls the sub class) More...
 
virtual ~MaterialEffectsSimulator ()
 

Private Member Functions

void compute (ParticlePropagator &Particle, RandomEngineAndDistribution const *)
 Generate a nuclear interaction according to the probability that it happens. More...
 
double distanceToPrimary (const RawParticle &Particle, const RawParticle &aDaughter) const
 Compute distance between secondary and primary. More...
 
unsigned index (int thePid)
 Return a hashed index for a given pid. More...
 

Private Attributes

bool currentValuesWereSet
 
unsigned ien4
 
unsigned myOutputBuffer
 
std::ofstream myOutputFile
 
std::vector< std::vector
< TBranch * > > 
theBranches
 
std::vector< std::vector
< unsigned > > 
theCurrentEntry
 
std::vector< std::vector
< unsigned > > 
theCurrentInteraction
 
unsigned theDistAlgo
 
double theDistCut
 
std::vector< std::vector
< std::string > > 
theFileNames
 
std::vector< std::vector
< TFile * > > 
theFiles
 
std::map< int, int > theIDMap
 
std::vector< double > theLengthRatio
 
std::vector< std::vector
< NUEvent * > > 
theNUEvents
 
std::vector< std::vector
< unsigned > > 
theNumberOfEntries
 
std::vector< std::vector
< unsigned > > 
theNumberOfInteractions
 
std::vector< std::vector
< double > > 
thePionCM
 
std::vector< double > thePionEN
 
double thePionEnergy
 
std::vector< int > thePionID
 
std::vector< double > thePionMA
 
std::vector< std::string > thePionNA
 
std::vector< double > thePionPMin
 
std::vector< std::vector
< double > > 
theRatios
 
std::vector< std::vector
< TTree * > > 
theTrees
 

Additional Inherited Members

- Public Types inherited from MaterialEffectsSimulator
typedef std::vector
< RawParticle >
::const_iterator 
RHEP_const_iter
 
- Protected Attributes inherited from MaterialEffectsSimulator
std::vector< RawParticle_theUpdatedState
 
double A
 
double density
 
double radLen
 
double radLengths
 
int theClosestChargedDaughterId
 
GlobalVector theNormalVector
 
double Z
 

Detailed Description

Definition at line 33 of file NuclearInteractionSimulator.h.

Constructor & Destructor Documentation

NuclearInteractionSimulator::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.

Definition at line 24 of file NuclearInteractionSimulator.cc.

References gather_cfg::cout, currentValuesWereSet, Exception, lut2db_cfg::filename, contentValuesFiles::fullPath, edm::FileInPath::fullPath(), ien4, myOutputBuffer, myOutputFile, read(), mathSSE::sqrt(), AlCaHLTBitMon_QueryRunRegistry::string, theBranches, theCurrentEntry, theCurrentInteraction, theFileNames, theFiles, theNUEvents, theNumberOfEntries, theNumberOfInteractions, thePionCM, thePionEN, thePionMA, thePionNA, and theTrees.

37  :
39  thePionEN(hadronEnergies),
40  thePionID(hadronTypes),
41  thePionNA(hadronNames),
42  thePionMA(hadronMasses),
43  thePionPMin(hadronPMin),
44  thePionEnergy(pionEnergy),
45  theLengthRatio(lengthRatio),
47  theIDMap(idMap),
48  theDistAlgo(distAlgo),
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)
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 }
dictionary ratios
Definition: plotFactory.py:68
std::vector< std::vector< TFile * > > theFiles
std::vector< std::vector< unsigned > > theNumberOfInteractions
std::vector< std::vector< TTree * > > theTrees
std::vector< std::vector< unsigned > > theCurrentEntry
std::vector< std::vector< std::string > > theFileNames
Definition: NUEvent.h:6
std::vector< std::vector< unsigned > > theNumberOfEntries
T sqrt(T t)
Definition: SSEVec.h:48
MaterialEffectsSimulator(double A=28.0855, double Z=14.0000, double density=2.329, double radLen=9.360)
std::vector< std::vector< double > > thePionCM
bool read(std::string inputFile)
Read former nuclear interaction (from previous run)
std::vector< std::string > thePionNA
std::vector< std::vector< TBranch * > > theBranches
std::vector< std::vector< unsigned > > theCurrentInteraction
std::vector< std::vector< double > > theRatios
tuple filename
Definition: lut2db_cfg.py:20
tuple cout
Definition: gather_cfg.py:121
std::vector< std::vector< NUEvent * > > theNUEvents
math::XYZTLorentzVector XYZTLorentzVector
Definition: RawParticle.h:15
NuclearInteractionSimulator::~NuclearInteractionSimulator ( )

Default Destructor.

Definition at line 179 of file NuclearInteractionSimulator.cc.

References compare_using_db::ifile, myOutputFile, and theFiles.

179  {
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  // std::cout << "Closing file " << iene << " with name " << theFileNames[ifile][iene] << std::endl;
187  theFiles[ifile][iene]->Close();
188  }
189  }
190 
191  // Close the output file
192  myOutputFile.close();
193 
194  // dbe->save("test.root");
195 
196 }
std::vector< std::vector< TFile * > > theFiles

Member Function Documentation

void NuclearInteractionSimulator::compute ( ParticlePropagator Particle,
RandomEngineAndDistribution const *  random 
)
privatevirtual

Generate a nuclear interaction according to the probability that it happens.

Implements MaterialEffectsSimulator.

Definition at line 198 of file NuclearInteractionSimulator.cc.

References MaterialEffectsSimulator::_theUpdatedState, funct::abs(), RawParticle::boost(), currentValuesWereSet, HLT_25ns14e33_v1_cff::distance, distanceToPrimary(), relval_parameters_module::energy, create_public_lumi_plots::exp, NUEvent::NUInteraction::first, RandomEngineAndDistribution::flatShoot(), NUEvent::NUParticle::id, ien4, index(), NUEvent::NUInteraction::last, cmsBatch::log, NUEvent::NUParticle::mass, RawParticle::mass(), MaterialEffectsSimulator::orthogonal(), phi, RawParticle::pid(), HLT_25ns14e33_v1_cff::pMin, funct::pow(), NUEvent::NUParticle::px, NUEvent::NUParticle::py, NUEvent::NUParticle::pz, MaterialEffectsSimulator::radLengths, RawParticle::rotate(), RawParticle::setID(), funct::sin(), slope, mathSSE::sqrt(), MaterialEffectsSimulator::theA(), MaterialEffectsSimulator::theClosestChargedDaughterId, theCurrentEntry, theCurrentInteraction, theDistCut, theIDMap, theLengthRatio, theNUEvents, theNumberOfEntries, theNumberOfInteractions, thePionCM, thePionEN, thePionEnergy, thePionNA, thePionPMin, theRatios, theta(), theTrees, and MetAnalyzer::zAxis.

199 {
200  if(!currentValuesWereSet) {
201  currentValuesWereSet = true;
202  for ( unsigned iname=0; iname<thePionNA.size(); ++iname ) {
203  for ( unsigned iene=0; iene<thePionEN.size(); ++iene ) {
204  theCurrentEntry[iname][iene] = (unsigned) (theNumberOfEntries[iname][iene] * random->flatShoot());
205 
206  theTrees[iname][iene]->GetEntry(theCurrentEntry[iname][iene]);
207  unsigned NInteractions = theNUEvents[iname][iene]->nInteractions();
208  theNumberOfInteractions[iname][iene] = NInteractions;
209 
210  theCurrentInteraction[iname][iene] = (unsigned) (theNumberOfInteractions[iname][iene] * random->flatShoot());
211  }
212  }
213  }
214 
215  // Read a Nuclear Interaction in a random manner
216 
217  double pHadron = std::sqrt(Particle.Vect().Mag2());
218  // htot->Fill(pHadron);
219 
220  // The hadron has enough momentum to create some relevant final state
221  if ( pHadron > thePionEnergy ) {
222 
223  // The particle type
224  std::map<int,int>::const_iterator thePit = theIDMap.find(Particle.pid());
225 
226  int thePid = thePit != theIDMap.end() ? thePit->second : Particle.pid();
227 
228  // Is this particle type foreseen?
229  unsigned fPid = abs(thePid);
230  if ( fPid != 211 && fPid != 130 && fPid != 321 && fPid != 2112 && fPid != 2212 ) {
231  return;
232  //std::cout << "Unknown particle type = " << thePid << std::endl;
233  //thePid = 211;
234  }
235 
236  // The inelastic interaction length at p(pion) = 5 GeV/c
237  unsigned thePidIndex = index(thePid);
238  double theInelasticLength = radLengths * theLengthRatio[thePidIndex];
239 
240  // The elastic interaction length
241  // The baroque parameterization is a fit to Fig. 40.13 of the PDG
242  double ee = pHadron > 0.6 ?
243  exp(-std::sqrt((pHadron-0.6)/1.122)) : exp(std::sqrt((0.6-pHadron)/1.122));
244  double theElasticLength = ( 0.8753 * ee + 0.15 )
245  // double theElasticLength = ( 0.15 + 0.195 / log(pHadron/0.4) )
246  // double theElasticLength = ( 0.15 + 0.305 / log(pHadron/0.35) )
247  * theInelasticLength;
248 
249  // The total interaction length
250  double theTotalInteractionLength = theInelasticLength + theElasticLength;
251 
252  // Probability to interact is dl/L0 (maximum for 4 GeV pion)
253  double aNuclInteraction = -std::log(random->flatShoot());
254  if ( aNuclInteraction < theTotalInteractionLength ) {
255 
256  // The elastic part
257  double elastic = random->flatShoot();
258  if ( elastic < theElasticLength/theTotalInteractionLength ) {
259 
260  // helas->Fill(pHadron);
261 
262  // Characteristic scattering angle for the elastic part
263  double theta0 = std::sqrt(3.)/ std::pow(theA(),1./3.) * Particle.mass()/pHadron;
264 
265  // Draw an angle with theta/theta0*exp[(-theta/2theta0)**2] shape
266  double theta = theta0 * std::sqrt(-2.*std::log(random->flatShoot()));
267  double phi = 2. * 3.14159265358979323 * random->flatShoot();
268 
269  // Rotate the particle accordingly
270  RawParticle::Rotation rotation1(orthogonal(Particle.Vect()),theta);
271  RawParticle::Rotation rotation2(Particle.Vect(),phi);
272  Particle.rotate(rotation1);
273  Particle.rotate(rotation2);
274 
275  // Distance
276  double distance = std::sin(theta);
277 
278  // Create a daughter if the kink is large engough
279  if ( distance > theDistCut ) {
280  _theUpdatedState.resize(1);
281  _theUpdatedState[0].SetXYZT(Particle.Px(), Particle.Py(),
282  Particle.Pz(), Particle.E());
283  _theUpdatedState[0].setID(Particle.pid());
284  }
285 
286  // hscatter->Fill(myTheta);
287  // hscatter2->Fill(pHadron,myTheta);
288 
289  }
290 
291  // The inelastic part
292  else {
293 
294  // Avoid multiple map access
295  const std::vector<double>& aPionCM = thePionCM[thePidIndex];
296  const std::vector<double>& aRatios = theRatios[thePidIndex];
297  // Find the file with the closest c.m energy
298  // The target nucleon
299  XYZTLorentzVector Proton(0.,0.,0.,0.939);
300  // The current particle
301  const XYZTLorentzVector& Hadron = (const XYZTLorentzVector&)Particle;
302  // The smallest momentum for inelastic interactions
303  double pMin = thePionPMin[thePidIndex];
304  // The correspong smallest four vector
305  XYZTLorentzVector Hadron0(0.,0.,pMin,std::sqrt(pMin*pMin+Particle.M2()));
306 
307  // The current centre-of-mass energy
308  double ecm = (Proton+Hadron).M();
309  //std::cout << "Proton = " << Proton << std::endl;
310  //std::cout << "Hadron = " << Hadron << std::endl;
311  //std::cout << "ecm = " << ecm << std::endl;
312  // Get the files of interest (closest c.m. energies)
313  unsigned ene1=0;
314  unsigned ene2=0;
315  // The smallest centre-of-mass energy
316  // double ecm1=1.63;
317  double ecm1= (Proton+Hadron0).M();
318  //std::cout << "ecm1 = " << ecm1 << std::endl;
319  //std::cout << "ecm[0] = " << aPionCM[0] << std::endl;
320  //std::cout << "ecm[11] = " << aPionCM [ aPionCM.size()-1 ] << std::endl;
321  double ecm2=aPionCM[0];
322  double ratio1=0.;
323  double ratio2=aRatios[0];
324  if ( ecm > aPionCM[0] && ecm < aPionCM [ aPionCM.size()-1 ] ) {
325  for ( unsigned ene=1;
326  ene < aPionCM.size() && ecm > aPionCM[ene-1];
327  ++ene ) {
328  if ( ecm<aPionCM[ene] ) {
329  ene2 = ene;
330  ene1 = ene2-1;
331  ecm1 = aPionCM[ene1];
332  ecm2 = aPionCM[ene2];
333  ratio1 = aRatios[ene1];
334  ratio2 = aRatios[ene2];
335  }
336  }
337  } else if ( ecm > aPionCM[ aPionCM.size()-1 ] ) {
338  ene1 = aPionCM.size()-1;
339  ene2 = aPionCM.size()-2;
340  ecm1 = aPionCM[ene1];
341  ecm2 = aPionCM[ene2];
342  ratio1 = aRatios[ene2];
343  ratio2 = aRatios[ene2];
344  }
345 
346 
347  // The inelastic part of the cross section depends cm energy
348  double slope = (std::log10(ecm )-std::log10(ecm1))
349  / (std::log10(ecm2)-std::log10(ecm1));
350  double inelastic = ratio1 + (ratio2-ratio1) * slope;
351  double inelastic4 = pHadron < 4. ? aRatios[ien4] : 1.;
352 
353  //std::cout << "Inelastic = " << ratio1 << " " << ratio2 << " " << inelastic << std::endl;
354  // std::cout << "Energy/Inelastic : "
355  // << Hadron.e() << " " << inelastic << std::endl;
356 
357  // Simulate an inelastic interaction
358  if ( elastic > 1.- (inelastic*theInelasticLength)
359  /theTotalInteractionLength ) {
360 
361  // Avoid mutliple map access
362  std::vector<unsigned>& aCurrentInteraction = theCurrentInteraction[thePidIndex];
363  std::vector<unsigned>& aNumberOfInteractions = theNumberOfInteractions[thePidIndex];
364  std::vector<NUEvent*>& aNUEvents = theNUEvents[thePidIndex];
365  // hinel->Fill(pHadron);
366  // std::cout << "INELASTIC INTERACTION ! "
367  // << pHadron << " " << theInelasticLength << " "
368  // << inelastic * theInelasticLength << std::endl;
369  // Choice of the file to read according the the log10(ecm) distance
370  // and protection against low momentum proton and neutron that never interacts
371  // (i.e., empty files)
372  unsigned ene;
373  if ( random->flatShoot() < slope || aNumberOfInteractions[ene1] == 0 )
374  ene = ene2;
375  else
376  ene = ene1;
377 
378  //std::cout << "Ecm1/2 = " << ecm1 << " " << ecm2 << std::endl;
379  //std::cout << "Ratio1/2 = " << ratio1 << " " << ratio2 << std::endl;
380  //std::cout << "Ene = " << ene << " slope = " << slope << std::endl;
381 
382  //std::cout << "Pion energy = " << Hadron.E()
383  // << "File chosen " << theFileNames[thePidIndex][ene]
384  // << std::endl;
385 
386  // The boost characteristics
387  XYZTLorentzVector theBoost = Proton + Hadron;
388  theBoost /= theBoost.e();
389 
390  // std::cout << "File chosen : " << thePid << "/" << ene
391  // << " Current interaction = " << aCurrentInteraction[ene]
392  // << " Total interactions = " << aNumberOfInteractions[ene]
393  // << std::endl;
394  // theFiles[thePidIndex][ene]->cd();
395  // gDirectory->ls();
396 
397  // Check we are not either at the end of an interaction bunch
398  // or at the end of a file
399  if ( aCurrentInteraction[ene] == aNumberOfInteractions[ene] ) {
400  // std::cout << "End of interaction bunch ! ";
401  std::vector<unsigned>& aCurrentEntry = theCurrentEntry[thePidIndex];
402  std::vector<unsigned>& aNumberOfEntries = theNumberOfEntries[thePidIndex];
403  std::vector<TTree*>& aTrees = theTrees[thePidIndex];
404  ++aCurrentEntry[ene];
405  // std::cerr << "Read the next entry "
406  // << aCurrentEntry[ene] << std::endl;
407  aCurrentInteraction[ene] = 0;
408  if ( aCurrentEntry[ene] == aNumberOfEntries[ene] ) {
409  aCurrentEntry[ene] = 0;
410  // std::cout << "End of file - Rewind! " << std::endl;
411  }
412  unsigned myEntry = aCurrentEntry[ene];
413  // std::cout << "The new entry " << myEntry
414  // << " is read ... in TTree " << aTrees[ene] << " ";
415  aTrees[ene]->GetEntry(myEntry);
416  // std::cout << "The number of interactions in the new entry is ... ";
417  aNumberOfInteractions[ene] = aNUEvents[ene]->nInteractions();
418  // std::cout << aNumberOfInteractions[ene] << std::endl;
419  }
420 
421  // Read the interaction
422  NUEvent::NUInteraction anInteraction
423  = aNUEvents[ene]->theNUInteractions()[aCurrentInteraction[ene]];
424 
425  unsigned firstTrack = anInteraction.first;
426  unsigned lastTrack = anInteraction.last;
427  // std::cout << "First and last tracks are " << firstTrack << " " << lastTrack << std::endl;
428 
429  _theUpdatedState.resize(lastTrack-firstTrack+1);
430 
431  double distMin = 1E99;
432 
433  // Some rotation around the boost axis, for more randomness
434  XYZVector theAxis = theBoost.Vect().Unit();
435  double theAngle = random->flatShoot() * 2. * 3.14159265358979323;
436  RawParticle::Rotation axisRotation(theAxis,theAngle);
437  RawParticle::Boost axisBoost(theBoost.x(),theBoost.y(),theBoost.z());
438 
439  // A rotation to bring the particles back to the pion direction
440  XYZVector zAxis(0.,0.,1.);
441  XYZVector orthAxis = (zAxis.Cross(theBoost.Vect())).Unit();
442  double orthAngle = acos(theBoost.Vect().Unit().Z());
443  RawParticle::Rotation orthRotation(orthAxis,orthAngle);
444 
445  // A few checks
446  // double eAfter = 0.;
447 
448  // Loop on the nuclear interaction products
449  for ( unsigned iTrack=firstTrack; iTrack<=lastTrack; ++iTrack ) {
450 
451  unsigned idaugh = iTrack - firstTrack;
452  NUEvent::NUParticle aParticle = aNUEvents[ene]->theNUParticles()[iTrack];
453  // std::cout << "Track " << iTrack
454  // << " id/px/py/pz/mass "
455  // << aParticle.id << " "
456  // << aParticle.px << " "
457  // << aParticle.py << " "
458  // << aParticle.pz << " "
459  // << aParticle.mass << " " << endl;
460 
461  // Add a RawParticle with the proper energy in the c.m frame of
462  // the nuclear interaction
463  double energy = std::sqrt( aParticle.px*aParticle.px
464  + aParticle.py*aParticle.py
465  + aParticle.pz*aParticle.pz
466  + aParticle.mass*aParticle.mass/(ecm*ecm) );
467 
468  RawParticle& aDaughter = _theUpdatedState[idaugh];
469  aDaughter.SetXYZT(aParticle.px*ecm,aParticle.py*ecm,
470  aParticle.pz*ecm,energy*ecm);
471  aDaughter.setID(aParticle.id);
472 
473  // Rotate to the collision axis
474  aDaughter.rotate(orthRotation);
475 
476  // Rotate around the boost axis for more randomness
477  aDaughter.rotate(axisRotation);
478 
479  // Boost it in the lab frame
480  aDaughter.boost(axisBoost);
481 
482  // Store the closest daughter index (for later tracking purposes, so charged particles only)
483  double distance = distanceToPrimary(Particle,aDaughter);
484  // Find the closest daughter, if closer than a given upper limit.
485  if ( distance < distMin && distance < theDistCut ) {
486  distMin = distance;
488  }
489 
490  // eAfter += aDaughter.E();
491 
492  }
493 
494  /*
495  double eBefore = Particle.E();
496  double rapid = Particle.momentum().Eta();
497  if ( eBefore > 0. ) {
498  hAfter->Fill(eAfter/eBefore);
499  hAfter2->Fill(rapid,eAfter/eBefore);
500  hAfter3->Fill(eBefore,eAfter/eBefore);
501  }
502  */
503 
504  // ERROR The way this loops through the events breaks
505  // replay. Which events are retrieved depends on
506  // which previous events were processed.
507 
508  // Increment for next time
509  ++aCurrentInteraction[ene];
510 
511  // Simulate a stopping hadron (low momentum)
512  } else if ( pHadron < 4. &&
513  elastic > 1.- (inelastic4*theInelasticLength)
514  /theTotalInteractionLength ) {
515  // A fake particle with 0 momentum as a daughter!
516  _theUpdatedState.resize(1);
517  _theUpdatedState[0].SetXYZT(0.,0.,0.,0.);
518  _theUpdatedState[0].setID(22);
519  }
520 
521  }
522 
523  }
524 
525  }
526 
527 }
void boost(double bx, double by, double bz)
Definition: RawParticle.cc:183
static const double slope[3]
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
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
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
unsigned index(int thePid)
Return a hashed index for a given pid.
std::vector< std::string > thePionNA
std::vector< std::vector< unsigned > > theCurrentInteraction
std::vector< std::vector< double > > theRatios
std::vector< RawParticle > _theUpdatedState
ROOT::Math::AxisAngle Rotation
Definition: RawParticle.h:35
tuple zAxis
Definition: MetAnalyzer.py:56
std::vector< std::vector< NUEvent * > > theNUEvents
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
double distanceToPrimary(const RawParticle &Particle, const RawParticle &aDaughter) const
Compute distance between secondary and primary.
math::XYZTLorentzVector XYZTLorentzVector
Definition: RawParticle.h:15
tuple log
Definition: cmsBatch.py:341
double NuclearInteractionSimulator::distanceToPrimary ( const RawParticle Particle,
const RawParticle aDaughter 
) const
private

Compute distance between secondary and primary.

Definition at line 530 of file NuclearInteractionSimulator.cc.

References RawParticle::charge(), HLT_25ns14e33_v1_cff::distance, and theDistAlgo.

Referenced by compute().

531  {
532 
533  double distance = 2E99;
534 
535  // Compute the distance only for charged primaries
536  if ( fabs(Particle.charge()) > 1E-12 ) {
537 
538  // The secondary must have the same charge
539  double chargeDiff = fabs(aDaughter.charge()-Particle.charge());
540  if ( fabs(chargeDiff) < 1E-12 ) {
541 
542  // Here are two distance definitions * to be tuned *
543  switch ( theDistAlgo ) {
544 
545  case 1:
546  // sin(theta12)
547  distance = (aDaughter.Vect().Unit().Cross(Particle.Vect().Unit())).R();
548  break;
549 
550  case 2:
551  // sin(theta12) * p1/p2
552  distance = (aDaughter.Vect().Cross(Particle.Vect())).R()
553  /aDaughter.Vect().Mag2();
554  break;
555 
556  default:
557  // Should not happen
558  distance = 2E99;
559  break;
560 
561  }
562 
563  }
564 
565  }
566 
567  return distance;
568 
569 }
double charge() const
get the MEASURED charge
Definition: RawParticle.h:282
unsigned NuclearInteractionSimulator::index ( int  thePid)
private

Return a hashed index for a given pid.

Definition at line 689 of file NuclearInteractionSimulator.cc.

References thePionID.

Referenced by BeautifulSoup.PageElement::_invert(), and compute().

689  {
690 
691  unsigned myIndex=0;
692  while ( thePid != thePionID[myIndex] ) ++myIndex;
693  // std::cout << "pid/index = " << thePid << " " << myIndex << std::endl;
694  return myIndex;
695 
696 }
bool NuclearInteractionSimulator::read ( std::string  inputFile)

Read former nuclear interaction (from previous run)

Definition at line 625 of file NuclearInteractionSimulator.cc.

References python.entryComment::results, findQualityFiles::size, theCurrentEntry, and theCurrentInteraction.

Referenced by edmIntegrityCheck.PublishToFileSystem::get(), Vispa.Plugins.EdmBrowser.EdmDataAccessor.EdmDataAccessor::goto(), NuclearInteractionSimulator(), and Vispa.Plugins.EdmBrowser.EdmDataAccessor.EdmDataAccessor::setFilterBranches().

625  {
626 
627  std::ifstream myInputFile;
628  struct stat results;
629  //
630  unsigned size1 =
631  theCurrentEntry.size()*
632  theCurrentEntry.begin()->size();
633  std::vector<unsigned> theCurrentEntries;
634  theCurrentEntries.resize(size1);
635  size1*=sizeof(unsigned);
636  //
637  unsigned size2 =
638  theCurrentInteraction.size()*
639  theCurrentInteraction.begin()->size();
640  std::vector<unsigned> theCurrentInteractions;
641  theCurrentInteractions.resize(size2);
642  size2 *= sizeof(unsigned);
643  //
644  unsigned size = 0;
645 
646 
647  // Open the file (if any)
648  myInputFile.open (inputFile.c_str());
649  if ( myInputFile.is_open() ) {
650 
651  // Get the size of the file
652  if ( stat(inputFile.c_str(), &results) == 0 ) size = results.st_size;
653  else return false; // Something is wrong with that file !
654 
655  // Position the pointer just before the last record
656  myInputFile.seekg(size-size1-size2);
657  myInputFile.read((char*)(&theCurrentEntries.front()),size1);
658  myInputFile.read((char*)(&theCurrentInteractions.front()),size2);
659  myInputFile.close();
660 
661  // Read the current entries
662  std::vector< std::vector<unsigned> >::iterator aCurrentEntry = theCurrentEntry.begin();
663  std::vector< std::vector<unsigned> >::iterator lastCurrentEntry = theCurrentEntry.end();
664  unsigned allEntries=0;
665  for ( ; aCurrentEntry!=lastCurrentEntry; ++aCurrentEntry ) {
666  unsigned size = aCurrentEntry->size();
667  for ( unsigned iene=0; iene<size; ++iene )
668  (*aCurrentEntry)[iene] = theCurrentEntries[allEntries++];
669  }
670 
671  // Read the current interactions
672  std::vector< std::vector<unsigned> >::iterator aCurrentInteraction = theCurrentInteraction.begin();
673  std::vector< std::vector<unsigned> >::iterator lastCurrentInteraction = theCurrentInteraction.end();
674  unsigned allInteractions=0;
675  for ( ; aCurrentInteraction!=lastCurrentInteraction; ++aCurrentInteraction ) {
676  unsigned size = aCurrentInteraction->size();
677  for ( unsigned iene=0; iene<size; ++iene )
678  (*aCurrentInteraction)[iene] = theCurrentInteractions[allInteractions++];
679  }
680 
681  return true;
682  }
683 
684  return false;
685 
686 }
std::vector< std::vector< unsigned > > theCurrentEntry
std::vector< std::vector< unsigned > > theCurrentInteraction
tuple size
Write out results.
void NuclearInteractionSimulator::save ( )
virtual

Save current nuclear interaction (for later use)

Reimplemented from MaterialEffectsSimulator.

Definition at line 572 of file NuclearInteractionSimulator.cc.

References myOutputBuffer, myOutputFile, findQualityFiles::size, theCurrentEntry, and theCurrentInteraction.

Referenced by Vispa.Main.TabController.TabController::allowClose(), and Vispa.Main.TabController.TabController::checkModificationTimestamp().

572  {
573 
574  // Size of buffer
575  ++myOutputBuffer;
576 
577  // Periodically close the current file and open a new one
578  if ( myOutputBuffer/1000*1000 == myOutputBuffer ) {
579  myOutputFile.close();
580  myOutputFile.open ("NuclearInteractionOutputFile.txt");
581  // myOutputFile.seekp(0); // No need to rewind in that case
582  }
583  //
584  unsigned size1 =
585  theCurrentEntry.size()*
586  theCurrentEntry.begin()->size();
587  std::vector<unsigned> theCurrentEntries;
588  theCurrentEntries.resize(size1);
589  size1*=sizeof(unsigned);
590  //
591  unsigned size2 =
592  theCurrentInteraction.size()*
593  theCurrentInteraction.begin()->size();
594  std::vector<unsigned> theCurrentInteractions;
595  theCurrentInteractions.resize(size2);
596  size2 *= sizeof(unsigned);
597 
598  // Save the current entries
599  std::vector< std::vector<unsigned> >::const_iterator aCurrentEntry = theCurrentEntry.begin();
600  std::vector< std::vector<unsigned> >::const_iterator lastCurrentEntry = theCurrentEntry.end();
601  unsigned allEntries=0;
602  for ( ; aCurrentEntry!=lastCurrentEntry; ++aCurrentEntry ) {
603  unsigned size = aCurrentEntry->size();
604  for ( unsigned iene=0; iene<size; ++iene )
605  theCurrentEntries[allEntries++] = (*aCurrentEntry)[iene];
606  }
607 
608  // Save the current interactions
609  std::vector< std::vector<unsigned> >::const_iterator aCurrentInteraction = theCurrentInteraction.begin();
610  std::vector< std::vector<unsigned> >::const_iterator lastCurrentInteraction = theCurrentInteraction.end();
611  unsigned allInteractions=0;
612  for ( ; aCurrentInteraction!=lastCurrentInteraction; ++aCurrentInteraction ) {
613  unsigned size = aCurrentInteraction->size();
614  for ( unsigned iene=0; iene<size; ++iene )
615  theCurrentInteractions[allInteractions++] = (*aCurrentInteraction)[iene];
616  }
617  //
618  myOutputFile.write((const char*)(&theCurrentEntries.front()), size1);
619  myOutputFile.write((const char*)(&theCurrentInteractions.front()), size2);
620  myOutputFile.flush();
621 
622 }
std::vector< std::vector< unsigned > > theCurrentEntry
std::vector< std::vector< unsigned > > theCurrentInteraction
tuple size
Write out results.

Member Data Documentation

bool NuclearInteractionSimulator::currentValuesWereSet
private

Definition at line 100 of file NuclearInteractionSimulator.h.

Referenced by compute(), and NuclearInteractionSimulator().

unsigned NuclearInteractionSimulator::ien4
private

Definition at line 93 of file NuclearInteractionSimulator.h.

Referenced by compute(), and NuclearInteractionSimulator().

unsigned NuclearInteractionSimulator::myOutputBuffer
private

Definition at line 98 of file NuclearInteractionSimulator.h.

Referenced by NuclearInteractionSimulator(), and save().

std::ofstream NuclearInteractionSimulator::myOutputFile
private
std::vector< std::vector<TBranch*> > NuclearInteractionSimulator::theBranches
private

Definition at line 83 of file NuclearInteractionSimulator.h.

Referenced by NuclearInteractionSimulator().

std::vector< std::vector<unsigned> > NuclearInteractionSimulator::theCurrentEntry
private

Definition at line 85 of file NuclearInteractionSimulator.h.

Referenced by compute(), NuclearInteractionSimulator(), read(), and save().

std::vector< std::vector<unsigned> > NuclearInteractionSimulator::theCurrentInteraction
private

Definition at line 86 of file NuclearInteractionSimulator.h.

Referenced by compute(), NuclearInteractionSimulator(), read(), and save().

unsigned NuclearInteractionSimulator::theDistAlgo
private

Definition at line 94 of file NuclearInteractionSimulator.h.

Referenced by distanceToPrimary().

double NuclearInteractionSimulator::theDistCut
private

Definition at line 95 of file NuclearInteractionSimulator.h.

Referenced by compute().

std::vector< std::vector<std::string> > NuclearInteractionSimulator::theFileNames
private

Definition at line 89 of file NuclearInteractionSimulator.h.

Referenced by NuclearInteractionSimulator().

std::vector< std::vector<TFile*> > NuclearInteractionSimulator::theFiles
private
std::map< int, int > NuclearInteractionSimulator::theIDMap
private

Definition at line 92 of file NuclearInteractionSimulator.h.

Referenced by compute().

std::vector<double> NuclearInteractionSimulator::theLengthRatio
private

Definition at line 78 of file NuclearInteractionSimulator.h.

Referenced by compute().

std::vector< std::vector<NUEvent*> > NuclearInteractionSimulator::theNUEvents
private

Definition at line 84 of file NuclearInteractionSimulator.h.

Referenced by compute(), and NuclearInteractionSimulator().

std::vector< std::vector<unsigned> > NuclearInteractionSimulator::theNumberOfEntries
private

Definition at line 87 of file NuclearInteractionSimulator.h.

Referenced by compute(), and NuclearInteractionSimulator().

std::vector< std::vector<unsigned> > NuclearInteractionSimulator::theNumberOfInteractions
private

Definition at line 88 of file NuclearInteractionSimulator.h.

Referenced by compute(), and NuclearInteractionSimulator().

std::vector< std::vector<double> > NuclearInteractionSimulator::thePionCM
private

Definition at line 90 of file NuclearInteractionSimulator.h.

Referenced by compute(), and NuclearInteractionSimulator().

std::vector<double> NuclearInteractionSimulator::thePionEN
private

Definition at line 72 of file NuclearInteractionSimulator.h.

Referenced by compute(), and NuclearInteractionSimulator().

double NuclearInteractionSimulator::thePionEnergy
private

Definition at line 77 of file NuclearInteractionSimulator.h.

Referenced by compute().

std::vector<int> NuclearInteractionSimulator::thePionID
private

Definition at line 73 of file NuclearInteractionSimulator.h.

Referenced by index().

std::vector<double> NuclearInteractionSimulator::thePionMA
private

Definition at line 75 of file NuclearInteractionSimulator.h.

Referenced by NuclearInteractionSimulator().

std::vector<std::string> NuclearInteractionSimulator::thePionNA
private

Definition at line 74 of file NuclearInteractionSimulator.h.

Referenced by compute(), and NuclearInteractionSimulator().

std::vector<double> NuclearInteractionSimulator::thePionPMin
private

Definition at line 76 of file NuclearInteractionSimulator.h.

Referenced by compute().

std::vector< std::vector<double> > NuclearInteractionSimulator::theRatios
private

Definition at line 79 of file NuclearInteractionSimulator.h.

Referenced by compute().

std::vector< std::vector<TTree*> > NuclearInteractionSimulator::theTrees
private

Definition at line 82 of file NuclearInteractionSimulator.h.

Referenced by compute(), and NuclearInteractionSimulator().