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

#include <G4DNABornExcitationModel2.hh>

+ Inheritance diagram for G4DNABornExcitationModel2:

Public Member Functions

 G4DNABornExcitationModel2 (const G4ParticleDefinition *p=0, const G4String &nam="DNABornExcitationModel")
 
virtual ~G4DNABornExcitationModel2 ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &= *(new G4DataVector()))
 
virtual G4double CrossSectionPerVolume (const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
 
virtual G4double GetPartialCrossSection (const G4Material *, G4int level, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
 
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 44 of file G4DNABornExcitationModel2.hh.

Constructor & Destructor Documentation

◆ G4DNABornExcitationModel2()

G4DNABornExcitationModel2::G4DNABornExcitationModel2 ( const G4ParticleDefinition p = 0,
const G4String nam = "DNABornExcitationModel" 
)

Definition at line 43 of file G4DNABornExcitationModel2.cc.

44 :
45G4VEmModel(nam), isInitialised(false), fTableData(0)
46{
47 fpMolWaterDensity = 0;
48 fHighEnergy = 0;
49 fLowEnergy = 0;
50 fParticleDefinition = 0;
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 if (verboseLevel > 0)
61 {
62 G4cout << "Born excitation model is constructed " << G4endl;
63 }
65 fLastBinCallForFinalXS = 0;
66 fTotalXS = 0;
67 fTableData = 0;
68
69 // Selection of stationary mode
70
71 statCode = false;
72}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4ParticleChangeForGamma * fParticleChangeForGamma

◆ ~G4DNABornExcitationModel2()

G4DNABornExcitationModel2::~G4DNABornExcitationModel2 ( )
virtual

Definition at line 76 of file G4DNABornExcitationModel2.cc.

77{
78 // Cross section
79 if (fTableData)
80 delete fTableData;
81}

Member Function Documentation

◆ CrossSectionPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 187 of file G4DNABornExcitationModel2.cc.

192{
193 if (verboseLevel > 3)
194 {
195 G4cout << "Calling CrossSectionPerVolume() of G4DNABornExcitationModel2"
196 << G4endl;
197 }
198
199 if(particleDefinition != fParticleDefinition) return 0;
200
201 // Calculate total cross section for model
202
203 G4double sigma=0;
204
205 G4double waterDensity = (*fpMolWaterDensity)[material->GetIndex()];
206
207 if (ekin >= fLowEnergy && ekin <= fHighEnergy)
208 {
209 sigma = fTotalXS->Value(ekin, fLastBinCallForFinalXS);
210
211 // for(size_t i = 0; i < 5; ++i)
212 // sigma += (*fTableData)[i]->Value(ekin);
213
214 if(sigma == 0)
215 {
216 G4cerr << "PROBLEM SIGMA = 0 at " << G4BestUnit(ekin, "Energy")<< G4endl;
217 }
218 }
219
220 if (verboseLevel > 2)
221 {
222 G4cout << "__________________________________" << G4endl;
223 G4cout << "G4DNABornExcitationModel2 - XS INFO START" << G4endl;
224 G4cout << "Kinetic energy(eV)=" << ekin/eV << " particle : " << particleDefinition->GetParticleName() << G4endl;
225 G4cout << "Cross section per water molecule (cm^2)=" << sigma/cm/cm << G4endl;
226 G4cout << "Cross section per water molecule (cm^-1)=" << sigma*waterDensity/(1./cm) << G4endl;
227 G4cout << "G4DNABornExcitationModel2 - XS INFO END" << G4endl;
228 }
229
230 return sigma*waterDensity;
231}
double G4double
Definition: G4Types.hh:83
G4GLOB_DLL std::ostream G4cerr
size_t GetIndex() const
Definition: G4Material.hh:258
G4double Value(G4double theEnergy, std::size_t &lastidx) const

◆ GetPartialCrossSection()

G4double G4DNABornExcitationModel2::GetPartialCrossSection ( const G4Material ,
G4int  level,
const G4ParticleDefinition particle,
G4double  kineticEnergy 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 272 of file G4DNABornExcitationModel2.cc.

276{
277 if (fParticleDefinition != particle)
278 {
279 G4Exception("G4DNABornExcitationModel2::GetPartialCrossSection",
280 "bornParticleType",
282 "Model initialized for another particle type.");
283 }
284
285 return (*fTableData)(level)->Value(kineticEnergy);
286}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
virtual G4double Value(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
Definition: G4VEmModel.cc:415

◆ Initialise()

void G4DNABornExcitationModel2::Initialise ( const G4ParticleDefinition particle,
const G4DataVector = *(new G4DataVector()) 
)
virtual

Implements G4VEmModel.

Definition at line 85 of file G4DNABornExcitationModel2.cc.

87{
88
89 if (verboseLevel > 3)
90 {
91 G4cout << "Calling G4DNABornExcitationModel2::Initialise()" << G4endl;
92 }
93
94 if(fParticleDefinition != 0 && fParticleDefinition != particle)
95 {
96 G4Exception("G4DNABornExcitationModel2::Initialise","em0001",
97 FatalException,"Model already initialized for another particle type.");
98 }
99
100 fParticleDefinition = particle;
101
102 std::ostringstream fullFileName;
103 char *path = getenv("G4LEDATA");
104
105 if(G4String(path) == "")
106 {
107 G4Exception("G4DNABornExcitationModel2::Initialise","G4LEDATA-CHECK",
108 FatalException, "G4LEDATA not defined in environment variables");
109 }
110
111 fullFileName << path;
112
113 if(particle->GetParticleName() == "e-")
114 {
115 fullFileName << "/dna/bornExcitation-e.dat";
116 fLowEnergy = 9*eV;
117 fHighEnergy = 1*MeV;
118 }
119 else if(particle->GetParticleName() == "proton")
120 {
121 fullFileName << "/dna/bornExcitation-p.dat";
122 fLowEnergy = 500. * keV;
123 fHighEnergy = 100. * MeV;
124 }
125
126 SetLowEnergyLimit(fLowEnergy);
127 SetHighEnergyLimit(fHighEnergy);
128
129 //G4double scaleFactor = (1.e-22 / 3.343) * m*m;
130
131 fTableData = new G4PhysicsTable();
132 fTableData->RetrievePhysicsTable(fullFileName.str().c_str(), true);
133 for(size_t level = 0; level<fTableData->size(); ++level)
134 {
135 //(*fTableData)(level)->ScaleVector(1,scaleFactor);
136 (*fTableData)(level)->SetSpline(true);
137 }
138
139 size_t finalBin_i = 2000;
140 G4double E_min = fLowEnergy;
141 G4double E_max = fHighEnergy;
142 fTotalXS = new G4PhysicsLogVector(E_min, E_max, finalBin_i);
143 fTotalXS->SetSpline(true);
145 G4double finalXS;
146
147 for(size_t energy_i = 0; energy_i < finalBin_i; ++energy_i)
148 {
149 energy = fTotalXS->Energy(energy_i);
150 finalXS = 0;
151
152 for(size_t level = 0; level<fTableData->size(); ++level)
153 {
154 finalXS += (*fTableData)(level)->Value(energy);
155 }
156 fTotalXS->PutValue(energy_i, finalXS);
157 //G4cout << "energy = " << energy << " " << fTotalXS->Value(energy)
158 // << " " << energy_i << " " << finalXS << G4endl;
159 }
160
161 // for(energy = LowEnergyLimit() ; energy < HighEnergyLimit() ; energy += 1*pow(10,log10(energy)))
162 // {
163 // G4cout << "energy = " << energy << " " << fTotalXS->Value(energy) << G4endl;
164 // }
165
166 if( verboseLevel>0 )
167 {
168 G4cout << "Born excitation model is initialized " << G4endl
169 << "Energy range: "
170 << LowEnergyLimit() / eV << " eV - "
171 << HighEnergyLimit() / keV << " keV for "
172 << particle->GetParticleName()
173 << G4endl;
174 }
175
176 // Initialize water density pointer
178
179 if (isInitialised)
180 { return;}
182 isInitialised = true;
183}
const std::vector< G4double > * GetNumMolPerVolTableFor(const G4Material *) const
Retrieve a table of molecular densities (number of molecules per unit volume) in the G4 unit system f...
static G4DNAMolecularMaterial * Instance()
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:651
const G4String & GetParticleName() const
G4bool RetrievePhysicsTable(const G4String &filename, G4bool ascii=false)
G4double Energy(std::size_t index) const
void PutValue(std::size_t index, G4double theValue)
void SetSpline(G4bool)
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:757
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:133
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:652
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:645
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:764
G4double energy(const ThreeVector &p, const G4double m)

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 235 of file G4DNABornExcitationModel2.cc.

240{
241
242 if (verboseLevel > 3)
243 {
244 G4cout << "Calling SampleSecondaries() of G4DNABornExcitationModel2"
245 << G4endl;
246 }
247
248 G4double k = aDynamicParticle->GetKineticEnergy();
249
250 G4int level = RandomSelect(k);
251 G4double excitationEnergy = waterStructure.ExcitationEnergy(level);
252 G4double newEnergy = k - excitationEnergy;
253
254 if (newEnergy > 0)
255 {
257
258 if (!statCode) fParticleChangeForGamma->SetProposedKineticEnergy(newEnergy);
260
262 }
263
264 const G4Track * theIncomingTrack = fParticleChangeForGamma->GetCurrentTrack();
266 level,
267 theIncomingTrack);
268}
@ eExcitedMolecule
int G4int
Definition: G4Types.hh:85
static G4DNAChemistryManager * Instance()
void CreateWaterMolecule(ElectronicModification, G4int, const G4Track *)
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
const G4Track * GetCurrentTrack() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)

◆ SelectStationary()

void G4DNABornExcitationModel2::SelectStationary ( G4bool  input)
inline

Definition at line 108 of file G4DNABornExcitationModel2.hh.

109{
110 statCode = input;
111}

Member Data Documentation

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4DNABornExcitationModel2::fParticleChangeForGamma
protected

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