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

#include <G4DNACPA100ExcitationModel.hh>

+ Inheritance diagram for G4DNACPA100ExcitationModel:

Public Member Functions

 G4DNACPA100ExcitationModel (const G4ParticleDefinition *p=nullptr, const G4String &nam="DNACPA100ExcitationModel")
 
 ~G4DNACPA100ExcitationModel () override=default
 
void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 Initialise Each model must implement an Initialize method.
 
G4double CrossSectionPerVolume (const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax) override
 CrossSectionPerVolume Every model must implement its own CrossSectionPerVolume method. It is used by the process to determine the step path and must return a cross section times a number of molecules per volume unit.
 
void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
 SampleSecondaries Each model must implement SampleSecondaries to decide if a particle will be created after the ModelInterface or if any charateristic of the incident particle will change.
 
void SelectStationary (G4bool input)
 
G4DNACPA100ExcitationModeloperator= (const G4DNACPA100ExcitationModel &right)=delete
 
 G4DNACPA100ExcitationModel (const G4DNACPA100ExcitationModel &)=delete
 
- Public Member Functions inherited from G4VDNAModel
 G4VDNAModel (const G4String &nam, const G4String &applyToMaterial="")
 G4VDNAModel Constructeur of the G4VDNAModel class.
 
 ~G4VDNAModel () override
 ~G4VDNAModel
 
G4bool IsMaterialDefine (const size_t &materialID)
 IsMaterialDefine Check if the given material is defined in the simulation.
 
G4bool IsParticleExistingInModelForMaterial (const G4ParticleDefinition *particleName, const size_t &materialID)
 IsParticleExistingInModelForMaterial To check two things: 1- is the material existing in model ? 2- if yes, is the particle defined for that material ?
 
G4String GetName ()
 GetName.
 
G4double GetHighELimit (const size_t &materialID, const G4ParticleDefinition *particle)
 GetHighEnergyLimit.
 
G4double GetLowELimit (const size_t &materialID, const G4ParticleDefinition *particle)
 GetLowEnergyLimit.
 
void SetHighELimit (const size_t &materialID, const G4ParticleDefinition *particle, G4double lim)
 SetHighEnergyLimit.
 
void SetLowELimit (const size_t &materialID, const G4ParticleDefinition *particle, G4double lim)
 SetLowEnergyLimit.
 
- 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
 

Public Attributes

G4int verboseLevel = 0
 

Additional Inherited Members

- Protected Types inherited from G4VDNAModel
using MaterialParticleMapData
 
- Protected Member Functions inherited from G4VDNAModel
MaterialParticleMapDataGetData ()
 GetTableData.
 
std::vector< G4StringBuildApplyToMatVect (const G4String &materials)
 BuildApplyToMatVect Build the material name vector which is used to know the materials the user want to include in the model.
 
void ReadAndSaveCSFile (const size_t &materialID, const G4ParticleDefinition *p, const G4String &file, const G4double &scaleFactor)
 ReadAndSaveCSFile Read and save a "simple" cross section file : use of G4DNACrossSectionDataSet->loadData()
 
G4int RandomSelectShell (const G4double &k, const G4ParticleDefinition *particle, const size_t &materialName)
 RandomSelectShell Method to randomely select a shell from the data table uploaded. The size of the table (number of columns) is used to determine the total number of possible shells.
 
void AddCrossSectionData (const size_t &materialName, const G4ParticleDefinition *particleName, const G4String &fileCS, const G4String &fileDiffCS, const G4double &scaleFactor)
 AddCrossSectionData Method used during the initialization of the model class to add a new material. It adds a material to the model and fills vectors with informations.
 
void AddCrossSectionData (const size_t &materialName, const G4ParticleDefinition *particleName, const G4String &fileCS, const G4double &scaleFactor)
 AddCrossSectionData Method used during the initialization of the model class to add a new material. It adds a material to the model and fills vectors with informations. Not every model needs differential cross sections.
 
void LoadCrossSectionData (const G4ParticleDefinition *particleName)
 LoadCrossSectionData Method to loop on all the registered materials in the model and load the corresponding data.
 
virtual void ReadDiffCSFile (const size_t &materialName, const G4ParticleDefinition *particleName, const G4String &path, const G4double &scaleFactor)
 ReadDiffCSFile Virtual method that need to be implemented if one wish to use the differential cross sections. The read method for that kind of information is not standardized yet.
 
