Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DNASancheExcitationModel Class Reference

#include <G4DNASancheExcitationModel.hh>

+ Inheritance 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 (G4int verbose)
 
void SelectStationary (G4bool input)
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)=0
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double tmax=DBL_MAX)=0
 
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 CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=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)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectTargetAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double logKineticEnergy, 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 *)
 
G4int SelectIsotopeNumber (const G4Element *)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
G4VEmModelGetTripletModel ()
 
void SetTripletModel (G4VEmModel *)
 
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) const
 
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
 
void SetUseBaseMaterials (G4bool val)
 
G4bool UseBaseMaterials () 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)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

Protected Attributes

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

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 *)
 

Detailed Description

Definition at line 43 of file G4DNASancheExcitationModel.hh.

Constructor & Destructor Documentation

◆ G4DNASancheExcitationModel()

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

Definition at line 42 of file G4DNASancheExcitationModel.cc.

43 :
44 G4VEmModel(nam), isInitialised(false)
45{
46 fpWaterDensity = 0;
47
48 SetLowEnergyLimit(2.*eV);
49 SetHighEnergyLimit(100*eV);
50 nLevels = 9;
51
52 verboseLevel = 0;
53 // Verbosity scale:
54 // 0 = nothing
55 // 1 = warning for energy non-conservation
56 // 2 = details of energy budget
57 // 3 = calculation of cross sections, file openings, sampling of atoms
58 // 4 = entering in methods
59
60#ifdef SANCHE_VERBOSE
61 if (verboseLevel > 0)
62 {
63 G4cout << "Sanche Excitation model is constructed "
64 << G4endl
65 << "Energy range: "
66 << LowEnergyLimit() / eV << " eV - "
67 << HighEnergyLimit() / eV << " eV"
68 << G4endl;
69 }
70#endif
71
73 fpWaterDensity = 0;
74
75 // Selection of stationary mode
76
77 statCode = false;
78}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4ParticleChangeForGamma * fParticleChangeForGamma
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:757
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:652
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:645
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:764

◆ ~G4DNASancheExcitationModel()

G4DNASancheExcitationModel::~G4DNASancheExcitationModel ( )
virtual

Definition at line 82 of file G4DNASancheExcitationModel.cc.

83{
84}

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 184 of file G4DNASancheExcitationModel.cc.

193{
194#ifdef SANCHE_VERBOSE
195 if (verboseLevel > 3)
196 {
197 G4cout << "Calling CrossSectionPerVolume() of G4DNASancheExcitationModel"
198 << G4endl;
199 }
200#endif
201
202 // Calculate total cross section for model
203
204 G4double sigma = 0.;
205
206 G4double waterDensity = (*fpWaterDensity)[material->GetIndex()];
207
208 if (ekin >= LowEnergyLimit() && ekin <= HighEnergyLimit())
209 sigma = TotalCrossSection(ekin);
210
211#ifdef SANCHE_VERBOSE
212 if (verboseLevel > 2)
213 {
214 G4cout << "__________________________________" << G4endl;
215 G4cout << "=== G4DNASancheExcitationModel - XS INFO START" << G4endl;
216 G4cout << "=== Kinetic energy(eV)=" << ekin/eV << " particle : " << particleDefinition->GetParticleName() << G4endl;
217 G4cout << "=== Cross section per water molecule (cm^2)=" << sigma/cm/cm << G4endl;
218 G4cout << "=== Cross section per water molecule (cm^-1)=" << sigma*waterDensity/(1./cm) << G4endl;
219 G4cout << "=== G4DNASancheExcitationModel - XS INFO END" << G4endl;
220 }
221#endif
222
223 return sigma*2.*waterDensity;
224 // see papers for factor 2 description (liquid phase)
225}
double G4double
Definition: G4Types.hh:83
size_t GetIndex() const
Definition: G4Material.hh:258

◆ ExtendLowEnergyLimit()

void G4DNASancheExcitationModel::ExtendLowEnergyLimit ( G4double  threshold)
inline

Definition at line 121 of file G4DNASancheExcitationModel.hh.

122{
123 if(threshold < 2 * CLHEP::eV)
124 G4Exception("*** WARNING : the G4DNASancheExcitationModel class is not "
125 "validated below 2 eV !",
126 "", JustWarning, "");
127
128 SetLowEnergyLimit(threshold);
129}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35

