Geant4 9.6.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)
 
void ExtendLowEnergyLimit (G4double)
 
void SetVerboseLevel (int verbose)
 
- 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 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 ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., 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 *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
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)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=0)
 
void SetCrossSectionTable (G4PhysicsTable *)
 
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
 
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 ForceBuildTable (G4bool val)
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
const G4ElementGetCurrentElement () const
 

Protected Attributes

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

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 45 of file G4DNASancheExcitationModel.hh.

Constructor & Destructor Documentation

◆ G4DNASancheExcitationModel()

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

Definition at line 41 of file G4DNASancheExcitationModel.cc.

43 :G4VEmModel(nam),isInitialised(false)
44{
45 // nistwater = G4NistManager::Instance()->FindOrBuildMaterial("G4_WATER");
46 fpWaterDensity = 0;
47
48 lowEnergyLimit = 2 * eV;
49 highEnergyLimit = 100 * eV;
50 SetLowEnergyLimit(lowEnergyLimit);
51 SetHighEnergyLimit(highEnergyLimit);
52 nLevels = 9;
53
54 verboseLevel= 0;
55 // Verbosity scale:
56 // 0 = nothing
57 // 1 = warning for energy non-conservation
58 // 2 = details of energy budget
59 // 3 = calculation of cross sections, file openings, sampling of atoms
60 // 4 = entering in methods
61
62 if (verboseLevel > 0)
63 {
64 G4cout << "Sanche Excitation model is constructed " << G4endl
65 << "Energy range: "
66 << lowEnergyLimit / eV << " eV - "
67 << highEnergyLimit / eV << " eV"
68 << G4endl;
69 }
71 fpWaterDensity = 0;
72}
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
G4ParticleChangeForGamma * fParticleChangeForGamma
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:585
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:592

◆ ~G4DNASancheExcitationModel()

G4DNASancheExcitationModel::~G4DNASancheExcitationModel ( )
virtual

Definition at line 76 of file G4DNASancheExcitationModel.cc.

77{}

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

151{
152 if (verboseLevel > 3)
153 G4cout << "Calling CrossSectionPerVolume() of G4DNASancheExcitationModel" << G4endl;
154
155 // Calculate total cross section for model
156
157 G4double sigma=0;
158
159 G4double waterDensity = (*fpWaterDensity)[material->GetIndex()];
160
161 if(waterDensity!= 0.0)
162 // if (material == nistwater || material->GetBaseMaterial() == nistwater)
163 {
164
165 if (particleDefinition == G4Electron::ElectronDefinition())
166 {
167 if (ekin >= lowEnergyLimit && ekin < highEnergyLimit)
168 {
169 sigma = Sum(ekin);
170 }
171 }
172
173 if (verboseLevel > 2)
174 {
175 G4cout << "__________________________________" << G4endl;
176 G4cout << "°°° G4DNASancheExcitationModel - XS INFO START" << G4endl;
177 G4cout << "°°° Kinetic energy(eV)=" << ekin/eV << " particle : " << particleDefinition->GetParticleName() << G4endl;
178 G4cout << "°°° Cross section per water molecule (cm^2)=" << sigma/cm/cm << G4endl;
179 G4cout << "°°° Cross section per water molecule (cm^-1)=" << sigma*waterDensity/(1./cm) << G4endl;
180 // G4cout << " - Cross section per water molecule (cm^-1)=" << sigma*material->GetAtomicNumDensityVector()[1]/(1./cm) << G4endl;
181 G4cout << "°°° G4DNASancheExcitationModel - XS INFO END" << G4endl;
182 }
183
184 } // if water
185
186
187 // return sigma*2*material->GetAtomicNumDensityVector()[1];
188 return sigma*2*waterDensity;
189 // see papers for factor 2 description
190
191}
double G4double
Definition: G4Types.hh:64
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:89
size_t GetIndex() const
Definition: G4Material.hh:261

◆ ExtendLowEnergyLimit()

void G4DNASancheExcitationModel::ExtendLowEnergyLimit ( G4double  threshold)
inline

Definition at line 113 of file G4DNASancheExcitationModel.hh.

114{
115 lowEnergyLimit = threshold;
116 if (lowEnergyLimit < 2*CLHEP::eV)
117 G4Exception ("*** WARNING : the G4DNASancheExcitationModel class is not validated below 2 eV !","",JustWarning,"") ;
118}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41

Referenced by G4EmDNAPhysicsChemistry::ConstructProcess().

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 81 of file G4DNASancheExcitationModel.cc.

