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

#include <G4DNAChampionElasticModel.hh>

+ Inheritance diagram for G4DNAChampionElasticModel:

Public Member Functions

 G4DNAChampionElasticModel (const G4ParticleDefinition *p=0, const G4String &nam="DNAChampionElasticModel")
 
virtual ~G4DNAChampionElasticModel ()
 
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)
 
void SetKillBelowThreshold (G4double threshold)
 
G4double GetKillBelowThreshold ()
 
- 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 40 of file G4DNAChampionElasticModel.hh.

Constructor & Destructor Documentation

◆ G4DNAChampionElasticModel()

G4DNAChampionElasticModel::G4DNAChampionElasticModel ( const G4ParticleDefinition p = 0,
const G4String nam = "DNAChampionElasticModel" 
)

Definition at line 42 of file G4DNAChampionElasticModel.cc.

43 :
44 G4VEmModel(nam), isInitialised(false)
45{
46 SetLowEnergyLimit(7.4 * eV);
47 SetHighEnergyLimit(1. * MeV);
48
49 verboseLevel = 0;
50 // Verbosity scale:
51 // 0 = nothing
52 // 1 = warning for energy non-conservation
53 // 2 = details of energy budget
54 // 3 = calculation of cross sections, file openings, sampling of atoms
55 // 4 = entering in methods
56
57#ifdef CHAMPION_VERBOSE
58 if (verboseLevel > 0)
59 {
60 G4cout << "Champion Elastic model is constructed "
61 << G4endl
62 << "Energy range: "
63 << LowEnergyLimit() / eV << " eV - "
64 << HighEnergyLimit() / MeV << " MeV"
65 << G4endl;
66 }
67#endif
68
70 fpMolWaterDensity = 0;
71 fpData = 0;
72}
#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

◆ ~G4DNAChampionElasticModel()

G4DNAChampionElasticModel::~G4DNAChampionElasticModel ( )
virtual

Definition at line 76 of file G4DNAChampionElasticModel.cc.

77{
78 // For total cross section
79 if(fpData) delete fpData;
80
81 // For final state
82 eVecm.clear();
83}

Member Function Documentation

◆ CrossSectionPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 229 of file G4DNAChampionElasticModel.cc.

239{
240#ifdef CHAMPION_VERBOSE
241 if (verboseLevel > 3)
242 {
243 G4cout << "Calling CrossSectionPerVolume() of G4DNAChampionElasticModel"
244 << G4endl;
245 }
246#endif
247
248 // Calculate total cross section for model
249
250 G4double sigma = 0.;
251 G4double waterDensity = (*fpMolWaterDensity)[material->GetIndex()];
252
253 if (ekin <= HighEnergyLimit() && ekin >= LowEnergyLimit())
254 {
255 //SI : XS must not be zero otherwise sampling of secondaries method ignored
256 //
257 sigma = fpData->FindValue(ekin);
258 }
259
260#ifdef CHAMPION_VERBOSE
261 if (verboseLevel > 2)
262 {
263 G4cout << "__________________________________" << G4endl;
264 G4cout << "=== G4DNAChampionElasticModel - XS INFO START" << G4endl;
265 G4cout << "=== Kinetic energy(eV)=" << ekin/eV << " particle : " << p->GetParticleName() << G4endl;
266 G4cout << "=== Cross section per water molecule (cm^2)=" << sigma/cm/cm << G4endl;
267 G4cout << "=== Cross section per water molecule (cm^-1)=" << sigma*waterDensity/(1./cm) << G4endl;
268 G4cout << "=== G4DNAChampionElasticModel - XS INFO END" << G4endl;
269 }
270#endif
271
272 return sigma*waterDensity;
273}
double G4double
Definition: G4Types.hh:83
virtual G4double FindValue(G4double e, G4int componentId=0) const
size_t GetIndex() const
Definition: G4Material.hh:258
const G4String & GetParticleName() const

◆ GetKillBelowThreshold()

G4double G4DNAChampionElasticModel::GetKillBelowThreshold ( )
inline

Definition at line 67 of file G4DNAChampionElasticModel.hh.

68 {
70 errMsg << "The method G4DNAChampionElasticModel::"
71 "GetKillBelowThreshold is deprecated";
72
73 G4Exception("G4DNAChampionElasticModel::GetKillBelowThreshold",
74 "deprecated",
76 errMsg);
77 return 0.;
78 }
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40

◆ Initialise()

void G4DNAChampionElasticModel::Initialise ( const G4ParticleDefinition particle,
const G4DataVector  
)
virtual

Implements G4VEmModel.

Definition at line 87 of file G4DNAChampionElasticModel.cc.

