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

#include <G4DNARPWBAExcitationModel.hh>

+ Inheritance diagram for G4DNARPWBAExcitationModel:

Public Member Functions

 G4DNARPWBAExcitationModel (const G4ParticleDefinition *p=nullptr, const G4String &nam="DNARPWBAExcitationModel")
 
 ~G4DNARPWBAExcitationModel () override
 
G4DNARPWBAExcitationModeloperator= (const G4DNARPWBAExcitationModel &right)=delete
 
 G4DNARPWBAExcitationModel (const G4DNARPWBAExcitationModel &)=delete
 
void Initialise (const G4ParticleDefinition *, const G4DataVector &= *(new G4DataVector())) override
 
G4double CrossSectionPerVolume (const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax) override
 
G4double GetPartialCrossSection (const G4Material *, G4int level, const G4ParticleDefinition *, G4double kineticEnergy) override
 
void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
 
void SelectStationary (const 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 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 = nullptr
 
- 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 52 of file G4DNARPWBAExcitationModel.hh.

Constructor & Destructor Documentation

◆ G4DNARPWBAExcitationModel() [1/2]

G4DNARPWBAExcitationModel::G4DNARPWBAExcitationModel ( const G4ParticleDefinition * p = nullptr,
const G4String & nam = "DNARPWBAExcitationModel" )
explicit

Definition at line 48 of file G4DNARPWBAExcitationModel.cc.

50 : G4VEmModel(nam)
51{
52 // Verbosity scale:
53 // 0 = nothing
54 // 1 = warning for energy non-conservation
55 // 2 = details of energy budget
56 // 3 = calculation of cross sections, file openings, sampling of atoms
57 // 4 = entering in methods
58 if(verboseLevel > 0)
59 {
60 G4cout << "RPWBA excitation model is constructed " << G4endl;
61 }
62}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
G4VEmModel(const G4String &nam)
Definition G4VEmModel.cc:67

◆ ~G4DNARPWBAExcitationModel()

G4DNARPWBAExcitationModel::~G4DNARPWBAExcitationModel ( )
overridedefault

◆ G4DNARPWBAExcitationModel() [2/2]

G4DNARPWBAExcitationModel::G4DNARPWBAExcitationModel ( const G4DNARPWBAExcitationModel & )
delete

Member Function Documentation

◆ CrossSectionPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 131 of file G4DNARPWBAExcitationModel.cc.

134{
135 if(verboseLevel > 3)
136 {
137 G4cout << "Calling CrossSectionPerVolume() of G4DNARPWBAExcitationModel"
138 << G4endl;
139 }
140
141 if(fTableData == nullptr)
142 {
143 G4ExceptionDescription exceptionDescription;
144 exceptionDescription << "No cross section data ";
145 G4Exception("G4DNARPWBAIonisationModel::CrossSectionPerVolume", "em00120",
146 FatalException, exceptionDescription);
147 }
148
149 if(particleDefinition != fParticleDefinition)
150 return 0;
151
152 // Calculate total cross section for model
153
154 G4double sigma = 0;
155
156 G4double waterDensity = (*fpMolWaterDensity)[material->GetIndex()];
157
158 if(ekin >= fLowEnergy && ekin <= fHighEnergy)
159 {
160 sigma = fTableData->FindValue(ekin);
161 }
162
163 if(verboseLevel > 2)
164 {
165 G4cout << "__________________________________" << G4endl;
166 G4cout << "G4DNARPWBAExcitationModel - XS INFO START" << G4endl;
167 G4cout << "Kinetic energy(eV)=" << ekin / eV
168 << " particle : " << particleDefinition->GetParticleName() << G4endl;
169 G4cout << "Cross section per water molecule (cm^2)=" << sigma / cm / cm
170 << G4endl;
171 G4cout << "Cross section per water molecule (cm^-1)="
172 << sigma * waterDensity / (1. / cm) << G4endl;
173 G4cout << "G4DNARPWBAExcitationModel - XS INFO END" << G4endl;
174 }
175
176 return sigma * waterDensity;
177}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
double G4double
Definition G4Types.hh:83
std::size_t GetIndex() const

◆ GetPartialCrossSection()

G4double G4DNARPWBAExcitationModel::GetPartialCrossSection ( const G4Material * ,
G4int level,
const G4ParticleDefinition * particle,
G4double kineticEnergy )
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 219 of file G4DNARPWBAExcitationModel.cc.

222{
223 if(fParticleDefinition != particle)
224 {
225 G4Exception("G4DNARPWBAExcitationModel::GetPartialCrossSection",
226 "RPWBAParticleType", FatalException,
227 "Model initialized for another particle type.");
228 }
229
230 return fTableData->GetComponent(level)->FindValue(kineticEnergy);
231}

◆ Initialise()

void G4DNARPWBAExcitationModel::Initialise ( const G4ParticleDefinition * particle,
const G4DataVector & = *(new G4DataVector()) )
overridevirtual

Implements G4VEmModel.

Definition at line 70 of file G4DNARPWBAExcitationModel.cc.

72{
73 if(isInitialised)
74 {
75 return;
76 }
77 if(verboseLevel > 3)
78 {
79 G4cout << "Calling G4DNARPWBAExcitationModel::Initialise()" << G4endl;
80 }
81
82 if(fParticleDefinition != nullptr && fParticleDefinition != particle)
83 {
84 G4Exception("G4DNARPWBAExcitationModel::Initialise", "em0001",
86 "Model already initialized for another particle type.");
87 }
88
89 fTableFile = "dna/sigma_excitation_p_RPWBA";
90 fLowEnergy = 100. * MeV;
91 fHighEnergy = 300. * MeV;
92
93 if(LowEnergyLimit() < fLowEnergy || HighEnergyLimit() > fHighEnergy)
94 {
96 ed << "Model is applicable from "<<fLowEnergy<<" to "<<fHighEnergy;
97 G4Exception("G4DNARPWBAExcitationModel::Initialise", "em0004",
98 FatalException, ed);
99 }
100
101 G4double scaleFactor = 1 * cm * cm;
102 fTableData = make_unique<G4DNACrossSectionDataSet>(new G4LogLogInterpolation,
103 eV, scaleFactor);
104 fTableData->LoadData(fTableFile);
105
106 if(verboseLevel > 0)
107 {
108 G4cout << "RPWBA excitation model is initialized " << G4endl
109 << "Energy range: " << LowEnergyLimit() / eV << " eV - "
110 << HighEnergyLimit() / keV << " keV for "
111 << particle->GetParticleName() << G4endl;
112 }
113
114 // Initialize water density pointer
115 if(G4Material::GetMaterial("G4_WATER") != nullptr){
116 fpMolWaterDensity =
118 G4Material::GetMaterial("G4_WATER"));
119 }else{
120 G4ExceptionDescription exceptionDescription;
121 exceptionDescription << "G4_WATER does not exist :";
122 G4Exception("G4DNARPWBAIonisationModel::Initialise", "em00020",
123 FatalException, exceptionDescription);
124 }
126 isInitialised = true;
127}
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()
G4ParticleChangeForGamma * fParticleChangeForGamma
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
const G4String & GetParticleName() const
G4ParticleChangeForGamma * GetParticleChangeForGamma()
G4double LowEnergyLimit() const
G4double HighEnergyLimit() const

