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 | Static Private Attributes
NuclearInteractionFTFSimulator Class Reference

#include <NuclearInteractionFTFSimulator.h>

Inheritance diagram for NuclearInteractionFTFSimulator:
MaterialEffectsSimulator

Public Member Functions

 NuclearInteractionFTFSimulator (unsigned int distAlgo, double distCut, double elimit, double eth)
 Constructor. More...
 
 ~NuclearInteractionFTFSimulator ()
 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...
 
virtual void save ()
 Used by NuclearInteractionSimulator to save last sampled event. 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
 
void saveDaughter (ParticlePropagator &Particle, const G4LorentzVector &lv, int pdgid)
 

Private Attributes

G4LorentzVector curr4Mom
 
int currIdx
 
const G4ParticleDefinition * currParticle
 
G4Track * currTrack
 
double distMin
 
G4Step * dummyStep
 
size_t index
 
double intLengthElastic
 
double intLengthInelastic
 
G4Nucleus targetNucleus
 
G4CascadeInterface * theBertiniCascade
 
double theBertiniLimit
 
G4ThreeVector theBoost
 
G4GeneratorPrecompoundInterface * theCascade
 
G4DiffuseElastic * theDiffuseElastic
 
unsigned int theDistAlgo
 
double theDistCut
 
double theEnergyLimit
 
G4TheoFSGenerator * theHadronicModel
 
G4LundStringFragmentation * theLund
 
G4HadProjectile theProjectile
 
G4ExcitedStringDecay * theStringDecay
 
G4FTFModel * theStringModel
 
G4PhysicsLogVector * vect
 
G4ThreeVector vectProj
 

Static Private Attributes

static const G4ParticleDefinition * theG4Hadron [numHadrons] = {0}
 
static int theId [numHadrons] = {0}
 

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 40 of file NuclearInteractionFTFSimulator.h.

Constructor & Destructor Documentation

NuclearInteractionFTFSimulator::NuclearInteractionFTFSimulator ( unsigned int  distAlgo,
double  distCut,
double  elimit,
double  eth 
)

Constructor.

Definition at line 153 of file NuclearInteractionFTFSimulator.cc.

References currIdx, currParticle, currTrack, dummyStep, i, index, initializeOnce, intLengthElastic, intLengthInelastic, MeV, npoints, numHadrons, targetNucleus, theBertiniCascade, theDiffuseElastic, theEnergyLimit, theG4Hadron, theHadronicModel, theId, theLund, theStringDecay, theStringModel, and vect.