void EnableForMaterialAndParticle (const size_t &materialID, const G4ParticleDefinition *p)
 EnableMaterialAndParticle.
 
- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 
- Protected Attributes inherited from G4VDNAModel
G4bool isInitialised = false
 
- 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
 

Detailed Description

Definition at line 56 of file G4DNACPA100ExcitationModel.hh.

Constructor & Destructor Documentation

◆ G4DNACPA100ExcitationModel() [1/2]

G4DNACPA100ExcitationModel::G4DNACPA100ExcitationModel ( const G4ParticleDefinition * p = nullptr,
const G4String & nam = "DNACPA100ExcitationModel" )
explicit

Definition at line 57 of file G4DNACPA100ExcitationModel.cc.

59 : G4VDNAModel(nam, "all")
60{
61 fpGuanine = G4Material::GetMaterial("G4_GUANINE", false);
62 fpG4_WATER = G4Material::GetMaterial("G4_WATER", false);
63 fpDeoxyribose = G4Material::GetMaterial("G4_DEOXYRIBOSE", false);
64 fpCytosine = G4Material::GetMaterial("G4_CYTOSINE", false);
65 fpThymine = G4Material::GetMaterial("G4_THYMINE", false);
66 fpAdenine = G4Material::GetMaterial("G4_ADENINE", false);
67 fpPhosphate = G4Material::GetMaterial("G4_PHOSPHORIC_ACID", false);
68 fpParticle = G4Electron::ElectronDefinition();
69}
static G4Electron * ElectronDefinition()
Definition G4Electron.cc:86
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
G4VDNAModel(const G4String &nam, const G4String &applyToMaterial="")
G4VDNAModel Constructeur of the G4VDNAModel class.

◆ ~G4DNACPA100ExcitationModel()

G4DNACPA100ExcitationModel::~G4DNACPA100ExcitationModel ( )
overridedefault

◆ G4DNACPA100ExcitationModel() [2/2]

G4DNACPA100ExcitationModel::G4DNACPA100ExcitationModel ( const G4DNACPA100ExcitationModel & )
delete

Member Function Documentation

◆ CrossSectionPerVolume()

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

CrossSectionPerVolume Every model must implement its own CrossSectionPerVolume method. It is used by the process to determine the step path and must return a cross section times a number of molecules per volume unit.

Parameters
material

param materialName

Parameters
p

param ekin

Parameters
emin
emax
Returns
crossSection*numberOfMoleculesPerVolumeUnit

Implements G4VDNAModel.

Definition at line 162 of file G4DNACPA100ExcitationModel.cc.

165{
166 // Get the name of the current particle
167 G4String particleName = p->GetParticleName();
168 auto MatID = material->GetIndex();
169 // initialise variables
170 G4double lowLim;
171 G4double highLim;
172 G4double sigma = 0;
173
174 // Get the low energy limit for the current particle
175 lowLim = fpModelData->GetLowELimit(MatID, p);
176
177 // Get the high energy limit for the current particle
178 highLim = fpModelData->GetHighELimit(MatID, p);
179
180 // Check that we are in the correct energy range
181 if (ekin >= lowLim && ekin < highLim) {
182 // Get the map with all the data tables
183 auto Data = fpModelData->GetData();
184
185 if ((*Data)[MatID][p] == nullptr) {
186 G4Exception("G4DNACPA100ExcitationModel::CrossSectionPerVolume", "em00236", FatalException,
187 "No model is registered");
188 }
189 // Retrieve the cross section value
190 sigma = (*Data)[MatID][p]->FindValue(ekin);
191
192 if (verboseLevel > 2) {
193 auto MolDensity =
195 G4cout << "__________________________________" << G4endl;
196 G4cout << "°°° G4DNACPA100ExcitationModel - XS INFO START" << G4endl;
197 G4cout << "°°° Kinetic energy(eV)=" << ekin / eV << " particle : " << particleName << G4endl;
198 G4cout << "°°° lowLim (eV) = " << lowLim / eV << " highLim (eV) : " << highLim / eV << G4endl;
199 G4cout << "°°° Materials = " << (*G4Material::GetMaterialTable())[MatID]->GetName() << G4endl;
200 G4cout << "°°° Cross section per " << MatID << " ID molecule (cm^2)=" << sigma / cm / cm
201 << G4endl;
202 G4cout << "°°° Cross section per Phosphate molecule (cm^-1)="
203 << sigma * MolDensity / (1. / cm) << G4endl;
204 G4cout << "°°° G4DNACPA100ExcitationModel - XS INFO END" << G4endl;
205 }
206 }
207
208 // Return the cross section value
209 auto MolDensity = (*G4DNAMolecularMaterial::Instance()->GetNumMolPerVolTableFor(material))[MatID];
210 return sigma * MolDensity;
211}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
double G4double
Definition G4Types.hh:83
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
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()
std::size_t GetIndex() const
static G4MaterialTable * GetMaterialTable()
const G4String & GetParticleName() const
G4String GetName()
GetName.
G4double GetHighELimit(const size_t &materialID, const G4ParticleDefinition *particle)
GetHighEnergyLimit.
G4double GetLowELimit(const size_t &materialID, const G4ParticleDefinition *particle)
GetLowEnergyLimit.
MaterialParticleMapData * GetData()
GetTableData.