◆ operator=()

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

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 181 of file G4DNARPWBAExcitationModel.cc.

185{
186 if(verboseLevel > 3)
187 {
188 G4cout << "Calling SampleSecondaries() of G4DNARPWBAExcitationModel"
189 << G4endl;
190 }
191
192 G4double k = aDynamicParticle->GetKineticEnergy();
193
194 G4int level = RandomSelect(k);
195 G4double excitationEnergy = waterStructure.ExcitationEnergy(level);
196 G4double newEnergy = k - excitationEnergy;
197
198 if(newEnergy > 0)
199 {
201 aDynamicParticle->GetMomentumDirection());
202
203 if(!statCode){
205 }
206 else{
208 }
210 }
211
212 const G4Track* theIncomingTrack = fParticleChangeForGamma->GetCurrentTrack();
214 eExcitedMolecule, level, theIncomingTrack);
215}
@ eExcitedMolecule
int G4int
Definition G4Types.hh:85
static G4DNAChemistryManager * Instance()
void CreateWaterMolecule(ElectronicModification, G4int, const G4Track *)
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(const G4ThreeVector &Pfinal)
const G4Track * GetCurrentTrack() const
void ProposeLocalEnergyDeposit(G4double anEnergyPart)

◆ SelectStationary()

void G4DNARPWBAExcitationModel::SelectStationary ( const G4bool & input)
inline

Definition at line 101 of file G4DNARPWBAExcitationModel.hh.

102{
103 statCode = input;
104}

Member Data Documentation

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4DNARPWBAExcitationModel::fParticleChangeForGamma = nullptr
protected

Definition at line 81 of file G4DNARPWBAExcitationModel.hh.

Referenced by Initialise(), and SampleSecondaries().


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