154  :
155  curr4Mom(0.,0.,0.,0.),
156  vectProj(0.,0.,1.),
157  theBoost(0.,0.,0.),
158  theBertiniLimit(elimit),
159  theEnergyLimit(eth),
161  distMin(1E99),
162  theDistAlgo(distAlgo)
163 {
164  // FTF model
165  theHadronicModel = new G4TheoFSGenerator("FTF");
166  theStringModel = new G4FTFModel();
167  G4GeneratorPrecompoundInterface* cascade
168  = new G4GeneratorPrecompoundInterface(new CMSDummyDeexcitation());
169  theLund = new G4LundStringFragmentation();
170  theStringDecay = new G4ExcitedStringDecay(theLund);
171  theStringModel->SetFragmentationModel(theStringDecay);
172 
173  theHadronicModel->SetTransport(cascade);
174  theHadronicModel->SetHighEnergyGenerator(theStringModel);
175  theHadronicModel->SetMinEnergy(theEnergyLimit);
176 
177  // Bertini Cascade
178  theBertiniCascade = new G4CascadeInterface();
179 
180  theDiffuseElastic = new G4DiffuseElastic();
181 
182  // Geant4 particles and cross sections
183  std::call_once(initializeOnce, [this] () {
184  theG4Hadron[0] = G4Proton::Proton();
185  theG4Hadron[1] = G4Neutron::Neutron();
186  theG4Hadron[2] = G4PionPlus::PionPlus();
187  theG4Hadron[3] = G4PionMinus::PionMinus();
188  theG4Hadron[4] = G4KaonPlus::KaonPlus();
189  theG4Hadron[5] = G4KaonMinus::KaonMinus();
190  theG4Hadron[6] = G4KaonZeroLong::KaonZeroLong();
191  theG4Hadron[7] = G4KaonZeroShort::KaonZeroShort();
192  theG4Hadron[8] = G4KaonZero::KaonZero();
193  theG4Hadron[9] = G4AntiKaonZero::AntiKaonZero();
194  theG4Hadron[10]= G4Lambda::Lambda();
195  theG4Hadron[11]= G4OmegaMinus::OmegaMinus();
196  theG4Hadron[12]= G4SigmaMinus::SigmaMinus();
197  theG4Hadron[13]= G4SigmaPlus::SigmaPlus();
198  theG4Hadron[14]= G4SigmaZero::SigmaZero();
199  theG4Hadron[15]= G4XiMinus::XiMinus();
200  theG4Hadron[16]= G4XiZero::XiZero();
201  theG4Hadron[17]= G4AntiProton::AntiProton();
202  theG4Hadron[18]= G4AntiNeutron::AntiNeutron();
203  theG4Hadron[19]= G4AntiLambda::AntiLambda();
204  theG4Hadron[20]= G4AntiOmegaMinus::AntiOmegaMinus();
205  theG4Hadron[21]= G4AntiSigmaMinus::AntiSigmaMinus();
206  theG4Hadron[22]= G4AntiSigmaPlus::AntiSigmaPlus();
207  theG4Hadron[23]= G4AntiSigmaZero::AntiSigmaZero();
208  theG4Hadron[24]= G4AntiXiMinus::AntiXiMinus();
209  theG4Hadron[25]= G4AntiXiZero::AntiXiZero();
210  theG4Hadron[26]= G4AntiAlpha::AntiAlpha();
211  theG4Hadron[27]= G4AntiDeuteron::AntiDeuteron();
212  theG4Hadron[28]= G4AntiTriton::AntiTriton();
213  theG4Hadron[29]= G4AntiHe3::AntiHe3();
214 
215  // other Geant4 particles
216  G4ParticleDefinition* ion = G4GenericIon::GenericIon();
217  ion->SetProcessManager(new G4ProcessManager(ion));
218  G4DecayPhysics decays;
219  decays.ConstructParticle();
220  G4ParticleTable* partTable = G4ParticleTable::GetParticleTable();
221  partTable->SetVerboseLevel(0);
222  partTable->SetReadiness();
223 
224  for(int i=0; i<numHadrons; ++i) {
225  theId[i] = theG4Hadron[i]->GetPDGEncoding();
226  }
227  });
228 
229  // local objects
230  vect = new G4PhysicsLogVector(npoints-1,100*MeV,TeV);
232  currIdx = 0;
233  index = 0;
234  currTrack = 0;
236 
237  // fill projectile particle definitions
238  dummyStep = new G4Step();
239  dummyStep->SetPreStepPoint(new G4StepPoint());
240 
241  // target is always Silicon
242  targetNucleus.SetParameters(28, 14);
243 }
int i
Definition: DBlmapReader.cc:9
const double MeV
static std::once_flag initializeOnce
static const int npoints
static const int numHadrons
const G4ParticleDefinition * currParticle
static const G4ParticleDefinition * theG4Hadron[numHadrons]
NuclearInteractionFTFSimulator::~NuclearInteractionFTFSimulator ( )

Default Destructor.

Definition at line 245 of file NuclearInteractionFTFSimulator.cc.

References theLund, theStringDecay, theStringModel, and vect.

245  {
246 
247  delete theStringDecay;
248  delete theStringModel;
249  delete theLund;
250  delete vect;
251 }

Member Function Documentation

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

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

Implements MaterialEffectsSimulator.

Definition at line 253 of file NuclearInteractionFTFSimulator.cc.