◆ Initialise()

void G4DNACPA100ExcitationModel::Initialise ( const G4ParticleDefinition * particle,
const G4DataVector & cuts )
overridevirtual

Initialise Each model must implement an Initialize method.

Parameters
particle
cuts

Implements G4VDNAModel.

Definition at line 73 of file G4DNACPA100ExcitationModel.cc.

75{
76 if (isInitialised) {
77 return;
78 }
79 if (verboseLevel > 3) {
80 G4cout << "Calling G4DNACPA100ExcitationModel::Initialise()" << G4endl;
81 }
82
84 if (p != fpParticle) {
85 std::ostringstream oss;
86 oss << " Model is not applied for this particle " << p->GetParticleName();
87 G4Exception("G4DNACPA100ExcitationModel::G4DNACPA100ExcitationModel", "CPA001",
88 FatalException, oss.str().c_str());
89 }
90
91 const char* path = G4FindDataDir("G4LEDATA");
92
93 if (path == nullptr) {
94 G4Exception("G4DNACPA100ExcitationModel::Initialise", "em0006", FatalException,
95 "G4LEDATA environment variable not set.");
96 return;
97 }
98
99 std::size_t index;
100 if (fpG4_WATER != nullptr) {
101 index = fpG4_WATER->GetIndex();
102 AddCrossSectionData(index, p, "dna/sigma_excitation_e_cpa100", 1.e-20 * m * m);
103 SetLowELimit(index, p, 11 * eV);
104 SetHighELimit(index, p, 255955 * eV);
105 }
106 if (fpGuanine != nullptr) {
107 index = fpGuanine->GetIndex();
108 AddCrossSectionData(index, p, "dna/sigma_excitation_e_cpa100_guanine", 1. * cm * cm);
109 SetLowELimit(index, p, 11 * eV);
110 SetHighELimit(index, p, 1 * MeV);
111 }
112 if (fpDeoxyribose != nullptr) {
113 index = fpDeoxyribose->GetIndex();
114 AddCrossSectionData(index, p, "dna/sigma_excitation_e_cpa100_deoxyribose", 1. * cm * cm);
115 SetLowELimit(index, p, 11 * eV);
116 SetHighELimit(index, p, 1 * MeV);
117 }
118 if (fpCytosine != nullptr) {
119 index = fpCytosine->GetIndex();
120 AddCrossSectionData(index, p, "dna/sigma_excitation_e_cpa100_cytosine", 1. * cm * cm);
121 SetLowELimit(index, p, 11 * eV);
122 SetHighELimit(index, p, 1 * MeV);
123 }
124 if (fpThymine != nullptr) {
125 index = fpThymine->GetIndex();
126 AddCrossSectionData(index, p, "dna/sigma_excitation_e_cpa100_thymine", 1. * cm * cm);
127 SetLowELimit(index, p, 11 * eV);
128 SetHighELimit(index, p, 1 * MeV);
129 }
130 if (fpAdenine != nullptr) {
131 index = fpAdenine->GetIndex();
132 AddCrossSectionData(index, p, "dna/sigma_excitation_e_cpa100_adenine", 1. * cm * cm);
133 SetLowELimit(index, p, 11 * eV);
134 SetHighELimit(index, p, 1 * MeV);
135 }
136 if (fpPhosphate != nullptr) {
137 index = fpPhosphate->GetIndex();
138 AddCrossSectionData(index, p, "dna/sigma_excitation_e_cpa100_phosphoric_acid", 1. * cm * cm);
139 SetLowELimit(index, p, 11 * eV);
140 SetHighELimit(index, p, 1 * MeV);
141 }
142
145 fpModelData = this;
146 }
147 else {
148 auto dataModel = dynamic_cast<G4DNACPA100ExcitationModel*>(
150 if (dataModel == nullptr) {
151 G4cout << "G4DNACPA100ExcitationModel::CrossSectionPerVolume:: not good modelData" << G4endl;
152 throw;
153 }
154 fpModelData = dataModel;
155 }
156 fParticleChangeForGamma = GetParticleChangeForGamma();
157 isInitialised = true;
158}
const char * G4FindDataDir(const char *)
static G4DNAMaterialManager * Instance()
void SetMasterDataModel(const DNAModelType &t, G4VEmModel *m)
G4VEmModel * GetModel(const DNAModelType &t)
void LoadCrossSectionData(const G4ParticleDefinition *particleName)
LoadCrossSectionData Method to loop on all the registered materials in the model and load the corresp...
void SetLowELimit(const size_t &materialID, const G4ParticleDefinition *particle, G4double lim)
SetLowEnergyLimit.
void SetHighELimit(const size_t &materialID, const G4ParticleDefinition *particle, G4double lim)
SetHighEnergyLimit.
void AddCrossSectionData(const size_t &materialName, const G4ParticleDefinition *particleName, const G4String &fileCS, const G4String &fileDiffCS, const G4double &scaleFactor)
AddCrossSectionData Method used during the initialization of the model class to add a new material....
G4ParticleChangeForGamma * GetParticleChangeForGamma()
G4bool IsLocked() const

