Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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)
 
void SelectStationary (G4bool input)
 
- 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 level, const G4ParticleDefinition *, G4double kineticEnergy)
 
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 * > *)
 
virtual 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=nullptr)
 
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 SetFluctuationFlag (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

G4ParticleChangeForGammafParticleChangeForGamma
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 
G4bool lossFlucFlag
 

Additional Inherited Members

- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
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::G4DNASancheExcitationModel ( const G4ParticleDefinition p = 0,
const G4String nam = "DNASancheExcitationModel" 
)

Definition at line 43 of file G4DNASancheExcitationModel.cc.

44  :
45  G4VEmModel(nam), isInitialised(false)
46 {
47  fpWaterDensity = 0;
48 
50  SetHighEnergyLimit(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 #ifdef SANCHE_VERBOSE
62  if (verboseLevel > 0)
63  {
64  G4cout << "Sanche Excitation model is constructed "
65  << G4endl
66  << "Energy range: "
67  << LowEnergyLimit() / eV << " eV - "
68  << HighEnergyLimit() / eV << " eV"
69  << G4endl;
70  }
71 #endif
72 
74  fpWaterDensity = 0;
75 
76  // Selection of stationary mode
77 
78  statCode = false;
79 }
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:643
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:636
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:68
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:732
G4GLOB_DLL std::ostream G4cout
G4ParticleChangeForGamma * fParticleChangeForGamma
static constexpr double eV
Definition: G4SIunits.hh:215
#define G4endl
Definition: G4ios.hh:61
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:739

Here is the call graph for this function:

G4DNASancheExcitationModel::~G4DNASancheExcitationModel ( )
virtual

Definition at line 83 of file G4DNASancheExcitationModel.cc.

84 {
85 }

Member Function Documentation

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

Reimplemented from G4VEmModel.

Definition at line 185 of file G4DNASancheExcitationModel.cc.

194 {
195 #ifdef SANCHE_VERBOSE
196  if (verboseLevel > 3)
197  {
198  G4cout << "Calling CrossSectionPerVolume() of G4DNASancheExcitationModel"
199  << G4endl;
200  }
201 #endif
202 
203  // Calculate total cross section for model
204 
205  G4double sigma = 0.;
206 
207  G4double waterDensity = (*fpWaterDensity)[material->GetIndex()];
208 
209  if(waterDensity!= 0.0)
210  {
211  if (ekin >= LowEnergyLimit() && ekin < HighEnergyLimit())
212  // for now this is necessary
213  {
214  //sigma = Sum(ekin);
215  sigma = TotalCrossSection(ekin);
216  }
217 
218 #ifdef SANCHE_VERBOSE
219  if (verboseLevel > 2)
220  {
221  G4cout << "__________________________________" << G4endl;
222  G4cout << "=== G4DNASancheExcitationModel - XS INFO START" << G4endl;
223  G4cout << "=== Kinetic energy(eV)=" << ekin/eV << " particle : " << particleDefinition->GetParticleName() << G4endl;
224  G4cout << "=== Cross section per water molecule (cm^2)=" << sigma/cm/cm << G4endl;
225  G4cout << "=== Cross section per water molecule (cm^-1)=" << sigma*waterDensity/(1./cm) << G4endl;
226  G4cout << "=== G4DNASancheExcitationModel - XS INFO END" << G4endl;
227  }
228 #endif
229  } // if water
230 
231  return sigma*2.*waterDensity;
232  // see papers for factor 2 description
233 
234 }
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:643
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:636
size_t GetIndex() const
Definition: G4Material.hh:262
G4GLOB_DLL std::ostream G4cout
static constexpr double cm
Definition: G4SIunits.hh:119
static constexpr double eV
Definition: G4SIunits.hh:215
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

void G4DNASancheExcitationModel::ExtendLowEnergyLimit ( G4double  threshold)
inline

Definition at line 123 of file G4DNASancheExcitationModel.hh.

124 {
125  if(threshold < 2 * CLHEP::eV)
126  G4Exception("*** WARNING : the G4DNASancheExcitationModel class is not "
127  "validated below 2 eV !",
128  "", JustWarning, "");
129 
130  SetLowEnergyLimit(threshold);
131 }
static constexpr double eV
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:739

Here is the call graph for this function:

Here is the caller graph for this function:

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

Implements G4VEmModel.

Definition at line 91 of file G4DNASancheExcitationModel.cc.

93 {
94 
95 #ifdef SANCHE_VERBOSE
96  if (verboseLevel > 3)
97  {
98  G4cout << "Calling G4DNASancheExcitationModel::Initialise()"
99  << G4endl;
100  }
101 #endif
102 
103  // Energy limits
104 
105  if (LowEnergyLimit() < 2.*eV)
106  {
107  G4Exception("*** WARNING : the G4DNASancheExcitationModel class is not "
108  "validated below 2 eV !",
109  "", JustWarning, "");
110  }
111 
112  if (HighEnergyLimit() > 100.*eV)
113  {
114  G4cout << "G4DNASancheExcitationModel: high energy limit decreased from " <<
115  HighEnergyLimit()/eV << " eV to " << 100. << " eV" << G4endl;
116  SetHighEnergyLimit(100.*eV);
117  }
118 
119  //
120 #ifdef SANCHE_VERBOSE
121  if (verboseLevel > 0)
122  {
123  G4cout << "Sanche Excitation model is initialized " << G4endl
124  << "Energy range: "
125  << LowEnergyLimit() / eV << " eV - "
126  << HighEnergyLimit() / eV << " eV"
127  << G4endl;
128  }
129 #endif
130 
131  // Initialize water density pointer
132  fpWaterDensity = G4DNAMolecularMaterial::Instance()->
133  GetNumMolPerVolTableFor(G4Material::GetMaterial("G4_WATER"));
134 
135  if (isInitialised) {return;} // RETURNS HERE
136 
138  isInitialised = true;
139 
140  char *path = getenv("G4LEDATA");
141  std::ostringstream eFullFileName;
142  eFullFileName << path << "/dna/sigma_excitationvib_e_sanche.dat";
143  std::ifstream input(eFullFileName.str().c_str());
144 
145  if (!input)
146  {
147  G4Exception("G4DNASancheExcitationModel::Initialise","em0003",
148  FatalException,"Missing data file:/dna/sigma_excitationvib_e_sanche.dat");
149  }
150 
151  // March 25th, 2014 - Vaclav Stepan, Sebastien Incerti
152  // Added clear for MT
153  tdummyVec.clear();
154  //
155 
156  double t;
157  double xs;
158 
159  while(!input.eof())
160  {
161  input>>t;
162  tdummyVec.push_back(t);
163 
164  fEnergyLevelXS.push_back(std::vector<double>());
165  fEnergyTotalXS.push_back(0);
166  std::vector<double>& levelXS = fEnergyLevelXS.back();
167  levelXS.reserve(9);
168 
169 // G4cout<<t;
170 
171  for(size_t i = 0 ; i < 9 ;++i)
172  {
173  input>>xs;
174  levelXS.push_back(xs);
175  fEnergyTotalXS.back() += xs;
176 // G4cout <<" " << levelXS[i];
177  }
178 
179 // G4cout << G4endl;
180  }
181 }
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:643
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:636
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:602
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:732
G4GLOB_DLL std::ostream G4cout
G4ParticleChangeForGamma * fParticleChangeForGamma
static constexpr double eV
Definition: G4SIunits.hh:215
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static G4DNAMolecularMaterial * Instance()
#define G4endl
Definition: G4ios.hh:61
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:132

Here is the call graph for this function:

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

Definition at line 300 of file G4DNASancheExcitationModel.cc.

302 {
303  std::vector<double>::iterator t2 = std::upper_bound(tdummyVec.begin(),
304  tdummyVec.end(), t / eV);
305  std::vector<double>::iterator t1 = t2 - 1;
306 
307  size_t i1 = t1 - tdummyVec.begin();
308  size_t i2 = t2 - tdummyVec.begin();
309 
310  double sigma = LinInterpolate((*t1), (*t2),
311  t / eV,
312  fEnergyLevelXS[i1][level],
313  fEnergyLevelXS[i2][level]);
314 
315  static const double conv_factor = 1e-16 * cm * cm;
316 
317  sigma *= conv_factor;
318  if (sigma == 0.) sigma = 1e-30;
319  return (sigma);
320 }
static constexpr double cm
Definition: G4SIunits.hh:119
static constexpr double eV
Definition: G4SIunits.hh:215
tuple t1
Definition: plottest35.py:33
void G4DNASancheExcitationModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  ,
const G4MaterialCutsCouple ,
const G4DynamicParticle aDynamicElectron,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Definition at line 238 of file G4DNASancheExcitationModel.cc.

244 {
245 #ifdef SANCHE_VERBOSE
246  if (verboseLevel > 3)
247  {
248  G4cout << "Calling SampleSecondaries() of G4DNASancheExcitationModel"
249  << G4endl;
250  }
251 #endif
252 
253  G4double electronEnergy0 = aDynamicElectron->GetKineticEnergy();
254  G4int level = RandomSelect(electronEnergy0);
255  G4double excitationEnergy = VibrationEnergy(level); // levels go from 0 to 8
256  G4double newEnergy = electronEnergy0 - excitationEnergy;
257 
258  /*
259  if (electronEnergy0 < highEnergyLimit)
260  {
261  if (newEnergy >= lowEnergyLimit)
262  {
263  fParticleChangeForGamma->ProposeMomentumDirection(aDynamicElectron->GetMomentumDirection());
264  fParticleChangeForGamma->SetProposedKineticEnergy(newEnergy);
265  fParticleChangeForGamma->ProposeLocalEnergyDeposit(excitationEnergy);
266  }
267 
268  else
269  {
270  fParticleChangeForGamma->ProposeTrackStatus(fStopAndKill);
271  fParticleChangeForGamma->ProposeLocalEnergyDeposit(electronEnergy0);
272  }
273  }
274  */
275 
276  if (electronEnergy0 < HighEnergyLimit() && newEnergy>0.)
277  {
278 
279  if (!statCode)
280  {
284  }
285 
286  else
287  {
291  }
292 
293  }
294 
295  //
296 }
G4double GetKineticEnergy() const
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:636
int G4int
Definition: G4Types.hh:78
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4GLOB_DLL std::ostream G4cout
const G4ThreeVector & GetMomentumDirection() const
G4ParticleChangeForGamma * fParticleChangeForGamma
void SetProposedKineticEnergy(G4double proposedKinEnergy)
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

void G4DNASancheExcitationModel::SelectStationary ( G4bool  input)
inline

Definition at line 135 of file G4DNASancheExcitationModel.hh.

136 {
137  statCode = input;
138 }
void G4DNASancheExcitationModel::SetVerboseLevel ( int  verbose)
inline

Definition at line 75 of file G4DNASancheExcitationModel.hh.

76  {
77  verboseLevel = verbose;
78  }
G4double G4DNASancheExcitationModel::TotalCrossSection ( G4double  t)

Definition at line 324 of file G4DNASancheExcitationModel.cc.

325 {
326  std::vector<double>::iterator t2 = std::upper_bound(tdummyVec.begin(),
327  tdummyVec.end(), t / eV);
328  std::vector<double>::iterator t1 = t2 - 1;
329 
330  size_t i1 = t1 - tdummyVec.begin();
331  size_t i2 = t2 - tdummyVec.begin();
332 
333  double sigma = LinInterpolate((*t1), (*t2),
334  t / eV,
335  fEnergyTotalXS[i1],
336  fEnergyTotalXS[i2]);
337 
338  static const double conv_factor = 1e-16 * cm * cm;
339 
340  sigma *= conv_factor;
341  if (sigma == 0.) sigma = 1e-30;
342  return (sigma);
343 }
static constexpr double cm
Definition: G4SIunits.hh:119
static constexpr double eV
Definition: G4SIunits.hh:215
tuple t1
Definition: plottest35.py:33

Here is the caller graph for this function:

Member Data Documentation

G4ParticleChangeForGamma* G4DNASancheExcitationModel::fParticleChangeForGamma
protected

Definition at line 84 of file G4DNASancheExcitationModel.hh.


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