References MaterialEffectsSimulator::_theUpdatedState, corrfactors_el, corrfactors_inel, funct::cos(), prof2calltree::cost, curr4Mom, currIdx, currParticle, currTrack, dir, distMin, dummyStep, alignCSCRings::e, RandomEngineAndDistribution::flatShoot(), GeV, index, intLengthElastic, intLengthInelastic, j, RawParticle::momentum(), npoints, nuclElLength, nuclInelLength, numHadrons, phi, RawParticle::pid(), MaterialEffectsSimulator::radLengths, query::result, saveDaughter(), funct::sin(), mathSSE::sqrt(), targetNucleus, theBertiniCascade, theBertiniLimit, theBoost, MaterialEffectsSimulator::theClosestChargedDaughterId, theDiffuseElastic, theDistCut, theEnergyLimit, theG4Hadron, theHadronicModel, theId, theProjectile, vect, and vectProj.

255 {
256  //std::cout << "#### Primary " << Particle.pid() << " E(GeV)= "
257  // << Particle.momentum().e() << std::endl;
258 
259  int thePid = Particle.pid();
260  if(thePid != theId[currIdx]) {
261  currParticle = 0;
262  currIdx = 0;
263  for(; currIdx<numHadrons; ++currIdx) {
264  if(theId[currIdx] == thePid) {
266  // neutral kaons
267  if(7 == currIdx || 8 == currIdx) {
269  if(random->flatShoot() > 0.5) { currParticle = theG4Hadron[10]; }
270  }
271  break;
272  }
273  }
274  }
275  if(!currParticle) { return; }
276 
277  // fill projectile for Geant4
278  double mass = currParticle->GetPDGMass();
279  double ekin = CLHEP::GeV*Particle.momentum().e() - mass;
280 
281  // check interaction length
284  if(0.0 == intLengthInelastic) { return; }
285 
286  // apply corrections
287  if(ekin <= vect->Energy(0)) {
290  } else if(ekin < vect->Energy(npoints-1)) {
291  index = vect->FindBin(ekin, index);
292  double e1 = vect->Energy(index);
293  double e2 = vect->Energy(index+1);
294  intLengthElastic *= ((corrfactors_el[currIdx][index]*(e2 - ekin) +
295  corrfactors_el[currIdx][index+1]*(ekin - e1))/(e2 - e1));
296  intLengthInelastic *= ((corrfactors_inel[currIdx][index]*(e2 - ekin) +
297  corrfactors_inel[currIdx][index+1]*(ekin - e1))/(e2 - e1));
298  }
299  /*
300  std::cout << " Primary " << currParticle->GetParticleName()
301  << " E(GeV)= " << e*fact << std::endl;
302  */
303 
304  double currInteractionLength = -G4Log(random->flatShoot())*intLengthElastic*intLengthInelastic
306  /*
307  std::cout << "*NuclearInteractionFTFSimulator::compute: R(X0)= " << radLengths
308  << " Rnuc(X0)= " << theNuclIntLength[currIdx] << " IntLength(X0)= "
309  << currInteractionLength << std::endl;
310  */
311  // Check position of nuclear interaction
312  if (currInteractionLength > radLengths) { return; }
313 
314  // fill projectile for Geant4
315  double px = Particle.momentum().px();
316  double py = Particle.momentum().py();
317  double pz = Particle.momentum().pz();
318  double ptot = sqrt(px*px + py*py + pz*pz);
319  double norm = 1./ptot;
320  G4ThreeVector dir(px*norm, py*norm, pz*norm);
321  /*
322  std::cout << " Primary " << currParticle->GetParticleName()
323  << " E(GeV)= " << e*fact << " P(GeV/c)= ("
324  << px << " " << py << " " << pz << ")" << std::endl;
325  */
326 
327  G4DynamicParticle* dynParticle = new G4DynamicParticle(theG4Hadron[currIdx],dir,ekin);
328  currTrack = new G4Track(dynParticle, 0.0, vectProj);
329  currTrack->SetStep(dummyStep);
330 
331  theProjectile.Initialise(*currTrack);
332  delete currTrack;
333 
334  G4HadFinalState* result;
335 
336  // elastic interaction
338 
339  result = theDiffuseElastic->ApplyYourself(theProjectile, targetNucleus);
340  G4ThreeVector dirnew = result->GetMomentumChange().unit();
341  double cost = (dir*dirnew);
342  double sint = std::sqrt((1. - cost)*(1. + cost));
343 
344  curr4Mom.set(ptot*dirnew.x(),ptot*dirnew.y(),ptot*dirnew.z(),Particle.momentum().e());
345 
346  // Always create a daughter if the kink is large engough
347  if (sint > theDistCut) {
348  saveDaughter(Particle, curr4Mom, thePid);
349  } else {
350  Particle.SetXYZT(curr4Mom.px(), curr4Mom.py(), curr4Mom.pz(), curr4Mom.e());
351  }
352 
353  // inelastic interaction
354  } else {
355 
356  // Bertini cascade for low-energy hadrons (except light anti-nuclei)
357  // FTFP is applied above energy limit and for all anti-hyperons and anti-ions
358  if(ekin <= theBertiniLimit && currIdx < 17) {
359  result = theBertiniCascade->ApplyYourself(theProjectile, targetNucleus);
360  } else {
361  result = theHadronicModel->ApplyYourself(theProjectile, targetNucleus);
362  }
363  if(result) {
364 
365  int nsec = result->GetNumberOfSecondaries();
366  if(0 < nsec) {
367 
368  result->SetTrafoToLab(theProjectile.GetTrafoToLab());
369  _theUpdatedState.clear();
370 
371  //std::cout << " " << nsec << " secondaries" << std::endl;
372  // Generate angle
373  double phi = random->flatShoot()*CLHEP::twopi;
375  distMin = 1e99;
376 
377  // rotate and store secondaries
378  for (int j=0; j<nsec; ++j) {
379 
380  const G4DynamicParticle* dp = result->GetSecondary(j)->GetParticle();
381  int thePid = dp->GetParticleDefinition()->GetPDGEncoding();
382 
383  // rotate around primary direction
384  curr4Mom = dp->Get4Momentum();
385  curr4Mom.rotate(phi, vectProj);
386  curr4Mom *= result->GetTrafoToLab();
387  /*
388  std::cout << j << ". " << dp->GetParticleDefinition()->GetParticleName()
389  << " " << thePid
390  << " " << curr4Mom*fact << std::endl;
391  */
392  // prompt 2-gamma decay for pi0, eta, eta'p
393  if(111 == thePid || 221 == thePid || 331 == thePid) {
394  theBoost = curr4Mom.boostVector();
395  double e = 0.5*dp->GetParticleDefinition()->GetPDGMass();
396  double fi = random->flatShoot()*CLHEP::twopi;
397  double cth = 2*random->flatShoot() - 1.0;
398  double sth = sqrt((1.0 - cth)*(1.0 + cth));
399  G4LorentzVector lv(e*sth*cos(fi),e*sth*sin(fi),e*cth,e);
400  lv.boost(theBoost);
401  saveDaughter(Particle, lv, 22);
402  curr4Mom -= lv;
403  if(curr4Mom.e() > theEnergyLimit) {
404  saveDaughter(Particle, curr4Mom, 22);
405  }
406  } else {
407  if(curr4Mom.e() > theEnergyLimit + dp->GetParticleDefinition()->GetPDGMass()) {
408  saveDaughter(Particle, curr4Mom, thePid);
409  }
410  }
411  }
412  result->Clear();
413  }
414  }
415  }
416 }
const double nuclInelLength[numHadrons]
const double GeV
Definition: MathUtil.h:16
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
int pid() const
get the HEP particle ID number
Definition: RawParticle.h:265
TRandom random
Definition: MVATrainer.cc:138
const double corrfactors_inel[numHadrons][npoints]
const XYZTLorentzVector & momentum() const
the momentum fourvector
Definition: RawParticle.h:286
void saveDaughter(ParticlePropagator &Particle, const G4LorentzVector &lv, int pdgid)
T sqrt(T t)
Definition: SSEVec.h:48
tuple result
Definition: query.py:137
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
int j
Definition: DBlmapReader.cc:9
const double nuclElLength[numHadrons]
static const int npoints
static const int numHadrons
const double corrfactors_el[numHadrons][npoints]
const G4ParticleDefinition * currParticle
static const G4ParticleDefinition * theG4Hadron[numHadrons]
std::vector< RawParticle > _theUpdatedState
dbl *** dir
Definition: mlp_gen.cc:35
double NuclearInteractionFTFSimulator::distanceToPrimary ( const RawParticle Particle,
const RawParticle aDaughter 
) const
private

