CMS 3D CMS Logo

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 () override
 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 *) override
 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 155 of file NuclearInteractionFTFSimulator.cc.

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

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

Default Destructor.

Definition at line 247 of file NuclearInteractionFTFSimulator.cc.

References theLund, theStringDecay, theStringModel, and vect.

247  {
248 
249  delete theStringDecay;
250  delete theStringModel;
251  delete theLund;
252  delete vect;
253 }

Member Function Documentation

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

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

Implements MaterialEffectsSimulator.

Definition at line 255 of file NuclearInteractionFTFSimulator.cc.

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

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

References RawParticle::charge(), SoftLeptonByDistance_cfi::distance, theDistAlgo, and RawParticle::Vect().

Referenced by saveDaughter().

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

Definition at line 420 of file NuclearInteractionFTFSimulator.cc.

References MaterialEffectsSimulator::_theUpdatedState, SoftLeptonByDistance_cfi::distance, distanceToPrimary(), distMin, fact, training_settings::idx, makeParticle(), BaseParticlePropagator::particle(), ParticlePropagator::particleDataTable(), BPhysicsValidation_cfi::pdgid, MaterialEffectsSimulator::theClosestChargedDaughterId, theDistCut, and RawParticle::vertex().

Referenced by compute().

422 {
423  unsigned int idx = _theUpdatedState.size();
424  _theUpdatedState.emplace_back(makeParticle(Particle.particleDataTable(),
425  pdgid,
426  XYZTLorentzVector{lv.px()*fact,lv.py()*fact,lv.pz()*fact,lv.e()*fact},
427  Particle.particle().vertex()));
428 
429  // Store the closest daughter index (for later tracking purposes, so charged particles only)
431  // Find the closest daughter, if closer than a given upper limit.
432  if ( distance < distMin && distance < theDistCut ) {
433  distMin = distance;
435  }
436  // std::cout << _theUpdatedState[idx] << std::endl;
437 }
const HepPDT::ParticleDataTable * particleDataTable() const
RawParticle const & particle() const
The particle being propagated.
RawParticle makeParticle(HepPDT::ParticleDataTable const *, int id, const math::XYZTLorentzVector &p)
Definition: makeParticle.cc:29
double distanceToPrimary(const RawParticle &Particle, const RawParticle &aDaughter) const
const XYZTLorentzVector & vertex() const
the vertex fourvector
Definition: RawParticle.h:339
const double fact
std::vector< RawParticle > _theUpdatedState
math::XYZTLorentzVector XYZTLorentzVector
Definition: RawParticle.h:27

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().