83{
84
85
86 if (verboseLevel > 3)
87 G4cout << "Calling G4DNASancheExcitationModel::Initialise()" << G4endl;
88
89 // Energy limits
90
91 if (LowEnergyLimit() < lowEnergyLimit)
92 {
93 G4cout << "G4DNASancheExcitationModel: low energy limit increased from " <<
94 LowEnergyLimit()/eV << " eV to " << lowEnergyLimit/eV << " eV" << G4endl;
95 SetLowEnergyLimit(lowEnergyLimit);
96 }
97
98 if (HighEnergyLimit() > highEnergyLimit)
99 {
100 G4cout << "G4DNASancheExcitationModel: high energy limit decreased from " <<
101 HighEnergyLimit()/eV << " eV to " << highEnergyLimit/eV << " eV" << G4endl;
102 SetHighEnergyLimit(highEnergyLimit);
103 }
104
105 //
106
107 if (verboseLevel > 0)
108 {
109 G4cout << "Sanche Excitation model is initialized " << G4endl
110 << "Energy range: "
111 << LowEnergyLimit() / eV << " eV - "
112 << HighEnergyLimit() / eV << " eV"
113 << G4endl;
114 }
115
116 // Initialize water density pointer
118
119 if (isInitialised) { return; }
121 isInitialised = true;
122
123 char *path = getenv("G4LEDATA");
124 std::ostringstream eFullFileName;
125 eFullFileName << path << "/dna/sigma_excitationvib_e_sanche.dat";
126 std::ifstream input(eFullFileName.str().c_str());
127
128 if (!input)
129 {
130 G4Exception("G4DNASancheExcitationModel::Initialise","em0003",
131 FatalException,"Missing data file:/dna/sigma_excitationvib_e_sanche.dat");
132 }
133
134 while(!input.eof())
135 {
136 double t;
137 input>>t;
138 tdummyVec.push_back(t);
139 input>>map1[t][0]>>map1[t][1]>>map1[t][2]>>map1[t][3]>>map1[t][4]>>map1[t][5]>>map1[t][6]>>map1[t][7]>>map1[t][8];
140 //G4cout<<t<<" "<<map1[t][0]<<map1[t][1]<<map1[t][2]<<map1[t][3]<<map1[t][4]<<map1[t][5]<<map1[t][6]<<map1[t][7]<<map1[t][8]<<G4endl;
141 }
142}
@ FatalException
static G4DNAMolecularMaterial * Instance()
const std::vector< double > * GetNumMolPerVolTableFor(const G4Material *) const
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:576
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:109
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:529
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:522

◆ PartialCrossSection()

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

Definition at line 240 of file G4DNASancheExcitationModel.cc.

241{
242 std::vector<double>::iterator t2 = std::upper_bound(tdummyVec.begin(),tdummyVec.end(), t/eV);
243 std::vector<double>::iterator t1 = t2-1;
244
245 double sigma = LinInterpolate((*t1), (*t2), t/eV, map1[*t1][level], map1[*t2][level]);
246 sigma*=1e-16*cm*cm;
247 if(sigma==0.)sigma=1e-30;
248 return (sigma);
249}

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 195 of file G4DNASancheExcitationModel.cc.

200{
201
202 if (verboseLevel > 3)
203 G4cout << "Calling SampleSecondaries() of G4DNASancheExcitationModel" << G4endl;
204
205 G4double electronEnergy0 = aDynamicElectron->GetKineticEnergy();
206 G4int level = RandomSelect(electronEnergy0);
207 G4double excitationEnergy = VibrationEnergy(level); // levels go from 0 to 8
208 G4double newEnergy = electronEnergy0 - excitationEnergy;
209
210 /*
211 if (electronEnergy0 < highEnergyLimit)
212 {
213 if (newEnergy >= lowEnergyLimit)
214 {
215 fParticleChangeForGamma->ProposeMomentumDirection(aDynamicElectron->GetMomentumDirection());
216 fParticleChangeForGamma->SetProposedKineticEnergy(newEnergy);
217 fParticleChangeForGamma->ProposeLocalEnergyDeposit(excitationEnergy);
218 }
219
220 else
221 {
222 fParticleChangeForGamma->ProposeTrackStatus(fStopAndKill);
223 fParticleChangeForGamma->ProposeLocalEnergyDeposit(electronEnergy0);
224 }
225 }
226*/
227
228 if (electronEnergy0 < highEnergyLimit && newEnergy>0.)
229 {
233 }
234
235 //
236}
int G4int
Definition: G4Types.hh:66
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)

◆ SetVerboseLevel()

void G4DNASancheExcitationModel::SetVerboseLevel ( int  verbose)
inline

Definition at line 75 of file G4DNASancheExcitationModel.hh.

75{verboseLevel = verbose;}

Member Data Documentation

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4DNASancheExcitationModel::fParticleChangeForGamma
protected

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