89{
90#ifdef CHAMPION_VERBOSE
91 if (verboseLevel > 3)
92 {
93 G4cout << "Calling G4DNAChampionElasticModel::Initialise()" << G4endl;
94 }
95#endif
96
97 if(particle->GetParticleName() != "e-")
98 {
99 G4Exception("G4DNAChampionElasticModel::Initialise",
100 "em0002",
102 "Model not applicable to particle type.");
103 }
104
105 // Energy limits
106
107 if (LowEnergyLimit() < 7.4*eV)
108 {
109 G4cout << "G4DNAChampionElasticModel: low energy limit increased from "
110 << LowEnergyLimit()/eV << " eV to " << 7.4 << " eV"
111 << G4endl;
112 SetLowEnergyLimit(7.4*eV);
113 }
114
115 if (HighEnergyLimit() > 1.*MeV)
116 {
117 G4cout << "G4DNAChampionElasticModel: high energy limit decreased from "
118 << HighEnergyLimit()/MeV << " MeV to " << 1. << " MeV"
119 << G4endl;
120 SetHighEnergyLimit(1.*MeV);
121 }
122
123 if (isInitialised) { return; }
124
125 // *** ELECTRON
126 // For total cross section
127 // Reading of data files
128
129 G4double scaleFactor = 1e-16*cm*cm;
130
131 G4String fileElectron("dna/sigma_elastic_e_champion");
132
134 eV,
135 scaleFactor );
136 fpData->LoadData(fileElectron);
137
138 // For final state
139
140 char *path = getenv("G4LEDATA");
141
142 if (!path)
143 {
144 G4Exception("G4ChampionElasticModel::Initialise",
145 "em0006",
147 "G4LEDATA environment variable not set.");
148 return;
149 }
150
151 std::ostringstream eFullFileName;
152 eFullFileName << path << "/dna/sigmadiff_cumulated_elastic_e_champion_hp.dat";
153 std::ifstream eDiffCrossSection(eFullFileName.str().c_str());
154
155 if (!eDiffCrossSection)
156 {
158 errMsg << "Missing data file:/dna/sigmadiff_cumulated_elastic_e_champion_hp.dat; "
159 << "please use G4EMLOW6.36 and above.";
160
161 G4Exception("G4DNAChampionElasticModel::Initialise",
162 "em0003",
164 errMsg);
165 }
166
167 // March 25th, 2014 - Vaclav Stepan, Sebastien Incerti
168 // Added clear for MT
169
170 eTdummyVec.clear();
171 eVecm.clear();
172 eDiffCrossSectionData.clear();
173
174 //
175
176 eTdummyVec.push_back(0.);
177
178 while(!eDiffCrossSection.eof())
179 {
180 G4double tDummy;
181 G4double eDummy;
182 eDiffCrossSection >> tDummy >> eDummy;
183
184 // SI : mandatory eVecm initialization
185
186 if (tDummy != eTdummyVec.back())
187 {
188 eTdummyVec.push_back(tDummy);
189 eVecm[tDummy].push_back(0.);
190 }
191
192 eDiffCrossSection >> eDiffCrossSectionData[tDummy][eDummy];
193
194 if (eDummy != eVecm[tDummy].back()) eVecm[tDummy].push_back(eDummy);
195 }
196
197 // End final state
198
199#ifdef CHAMPION_VERBOSE
200 if (verboseLevel>0)
201 {
202 if (verboseLevel > 2)
203 {
204 G4cout << "Loaded cross section files for Champion Elastic model" << G4endl;
205 }
206
207 G4cout << "Champion Elastic model is initialized " << G4endl
208 << "Energy range: "
209 << LowEnergyLimit() / eV << " eV - "
210 << HighEnergyLimit() / MeV << " MeV"
211 << G4endl;
212 }
213#endif
214
215 // Initialize water density pointer
217
218 fpMolWaterDensity = G4DNAMolecularMaterial::Instance()->
219 GetNumMolPerVolTableFor(G4Material::GetMaterial("G4_WATER"));
220
222 isInitialised = true;
223
224}
@ FatalException
virtual G4bool LoadData(const G4String &argFileName)
static G4DNAMolecularMaterial * Instance()
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:651
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:133

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 277 of file G4DNAChampionElasticModel.cc.

282{
283
284#ifdef CHAMPION_VERBOSE
285 if (verboseLevel > 3)
286 {
287 G4cout << "Calling SampleSecondaries() of G4DNAChampionElasticModel" << G4endl;
288 }
289#endif
290
291 G4double electronEnergy0 = aDynamicElectron->GetKineticEnergy();
292
293 G4double cosTheta = RandomizeCosTheta(electronEnergy0);
294
295 G4double phi = 2. * pi * G4UniformRand();
296
297 G4ThreeVector zVers = aDynamicElectron->GetMomentumDirection();
298 G4ThreeVector xVers = zVers.orthogonal();
299 G4ThreeVector yVers = zVers.cross(xVers);
300
301 G4double xDir = std::sqrt(1. - cosTheta*cosTheta);
302 G4double yDir = xDir;
303 xDir *= std::cos(phi);
304 yDir *= std::sin(phi);
305
306 G4ThreeVector zPrimeVers((xDir*xVers + yDir*yVers + cosTheta*zVers));
307
309
311
312}
#define G4UniformRand()
Definition: Randomize.hh:52
Hep3Vector orthogonal() const
Hep3Vector cross(const Hep3Vector &) const
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
const G4double pi

◆ SetKillBelowThreshold()

void G4DNAChampionElasticModel::SetKillBelowThreshold ( G4double  threshold)

Definition at line 466 of file G4DNAChampionElasticModel.cc.

467{
469 errMsg << "The method G4DNAChampionElasticModel::SetKillBelowThreshold is deprecated";
470
471 G4Exception("G4DNAChampionElasticModel::SetKillBelowThreshold",
472 "deprecated",
474 errMsg);
475}

Member Data Documentation

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4DNAChampionElasticModel::fParticleChangeForGamma
protected

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