Definition at line 437 of file NuclearInteractionFTFSimulator.cc.

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

Referenced by saveDaughter().

439 {
440  double distance = 2E99;
441  // Compute the distance only for charged primaries
442  if ( fabs(Particle.charge()) > 1E-6 ) {
443 
444  // The secondary must have the same charge
445  double chargeDiff = fabs(aDaughter.charge()-Particle.charge());
446  if ( fabs(chargeDiff) < 1E-6 ) {
447 
448  // Here are two distance definitions * to be tuned *
449  switch ( theDistAlgo ) {
450 
451  case 1:
452  // sin(theta12)
453  distance = (aDaughter.Vect().Unit().Cross(Particle.Vect().Unit())).R();
454  break;
455 
456  case 2:
457  // sin(theta12) * p1/p2
458  distance = (aDaughter.Vect().Cross(Particle.Vect())).R()
459  /aDaughter.Vect().Mag2();
460  break;
461 
462  default:
463  // Should not happen
464  break;
465  }
466  }
467  }
468  return distance;
469 }
double charge() const
get the MEASURED charge
Definition: RawParticle.h:282
void NuclearInteractionFTFSimulator::saveDaughter ( ParticlePropagator Particle,
const G4LorentzVector &  lv,
int  pdgid 
)
private

Definition at line 418 of file NuclearInteractionFTFSimulator.cc.