Referenced by G4EmDNAChemistry::ConstructProcess(), G4EmDNAChemistry_option1::ConstructProcess(), G4EmDNAChemistry_option3::ConstructProcess(), and G4EmDNAChemistry_option2::ConstructProcess().

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 89 of file G4DNASancheExcitationModel.cc.

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

◆ PartialCrossSection()

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

Definition at line 291 of file G4DNASancheExcitationModel.cc.

293{
294 // Protection against out of boundary access
295 if (t/eV==tdummyVec.back()) t=t*(1.-1e-12);
296 //
297
298 std::vector<G4double>::iterator t2 = std::upper_bound(tdummyVec.begin(),
299 tdummyVec.end(), t / eV);
300 std::vector<G4double>::iterator t1 = t2 - 1;
301
302 size_t i1 = t1 - tdummyVec.begin();
303 size_t i2 = t2 - tdummyVec.begin();
304
305 G4double sigma = LinInterpolate((*t1), (*t2),
306 t / eV,
307 fEnergyLevelXS[i1][level],
308 fEnergyLevelXS[i2][level]);
309
310 static const G4double conv_factor = 1e-16 * cm * cm;
311
312 sigma *= conv_factor;
313 if (sigma == 0.) sigma = 1e-30;
314 return (sigma);
315}

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 229 of file G4DNASancheExcitationModel.cc.

235{
236#ifdef SANCHE_VERBOSE
237 if (verboseLevel > 3)
238 {
239 G4cout << "Calling SampleSecondaries() of G4DNASancheExcitationModel"
240 << G4endl;
241 }
242#endif
243
244 G4double electronEnergy0 = aDynamicElectron->GetKineticEnergy();
245 G4int level = RandomSelect(electronEnergy0);
246 G4double excitationEnergy = VibrationEnergy(level); // levels go from 0 to 8
247 G4double newEnergy = electronEnergy0 - excitationEnergy;
248
249 /*
250 if (electronEnergy0 < highEnergyLimit)
251 {
252 if (newEnergy >= lowEnergyLimit)
253 {
254 fParticleChangeForGamma->ProposeMomentumDirection(aDynamicElectron->GetMomentumDirection());
255 fParticleChangeForGamma->SetProposedKineticEnergy(newEnergy);
256 fParticleChangeForGamma->ProposeLocalEnergyDeposit(excitationEnergy);
257 }
258
259 else
260 {
261 fParticleChangeForGamma->ProposeTrackStatus(fStopAndKill);
262 fParticleChangeForGamma->ProposeLocalEnergyDeposit(electronEnergy0);
263 }
264 }
265 */
266
267 if (electronEnergy0 <= HighEnergyLimit() && newEnergy>0.)
268 {
269
270 if (!statCode)
271 {
275 }
276
277 else
278 {
282 }
283
284 }
285
286 //
287}
int G4int
Definition: G4Types.hh:85
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)

◆ SelectStationary()

void G4DNASancheExcitationModel::SelectStationary ( G4bool  input)
inline

Definition at line 133 of file G4DNASancheExcitationModel.hh.

134{
135 statCode = input;
136}

Referenced by G4EmDNAPhysics_stationary::ConstructProcess().

◆ SetVerboseLevel()

void G4DNASancheExcitationModel::SetVerboseLevel ( G4int  verbose)
inline

Definition at line 73 of file G4DNASancheExcitationModel.hh.

74 {
75 verboseLevel = verbose;
76 }

◆ TotalCrossSection()

G4double G4DNASancheExcitationModel::TotalCrossSection ( G4double  t)

Definition at line 319 of file G4DNASancheExcitationModel.cc.

320{
321 // Protection against out of boundary access
322 if (t/eV==tdummyVec.back()) t=t*(1.-1e-12);
323 //
324
325 std::vector<G4double>::iterator t2 = std::upper_bound(tdummyVec.begin(),
326 tdummyVec.end(), t / eV);
327 std::vector<G4double>::iterator t1 = t2 - 1;
328
329 size_t i1 = t1 - tdummyVec.begin();
330 size_t i2 = t2 - tdummyVec.begin();
331
332 G4double sigma = LinInterpolate((*t1), (*t2),
333 t / eV,
334 fEnergyTotalXS[i1],
335 fEnergyTotalXS[i2]);
336
337 static const G4double conv_factor = 1e-16 * cm * cm;
338
339 sigma *= conv_factor;
340 if (sigma == 0.) sigma = 1e-30;
341 return (sigma);
342}

Referenced by CrossSectionPerVolume().

Member Data Documentation

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4DNASancheExcitationModel::fParticleChangeForGamma
protected

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