◆ operator=()

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

◆ SampleSecondaries()

void G4DNACPA100ExcitationModel::SampleSecondaries ( std::vector< G4DynamicParticle * > * ,
const G4MaterialCutsCouple * ,
const G4DynamicParticle * ,
G4double tmin,
G4double tmax )
overridevirtual

SampleSecondaries Each model must implement SampleSecondaries to decide if a particle will be created after the ModelInterface or if any charateristic of the incident particle will change.

Parameters
materialName

param particleChangeForGamma

Parameters
tmin

param tmax

Implements G4VDNAModel.

Definition at line 215 of file G4DNACPA100ExcitationModel.cc.

219{
220 auto materialID = couple->GetMaterial()->GetIndex();
221 G4double k = aDynamicParticle->GetKineticEnergy();
222 const auto& particle = aDynamicParticle->GetDefinition();
223 G4double lowLim = fpModelData->GetLowELimit(materialID, particle);
224 G4double highLim = fpModelData->GetHighELimit(materialID, particle);
225
226 // Check if we are in the correct energy range
227 if (k >= lowLim && k < highLim) {
228 G4int level;
229 G4double excitationEnergy;
230 G4double newEnergy;
231 if (materialID == fpG4_WATER->GetIndex()) {
232 level = fpModelData->RandomSelectShell(k, particle, materialID);
233 excitationEnergy = eStructure.ExcitationEnergy(level, materialID);
234 }
235 else {
236 do {
237 level = eStructure.NumberOfLevels(materialID) * G4UniformRand();
238 excitationEnergy = eStructure.ExcitationEnergy(level, materialID);
239 } while ((k - eStructure.ExcitationEnergy(level, materialID)) < 0);
240 }
241 newEnergy = k - excitationEnergy;
242
243 if (k - newEnergy <= 0) {
244 G4cout << "k : " << k << " newEnergy : " << newEnergy << G4endl;
245 G4cout << "newEnergy : " << newEnergy << " k : " << k
246 << " excitationEnergy: " << excitationEnergy << G4endl;
247 G4cout << "G4DNACPA100ExcitationModel::level : " << eStructure.NumberOfLevels(materialID)
248 << " excitationEnergy : " << excitationEnergy << G4endl;
249 G4cout << "°°° Materials = " << (*G4Material::GetMaterialTable())[materialID]->GetName()
250 << G4endl;
251 G4cout << "Attention an error occured !!!" << G4endl;
252 abort();
253 }
254
255 if (newEnergy >= 0) {
256 // We take into account direction change as described page 87 (II.92) in thesis by S. Edel
257
258 G4double cosTheta =
259 (excitationEnergy / k) / (1. + (k / (2 * electron_mass_c2)) * (1. - excitationEnergy / k));
260
261 cosTheta = std::sqrt(1. - cosTheta);
262 G4double phi = 2. * pi * G4UniformRand();
263 const G4ThreeVector& zVers = aDynamicParticle->GetMomentumDirection();
264 // Computation of scattering angles (from Subroutine DIRAN in CPA100)
265
266 G4double CT1, ST1, CF1, SF1, CT2, ST2, CF2, SF2;
267 G4double sinTheta = std::sqrt(1 - cosTheta * cosTheta);
268 CT1 = zVers.z();
269 ST1 = std::sqrt(1. - CT1 * CT1);
270
271 ST1 != 0 ? CF1 = zVers.x() / ST1 : CF1 = std::cos(2. * pi * G4UniformRand());
272 ST1 != 0 ? SF1 = zVers.y() / ST1 : SF1 = std::sqrt(1. - CF1 * CF1);
273 G4double A3, A4, A5, A2, A1;
274 A3 = sinTheta * std::cos(phi);
275 A4 = A3 * CT1 + ST1 * cosTheta;
276 A5 = sinTheta * std::sin(phi);
277 A2 = A4 * SF1 + A5 * CF1;
278 A1 = A4 * CF1 - A5 * SF1;
279
280 CT2 = CT1 * cosTheta - ST1 * A3;
281 ST2 = std::sqrt(1. - CT2 * CT2);
282
283 if (ST2 == 0) {
284 ST2 = 1E-6;
285 }
286 CF2 = A1 / ST2;
287 SF2 = A2 / ST2;
288
289 G4ThreeVector zPrimeVers(ST2 * CF2, ST2 * SF2, CT2);
290 fParticleChangeForGamma->ProposeMomentumDirection(zPrimeVers.unit());
291 if (!statCode) {
292 fParticleChangeForGamma->SetProposedKineticEnergy(newEnergy);
293 }
294 else {
295 fParticleChangeForGamma->SetProposedKineticEnergy(k);
296 }
297
298 fParticleChangeForGamma->ProposeLocalEnergyDeposit(excitationEnergy);
299
300 // Chemistry only for water;
301 if (materialID == fpG4_WATER->GetIndex()) {
302 const G4Track* theIncomingTrack = fParticleChangeForGamma->GetCurrentTrack();
304 theIncomingTrack);
305 }
306 }
307 else {
308 G4cerr << "newEnergy : " << newEnergy << " k : " << k
309 << " excitationEnergy: " << excitationEnergy << G4endl;
310 G4cerr << "G4DNACPA100ExcitationModel::level : " << eStructure.NumberOfLevels(materialID)
311 << " excitationEnergy : " << excitationEnergy << G4endl;
312 G4cerr << "°°° Materials = " << (*G4Material::GetMaterialTable())[materialID]->GetName()
313 << G4endl;
314 G4cerr << "Attention an error occured !!!" << G4endl;
315 G4Exception("G4DNACPA100ExcitationModel::SampleSecondaries", "em00236", FatalException,
316 "model is not registered for this energy");
317 }
318 }
319 else {
320 G4cerr << "k : " << k << " lowLim : " << lowLim << " highLim : " << highLim << G4endl;
321 G4Exception("G4DNACPA100ExcitationModel::SampleSecondaries", "em00236", FatalException,
322 "model is not registered for this energy");
323 }
324}
@ eExcitedMolecule
int G4int
Definition G4Types.hh:85
G4GLOB_DLL std::ostream G4cerr
#define G4UniformRand()
Definition Randomize.hh:52
double z() const
double x() const
double y() const
std::size_t NumberOfLevels(const std::size_t &MatID)
G4double ExcitationEnergy(const std::size_t &level, const std::size_t &MatID)
static G4DNAChemistryManager * Instance()
void CreateWaterMolecule(ElectronicModification, G4int, const G4Track *)
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(const G4ThreeVector &Pfinal)
G4int RandomSelectShell(const G4double &k, const G4ParticleDefinition *particle, const size_t &materialName)
RandomSelectShell Method to randomely select a shell from the data table uploaded....
const G4Track * GetCurrentTrack() const
void ProposeLocalEnergyDeposit(G4double anEnergyPart)

◆ SelectStationary()

void G4DNACPA100ExcitationModel::SelectStationary ( G4bool input)
inline

Definition at line 97 of file G4DNACPA100ExcitationModel.hh.

98{
99 statCode = input;
100}

Member Data Documentation

◆ verboseLevel

G4int G4DNACPA100ExcitationModel::verboseLevel = 0

Definition at line 77 of file G4DNACPA100ExcitationModel.hh.

Referenced by CrossSectionPerVolume(), and Initialise().


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