References MaterialEffectsSimulator::_theUpdatedState, HLT_25ns14e33_v1_cff::distance, distanceToPrimary(), distMin, fact, customizeTrackingMonitorSeedNumber::idx, MaterialEffectsSimulator::theClosestChargedDaughterId, and theDistCut.

Referenced by compute().

420 {
421  unsigned int idx = _theUpdatedState.size();
422  _theUpdatedState.push_back(Particle);
423  _theUpdatedState[idx].SetXYZT(lv.px()*fact,lv.py()*fact,lv.pz()*fact,lv.e()*fact);
424  _theUpdatedState[idx].setID(pdgid);
425 
426  // Store the closest daughter index (for later tracking purposes, so charged particles only)
427  double distance = distanceToPrimary(Particle,_theUpdatedState[idx]);
428  // Find the closest daughter, if closer than a given upper limit.
429  if ( distance < distMin && distance < theDistCut ) {
430  distMin = distance;
432  }
433  // std::cout << _theUpdatedState[idx] << std::endl;
434 }
double distanceToPrimary(const RawParticle &Particle, const RawParticle &aDaughter) const
const double fact
tuple idx
DEBUGGING if hasattr(process,&quot;trackMonIterativeTracking2012&quot;): print &quot;trackMonIterativeTracking2012 D...
std::vector< RawParticle > _theUpdatedState

Member Data Documentation

G4LorentzVector NuclearInteractionFTFSimulator::curr4Mom
private

Definition at line 80 of file NuclearInteractionFTFSimulator.h.

Referenced by compute().

int NuclearInteractionFTFSimulator::currIdx
private

Definition at line 93 of file NuclearInteractionFTFSimulator.h.

Referenced by compute(), and NuclearInteractionFTFSimulator().

const G4ParticleDefinition* NuclearInteractionFTFSimulator::currParticle
private

Definition at line 76 of file NuclearInteractionFTFSimulator.h.

Referenced by compute(), and NuclearInteractionFTFSimulator().

G4Track* NuclearInteractionFTFSimulator::currTrack
private

Definition at line 75 of file NuclearInteractionFTFSimulator.h.

Referenced by compute(), and NuclearInteractionFTFSimulator().

