Geant4  10.02.p03
G4DNASancheExcitationModel Class Reference

#include <G4DNASancheExcitationModel.hh>

Inheritance diagram for G4DNASancheExcitationModel:
Collaboration diagram for G4DNASancheExcitationModel:

Public Member Functions

 G4DNASancheExcitationModel (const G4ParticleDefinition *p=0, const G4String &nam="DNASancheExcitationModel")
 
virtual ~G4DNASancheExcitationModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
 
virtual G4double CrossSectionPerVolume (const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle *> *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
 
G4double PartialCrossSection (G4double energy, G4int level)
 
G4double TotalCrossSection (G4double t)
 
void ExtendLowEnergyLimit (G4double)
 
void SetVerboseLevel (int verbose)
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel)
 
virtual void InitialiseForMaterial (const G4ParticleDefinition *, const G4Material *)
 
virtual void InitialiseForElement (const G4ParticleDefinition *, G4int Z)
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
virtual G4double GetPartialCrossSection (const G4Material *, G4int, const G4ParticleDefinition *, G4double)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ComputeCrossSectionPerShell (const G4ParticleDefinition *, G4int Z, G4int shellIdx, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ChargeSquareRatio (const G4Track &)
 
virtual G4double GetChargeSquareRatio (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual G4double GetParticleCharge (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void StartTracking (G4Track *)
 
virtual void CorrectionsAlongStep (const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double &eloss, G4double &niel, G4double length)
 
virtual G4double Value (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *, G4double cut=0.0)
 
virtual G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
virtual void ModelDescription (std::ostream &outFile) const
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
std::vector< G4EmElementSelector * > * GetElementSelectors ()
 
void SetElementSelectors (std::vector< G4EmElementSelector *> *)
 
G4double ComputeDEDX (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
G4double CrossSection (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeMeanFreePath (const G4ParticleDefinition *, G4double kineticEnergy, const G4Material *, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, const G4Element *, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4int SelectIsotopeNumber (const G4Element *)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4int SelectRandomAtomNumber (const G4Material *)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=0)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
void SetAngularDistribution (G4VEmAngularDistribution *)
 
G4double HighEnergyLimit () const
 
G4double LowEnergyLimit () const
 
G4double HighEnergyActivationLimit () const
 
G4double LowEnergyActivationLimit () const
 
G4double PolarAngleLimit () const
 
G4double SecondaryThreshold () const
 
G4bool LPMFlag () const
 
G4bool DeexcitationFlag () const
 
G4bool ForceBuildTableFlag () const
 
G4bool UseAngularGeneratorFlag () const
 
void SetAngularGeneratorFlag (G4bool)
 
void SetHighEnergyLimit (G4double)
 
void SetLowEnergyLimit (G4double)
 
void SetActivationHighEnergyLimit (G4double)
 
void SetActivationLowEnergyLimit (G4double)
 
G4bool IsActive (G4double kinEnergy)
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetLPMFlag (G4bool val)
 
void SetDeexcitationFlag (G4bool val)
 
void SetForceBuildTable (G4bool val)
 
void SetMasterThread (G4bool val)
 
G4bool IsMaster () const
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
const G4ElementGetCurrentElement () const
 
const G4IsotopeGetCurrentIsotope () const
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 

Protected Attributes

G4ParticleChangeForGamma * fParticleChangeForGamma
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChange * pParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 

Private Member Functions

G4int RandomSelect (G4double energy)
 
G4double VibrationEnergy (G4int level)
 
G4double Sum (G4double k)
 
G4double LinInterpolate (G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2)
 
G4DNASancheExcitationModeloperator= (const G4DNASancheExcitationModel &right)
 
 G4DNASancheExcitationModel (const G4DNASancheExcitationModel &)
 

Private Attributes

const std::vector< G4double > * fpWaterDensity
 
G4double lowEnergyLimit
 
G4double highEnergyLimit
 
G4bool isInitialised
 
G4int verboseLevel
 
G4int nLevels
 
std::vector< double > tdummyVec
 
std::vector< std::vector< double > > fEnergyLevelXS
 
std::vector< double > fEnergyTotalXS
 

Additional Inherited Members

- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLoss * GetParticleChangeForLoss ()
 
G4ParticleChangeForGamma * GetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 
- Static Protected Attributes inherited from G4VEmModel
static const G4double inveplus = 1.0/CLHEP::eplus
 

Detailed Description

Definition at line 45 of file G4DNASancheExcitationModel.hh.

Constructor & Destructor Documentation

◆ G4DNASancheExcitationModel() [1/2]

G4DNASancheExcitationModel::G4DNASancheExcitationModel ( const G4ParticleDefinition p = 0,
const G4String nam = "DNASancheExcitationModel" 
)

Definition at line 41 of file G4DNASancheExcitationModel.cc.

42  :
43  G4VEmModel(nam), isInitialised(false)
44 {
45  fpWaterDensity = 0;
46 
47  lowEnergyLimit = 2 * eV;
48  highEnergyLimit = 100 * eV;
51  nLevels = 9;
52 
53  verboseLevel = 0;
54  // Verbosity scale:
55  // 0 = nothing
56  // 1 = warning for energy non-conservation
57  // 2 = details of energy budget
58  // 3 = calculation of cross sections, file openings, sampling of atoms
59  // 4 = entering in methods
60 
61  if (verboseLevel > 0)
62  {
63  G4cout << "Sanche Excitation model is constructed " << G4endl<< "Energy range: "
64  << lowEnergyLimit / eV << " eV - "
65  << highEnergyLimit / eV << " eV"
66  << G4endl;
67  }
69  fpWaterDensity = 0;
70 }
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:69
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:725
G4GLOB_DLL std::ostream G4cout
const std::vector< G4double > * fpWaterDensity
G4ParticleChangeForGamma * fParticleChangeForGamma
static const double eV
Definition: G4SIunits.hh:212
#define G4endl
Definition: G4ios.hh:61
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:732
Here is the call graph for this function:

◆ ~G4DNASancheExcitationModel()

G4DNASancheExcitationModel::~G4DNASancheExcitationModel ( )
virtual

Definition at line 74 of file G4DNASancheExcitationModel.cc.

75 {
76 }

◆ G4DNASancheExcitationModel() [2/2]

G4DNASancheExcitationModel::G4DNASancheExcitationModel ( const G4DNASancheExcitationModel )
private

Member Function Documentation

◆ CrossSectionPerVolume()

G4double G4DNASancheExcitationModel::CrossSectionPerVolume ( const G4Material material,
const G4ParticleDefinition p,
G4double  ekin,
G4double  emin,
G4double  emax 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 168 of file G4DNASancheExcitationModel.cc.

173 {
174  if (verboseLevel > 3)
175  G4cout << "Calling CrossSectionPerVolume() of G4DNASancheExcitationModel"
176  << G4endl;
177 
178  // Calculate total cross section for model
179 
180  G4double sigma=0;
181 
182  G4double waterDensity = (*fpWaterDensity)[material->GetIndex()];
183 
184  if(waterDensity!= 0.0)
185  {
186 
187 // if (particleDefinition == G4Electron::ElectronDefinition()) // pas besoin
188  {
189  if (ekin >= lowEnergyLimit && ekin < highEnergyLimit)
190  // for now this is necessary
191  {
192  //sigma = Sum(ekin);
193  sigma = TotalCrossSection(ekin);
194  }
195  }
196 
197  if (verboseLevel > 2)
198  {
199  G4cout << "__________________________________" << G4endl;
200  G4cout << "=== G4DNASancheExcitationModel - XS INFO START" << G4endl;
201  G4cout << "=== Kinetic energy(eV)=" << ekin/eV << " particle : " << particleDefinition->GetParticleName() << G4endl;
202  G4cout << "=== Cross section per water molecule (cm^2)=" << sigma/cm/cm << G4endl;
203  G4cout << "=== Cross section per water molecule (cm^-1)=" << sigma*waterDensity/(1./cm) << G4endl;
204  G4cout << "=== G4DNASancheExcitationModel - XS INFO END" << G4endl;
205  }
206 
207  } // if water
208 
209  return sigma*2*waterDensity;
210  // see papers for factor 2 description
211 
212 }
static const double cm
Definition: G4SIunits.hh:118
size_t GetIndex() const
Definition: G4Material.hh:262
G4GLOB_DLL std::ostream G4cout
static const double eV
Definition: G4SIunits.hh:212
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:

◆ ExtendLowEnergyLimit()

void G4DNASancheExcitationModel::ExtendLowEnergyLimit ( G4double  threshold)
inline

Definition at line 121 of file G4DNASancheExcitationModel.hh.

122 {
123  lowEnergyLimit = threshold;
124  if(lowEnergyLimit < 2 * CLHEP::eV)
125  G4Exception("*** WARNING : the G4DNASancheExcitationModel class is not "
126  "validated below 2 eV !",
127  "", JustWarning, "");
128 }
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static const double eV
Here is the call graph for this function:
Here is the caller graph for this function:

◆ Initialise()

void G4DNASancheExcitationModel::Initialise ( const G4ParticleDefinition ,
const G4DataVector  
)
virtual

Implements G4VEmModel.

Definition at line 80 of file G4DNASancheExcitationModel.cc.

82 {
83 
84  if (verboseLevel > 3)
85  G4cout << "Calling G4DNASancheExcitationModel::Initialise()" << G4endl;
86 
87  // Energy limits
88 
90  {
91  G4cout << "G4DNASancheExcitationModel: low energy limit increased from " <<
92  LowEnergyLimit()/eV << " eV to " << lowEnergyLimit/eV << " eV" << G4endl;
94  }
95 
97  {
98  G4cout << "G4DNASancheExcitationModel: high energy limit decreased from " <<
99  HighEnergyLimit()/eV << " eV to " << highEnergyLimit/eV << " eV" << G4endl;
101  }
102 
103  //
104 
105  if (verboseLevel > 0)
106  {
107  G4cout << "Sanche Excitation model is initialized " << G4endl
108  << "Energy range: "
109  << LowEnergyLimit() / eV << " eV - "
110  << HighEnergyLimit() / eV << " eV"
111  << G4endl;
112  }
113 
114  // Initialize water density pointer
116  GetNumMolPerVolTableFor(G4Material::GetMaterial("G4_WATER"));
117 
118  if (isInitialised) {return;} // RETURNS HERE
119 
121  isInitialised = true;
122 
123  char *path = getenv("G4LEDATA");
124  std::ostringstream eFullFileName;
125  eFullFileName << path << "/dna/sigma_excitationvib_e_sanche.dat";
126  std::ifstream input(eFullFileName.str().c_str());
127 
128  if (!input)
129  {
130  G4Exception("G4DNASancheExcitationModel::Initialise","em0003",
131  FatalException,"Missing data file:/dna/sigma_excitationvib_e_sanche.dat");
132  }
133 
134  // March 25th, 2014 - Vaclav Stepan, Sebastien Incerti
135  // Added clear for MT
136  tdummyVec.clear();
137  //
138 
139  double t;
140  double xs;
141 
142  while(!input.eof())
143  {
144  input>>t;
145  tdummyVec.push_back(t);
146 
147  fEnergyLevelXS.push_back(std::vector<double>());
148  fEnergyTotalXS.push_back(0);
149  std::vector<double>& levelXS = fEnergyLevelXS.back();
150  levelXS.reserve(9);
151 
152 // G4cout<<t;
153 
154  for(size_t i = 0 ; i < 9 ;++i)
155  {
156  input>>xs;
157  levelXS.push_back(xs);
158  fEnergyTotalXS.back() += xs;
159 // G4cout <<" " << levelXS[i];
160  }
161 
162 // G4cout << G4endl;
163  }
164 }
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:641
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:604
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:725
G4GLOB_DLL std::ostream G4cout
const std::vector< G4double > * fpWaterDensity
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:634
G4ParticleChangeForGamma * fParticleChangeForGamma
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static const double eV
Definition: G4SIunits.hh:212
static G4DNAMolecularMaterial * Instance()
#define G4endl
Definition: G4ios.hh:61
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:732
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:134
std::vector< std::vector< double > > fEnergyLevelXS
Here is the call graph for this function:

◆ LinInterpolate()

G4double G4DNASancheExcitationModel::LinInterpolate ( G4double  e1,
G4double  e2,
G4double  e,
G4double  xs1,
G4double  xs2 
)
private

Definition at line 372 of file G4DNASancheExcitationModel.cc.

377 {
378  G4double a = (xs2 - xs1) / (e2 - e1);
379  G4double b = xs2 - a * e2;
380  G4double value = a * e + b;
381  // G4cout<<"interP >> "<<e1<<" "<<e2<<" "<<e<<" "<<xs1<<" "<<xs2<<" "<<a<<" "<<b<<" "<<value<<G4endl;
382 
383  return value;
384 }
static const G4double e2
static const G4double e1
double G4double
Definition: G4Types.hh:76
Here is the caller graph for this function:

◆ operator=()

G4DNASancheExcitationModel& G4DNASancheExcitationModel::operator= ( const G4DNASancheExcitationModel right)
private

◆ PartialCrossSection()

G4double G4DNASancheExcitationModel::PartialCrossSection ( G4double  energy,
G4int  level 
)

Definition at line 263 of file G4DNASancheExcitationModel.cc.

265 {
266  std::vector<double>::iterator t2 = std::upper_bound(tdummyVec.begin(),
267  tdummyVec.end(), t / eV);
268  std::vector<double>::iterator t1 = t2 - 1;
269 
270  size_t i1 = t1 - tdummyVec.begin();
271  size_t i2 = t2 - tdummyVec.begin();
272 
273  double sigma = LinInterpolate((*t1), (*t2),
274  t / eV,
275  fEnergyLevelXS[i1][level],
276  fEnergyLevelXS[i2][level]);
277 
278  static const double conv_factor = 1e-16 * cm * cm;
279 
280  sigma *= conv_factor;
281  if (sigma == 0.) sigma = 1e-30;
282  return (sigma);
283 }
static const double cm
Definition: G4SIunits.hh:118
TTree * t1
Definition: plottest35.C:26
G4double LinInterpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2)
static const double eV
Definition: G4SIunits.hh:212
TTree * t2
Definition: plottest35.C:36
std::vector< std::vector< double > > fEnergyLevelXS
Here is the call graph for this function:
Here is the caller graph for this function:

◆ RandomSelect()

G4int G4DNASancheExcitationModel::RandomSelect ( G4double  energy)
private

Definition at line 319 of file G4DNASancheExcitationModel.cc.

320 {
321 
322  // Level Selection Counting can be done here !
323 
324  G4int i = nLevels;
325  G4double value = 0.;
326  std::deque<double> values;
327 
328  while (i > 0)
329  {
330  i--;
331  G4double partial = PartialCrossSection(k, i);
332  values.push_front(partial);
333  value += partial;
334  }
335 
336  value *= G4UniformRand();
337 
338  i = nLevels;
339 
340  while (i > 0)
341  {
342  i--;
343  if (values[i] > value)
344  {
345  //outcount<<i<<" "<<VibrationEnergy(i)<<G4endl;
346  return i;
347  }
348  value -= values[i];
349  }
350 
351  //outcount<<0<<" "<<VibrationEnergy(0)<<G4endl;
352 
353  return 0;
354 }
G4double PartialCrossSection(G4double energy, G4int level)
int G4int
Definition: G4Types.hh:78
#define G4UniformRand()
Definition: Randomize.hh:97
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SampleSecondaries()

void G4DNASancheExcitationModel::SampleSecondaries ( std::vector< G4DynamicParticle *> *  ,
const G4MaterialCutsCouple ,
const G4DynamicParticle aDynamicElectron,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Definition at line 216 of file G4DNASancheExcitationModel.cc.

222 {
223 
224  if (verboseLevel > 3)
225  G4cout << "Calling SampleSecondaries() of G4DNASancheExcitationModel"
226  << G4endl;
227 
228  G4double electronEnergy0 = aDynamicElectron->GetKineticEnergy();
229  G4int level = RandomSelect(electronEnergy0);
230  G4double excitationEnergy = VibrationEnergy(level); // levels go from 0 to 8
231  G4double newEnergy = electronEnergy0 - excitationEnergy;
232 
233  /*
234  if (electronEnergy0 < highEnergyLimit)
235  {
236  if (newEnergy >= lowEnergyLimit)
237  {
238  fParticleChangeForGamma->ProposeMomentumDirection(aDynamicElectron->GetMomentumDirection());
239  fParticleChangeForGamma->SetProposedKineticEnergy(newEnergy);
240  fParticleChangeForGamma->ProposeLocalEnergyDeposit(excitationEnergy);
241  }
242 
243  else
244  {
245  fParticleChangeForGamma->ProposeTrackStatus(fStopAndKill);
246  fParticleChangeForGamma->ProposeLocalEnergyDeposit(electronEnergy0);
247  }
248  }
249  */
250 
251  if (electronEnergy0 < highEnergyLimit && newEnergy>0.)
252  {
253  fParticleChangeForGamma->ProposeMomentumDirection(aDynamicElectron->GetMomentumDirection());
254  fParticleChangeForGamma->SetProposedKineticEnergy(newEnergy);
255  fParticleChangeForGamma->ProposeLocalEnergyDeposit(excitationEnergy);
256  }
257 
258  //
259 }
int G4int
Definition: G4Types.hh:78
G4double GetKineticEnergy() const
G4GLOB_DLL std::ostream G4cout
G4ParticleChangeForGamma * fParticleChangeForGamma
const G4ThreeVector & GetMomentumDirection() const
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:

◆ SetVerboseLevel()

void G4DNASancheExcitationModel::SetVerboseLevel ( int  verbose)
inline

Definition at line 76 of file G4DNASancheExcitationModel.hh.

77  {
78  verboseLevel = verbose;
79  }

◆ Sum()

G4double G4DNASancheExcitationModel::Sum ( G4double  k)
private

Definition at line 358 of file G4DNASancheExcitationModel.cc.

359 {
360  G4double totalCrossSection = 0.;
361 
362  for (G4int i = 0; i < nLevels; i++)
363  {
364  totalCrossSection += PartialCrossSection(k, i);
365  }
366 
367  return totalCrossSection;
368 }
G4double PartialCrossSection(G4double energy, G4int level)
int G4int
Definition: G4Types.hh:78
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:

◆ TotalCrossSection()

G4double G4DNASancheExcitationModel::TotalCrossSection ( G4double  t)

Definition at line 287 of file G4DNASancheExcitationModel.cc.

288 {
289  std::vector<double>::iterator t2 = std::upper_bound(tdummyVec.begin(),
290  tdummyVec.end(), t / eV);
291  std::vector<double>::iterator t1 = t2 - 1;
292 
293  size_t i1 = t1 - tdummyVec.begin();
294  size_t i2 = t2 - tdummyVec.begin();
295 
296  double sigma = LinInterpolate((*t1), (*t2),
297  t / eV,
298  fEnergyTotalXS[i1],
299  fEnergyTotalXS[i2]);
300 
301  static const double conv_factor = 1e-16 * cm * cm;
302 
303  sigma *= conv_factor;
304  if (sigma == 0.) sigma = 1e-30;
305  return (sigma);
306 }
static const double cm
Definition: G4SIunits.hh:118
TTree * t1
Definition: plottest35.C:26
G4double LinInterpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2)
static const double eV
Definition: G4SIunits.hh:212
TTree * t2
Definition: plottest35.C:36
Here is the call graph for this function:
Here is the caller graph for this function:

◆ VibrationEnergy()

G4double G4DNASancheExcitationModel::VibrationEnergy ( G4int  level)
private

Definition at line 310 of file G4DNASancheExcitationModel.cc.

311 {
312  static G4double energies[9] = { 0.01, 0.024, 0.061, 0.092, 0.204, 0.417, 0.460,
313  0.500, 0.835 };
314  return (energies[level] * eV);
315 }
static const double eV
Definition: G4SIunits.hh:212
double G4double
Definition: G4Types.hh:76
Here is the caller graph for this function:

Member Data Documentation

◆ fEnergyLevelXS

std::vector<std::vector<double> > G4DNASancheExcitationModel::fEnergyLevelXS
private

Definition at line 110 of file G4DNASancheExcitationModel.hh.

◆ fEnergyTotalXS

std::vector<double> G4DNASancheExcitationModel::fEnergyTotalXS
private

Definition at line 111 of file G4DNASancheExcitationModel.hh.

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4DNASancheExcitationModel::fParticleChangeForGamma
protected

Definition at line 83 of file G4DNASancheExcitationModel.hh.

◆ fpWaterDensity

const std::vector<G4double>* G4DNASancheExcitationModel::fpWaterDensity
private

Definition at line 87 of file G4DNASancheExcitationModel.hh.

◆ highEnergyLimit

G4double G4DNASancheExcitationModel::highEnergyLimit
private

Definition at line 90 of file G4DNASancheExcitationModel.hh.

◆ isInitialised

G4bool G4DNASancheExcitationModel::isInitialised
private

Definition at line 91 of file G4DNASancheExcitationModel.hh.

◆ lowEnergyLimit

G4double G4DNASancheExcitationModel::lowEnergyLimit
private

Definition at line 89 of file G4DNASancheExcitationModel.hh.

◆ nLevels

G4int G4DNASancheExcitationModel::nLevels
private

Definition at line 97 of file G4DNASancheExcitationModel.hh.

◆ tdummyVec

std::vector<double> G4DNASancheExcitationModel::tdummyVec
private

Definition at line 109 of file G4DNASancheExcitationModel.hh.

◆ verboseLevel

G4int G4DNASancheExcitationModel::verboseLevel
private

Definition at line 92 of file G4DNASancheExcitationModel.hh.


The documentation for this class was generated from the following files: