Geant4 11.2.2
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=nullptr, const G4String &nam="DNASancheExcitationModel")
 
 ~G4DNASancheExcitationModel () override
 
G4DNASancheExcitationModeloperator= (const G4DNASancheExcitationModel &right)=delete
 
 G4DNASancheExcitationModel (const G4DNASancheExcitationModel &)=delete
 
void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
G4double CrossSectionPerVolume (const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax) override
 
void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
 
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 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 *, const G4double &length, G4double &eloss)
 
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 FillNumberOfSecondaries (G4int &numberOfTriplets, G4int &numberOfRecoil)
 
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)
 
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)
 
const G4ElementGetCurrentElement (const G4Material *mat=nullptr) const
 
G4int SelectRandomAtomNumber (const G4Material *) const
 
const G4IsotopeGetCurrentIsotope (const G4Element *elm=nullptr) const
 
G4int SelectIsotopeNumber (const G4Element *) const
 
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 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 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 *)
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 
void SetLPMFlag (G4bool)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

Protected Attributes

G4ParticleChangeForGammafParticleChangeForGamma
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData = nullptr
 
G4VParticleChangepParticleChange = nullptr
 
G4PhysicsTablexSectionTable = nullptr
 
const G4MaterialpBaseMaterial = nullptr
 
const std::vector< G4double > * theDensityFactor = nullptr
 
const std::vector< G4int > * theDensityIdx = nullptr
 
G4double inveplus
 
G4double pFactor = 1.0
 
std::size_t currentCoupleIndex = 0
 
std::size_t basedCoupleIndex = 0
 
G4bool lossFlucFlag = true
 

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() [1/2]

G4DNASancheExcitationModel::G4DNASancheExcitationModel ( const G4ParticleDefinition * p = nullptr,
const G4String & nam = "DNASancheExcitationModel" )

Definition at line 42 of file G4DNASancheExcitationModel.cc.

43 :
44 G4VEmModel(nam)
45{
46 fpWaterDensity = nullptr;
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 = nullptr;
74
75 // Selection of stationary mode
76
77 statCode = false;
78}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
G4ParticleChangeForGamma * fParticleChangeForGamma
void SetHighEnergyLimit(G4double)
G4double LowEnergyLimit() const
G4double HighEnergyLimit() const
void SetLowEnergyLimit(G4double)
G4VEmModel(const G4String &nam)
Definition G4VEmModel.cc:67

◆ ~G4DNASancheExcitationModel()

G4DNASancheExcitationModel::~G4DNASancheExcitationModel ( )
overridedefault

◆ G4DNASancheExcitationModel() [2/2]

G4DNASancheExcitationModel::G4DNASancheExcitationModel ( const G4DNASancheExcitationModel & )
delete

Member Function Documentation

◆ CrossSectionPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 183 of file G4DNASancheExcitationModel.cc.

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

◆ ExtendLowEnergyLimit()

void G4DNASancheExcitationModel::ExtendLowEnergyLimit ( G4double threshold)
inline

Definition at line 120 of file G4DNASancheExcitationModel.hh.

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

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

◆ Initialise()

void G4DNASancheExcitationModel::Initialise ( const G4ParticleDefinition * ,
const G4DataVector &  )
overridevirtual

Implements G4VEmModel.

Definition at line 88 of file G4DNASancheExcitationModel.cc.

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

◆ operator=()

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

◆ PartialCrossSection()

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

Definition at line 290 of file G4DNASancheExcitationModel.cc.

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

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 228 of file G4DNASancheExcitationModel.cc.

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

◆ SelectStationary()

void G4DNASancheExcitationModel::SelectStationary ( G4bool input)
inline

Definition at line 132 of file G4DNASancheExcitationModel.hh.

133{
134 statCode = input;
135}

◆ SetVerboseLevel()

void G4DNASancheExcitationModel::SetVerboseLevel ( G4int verbose)
inline

Definition at line 76 of file G4DNASancheExcitationModel.hh.

77 {
78 verboseLevel = verbose;
79 }

◆ TotalCrossSection()

G4double G4DNASancheExcitationModel::TotalCrossSection ( G4double t)

Definition at line 318 of file G4DNASancheExcitationModel.cc.

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

Referenced by CrossSectionPerVolume().

Member Data Documentation

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4DNASancheExcitationModel::fParticleChangeForGamma
protected

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