double NuclearInteractionFTFSimulator::distMin
private

Definition at line 88 of file NuclearInteractionFTFSimulator.h.

Referenced by compute(), and saveDaughter().

G4Step* NuclearInteractionFTFSimulator::dummyStep
private

Definition at line 74 of file NuclearInteractionFTFSimulator.h.

Referenced by compute(), and NuclearInteractionFTFSimulator().

size_t NuclearInteractionFTFSimulator::index
private
double NuclearInteractionFTFSimulator::intLengthElastic
private

Definition at line 90 of file NuclearInteractionFTFSimulator.h.

Referenced by compute(), and NuclearInteractionFTFSimulator().

double NuclearInteractionFTFSimulator::intLengthInelastic
private

Definition at line 91 of file NuclearInteractionFTFSimulator.h.

Referenced by compute(), and NuclearInteractionFTFSimulator().

G4Nucleus NuclearInteractionFTFSimulator::targetNucleus
private

Definition at line 78 of file NuclearInteractionFTFSimulator.h.

Referenced by compute(), and NuclearInteractionFTFSimulator().

G4CascadeInterface* NuclearInteractionFTFSimulator::theBertiniCascade
private

Definition at line 71 of file NuclearInteractionFTFSimulator.h.

Referenced by compute(), and NuclearInteractionFTFSimulator().

double NuclearInteractionFTFSimulator::theBertiniLimit
private

Definition at line 84 of file NuclearInteractionFTFSimulator.h.

Referenced by compute().

G4ThreeVector NuclearInteractionFTFSimulator::theBoost
private

Definition at line 82 of file NuclearInteractionFTFSimulator.h.

Referenced by compute().

G4GeneratorPrecompoundInterface* NuclearInteractionFTFSimulator::theCascade
private

Definition at line 69 of file NuclearInteractionFTFSimulator.h.

G4DiffuseElastic* NuclearInteractionFTFSimulator::theDiffuseElastic
private

Definition at line 72 of file NuclearInteractionFTFSimulator.h.

Referenced by compute(), and NuclearInteractionFTFSimulator().

unsigned int NuclearInteractionFTFSimulator::theDistAlgo
private

Definition at line 95 of file NuclearInteractionFTFSimulator.h.

Referenced by distanceToPrimary().

double NuclearInteractionFTFSimulator::theDistCut
private

Definition at line 87 of file NuclearInteractionFTFSimulator.h.

Referenced by compute(), and saveDaughter().

double NuclearInteractionFTFSimulator::theEnergyLimit
private

Definition at line 85 of file NuclearInteractionFTFSimulator.h.

Referenced by compute(), and NuclearInteractionFTFSimulator().

const G4ParticleDefinition * NuclearInteractionFTFSimulator::theG4Hadron = {0}
staticprivate

Definition at line 61 of file NuclearInteractionFTFSimulator.h.

Referenced by compute(), and NuclearInteractionFTFSimulator().

G4TheoFSGenerator* NuclearInteractionFTFSimulator::theHadronicModel
private

Definition at line 65 of file NuclearInteractionFTFSimulator.h.

Referenced by compute(), and NuclearInteractionFTFSimulator().

int NuclearInteractionFTFSimulator::theId = {0}
staticprivate

Definition at line 62 of file NuclearInteractionFTFSimulator.h.

Referenced by compute(), and NuclearInteractionFTFSimulator().

G4LundStringFragmentation* NuclearInteractionFTFSimulator::theLund
private
G4HadProjectile NuclearInteractionFTFSimulator::theProjectile
private

Definition at line 79 of file NuclearInteractionFTFSimulator.h.

Referenced by compute().

G4ExcitedStringDecay* NuclearInteractionFTFSimulator::theStringDecay
private
G4FTFModel* NuclearInteractionFTFSimulator::theStringModel
private
G4PhysicsLogVector* NuclearInteractionFTFSimulator::vect
private
G4ThreeVector NuclearInteractionFTFSimulator::vectProj
private

Definition at line 81 of file NuclearInteractionFTFSimulator.h.

Referenced by compute().