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

#include <G4DNACPA100ElasticModel.hh>

+ Inheritance diagram for G4DNACPA100ElasticModel:

Public Member Functions

 G4DNACPA100ElasticModel (const G4ParticleDefinition *p=nullptr, const G4String &nam="DNACPA100ElasticModel")
 
 ~G4DNACPA100ElasticModel () 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)
 
G4DNACPA100ElasticModeloperator= (const G4DNACPA100ElasticModel &right)=delete
 
 G4DNACPA100ElasticModel (const G4DNACPA100ElasticModel &)=delete
 
G4double GetElasticLevel (const std::size_t &l)
 
- 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
 

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 55 of file G4DNACPA100ElasticModel.hh.

Constructor & Destructor Documentation

◆ G4DNACPA100ElasticModel() [1/2]

G4DNACPA100ElasticModel::G4DNACPA100ElasticModel ( const G4ParticleDefinition * p = nullptr,
const G4String & nam = "DNACPA100ElasticModel" )
explicit

Definition at line 50 of file G4DNACPA100ElasticModel.cc.

51 : G4VDNAModel(nam, "all")
52{
53 fpGuanine = G4Material::GetMaterial("G4_GUANINE", false);
54 fpG4_WATER = G4Material::GetMaterial("G4_WATER", false);
55 fpDeoxyribose = G4Material::GetMaterial("G4_DEOXYRIBOSE", false);
56 fpCytosine = G4Material::GetMaterial("G4_CYTOSINE", false);
57 fpThymine = G4Material::GetMaterial("G4_THYMINE", false);
58 fpAdenine = G4Material::GetMaterial("G4_ADENINE", false);
59 fpPhosphate = G4Material::GetMaterial("G4_PHOSPHORIC_ACID", false);
60 fpParticle = G4Electron::ElectronDefinition();
61}
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.

◆ ~G4DNACPA100ElasticModel()

G4DNACPA100ElasticModel::~G4DNACPA100ElasticModel ( )
overridedefault

◆ G4DNACPA100ElasticModel() [2/2]

G4DNACPA100ElasticModel::G4DNACPA100ElasticModel ( const G4DNACPA100ElasticModel & )
delete

Member Function Documentation

◆ CrossSectionPerVolume()

G4double G4DNACPA100ElasticModel::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 171 of file G4DNACPA100ElasticModel.cc.

174{
175 // Get the name of the current particle
176 const G4String& particleName = p->GetParticleName();
177 auto materialID = pMaterial->GetIndex();
178
179 // set killBelowEnergy value for current material
180 fKillBelowEnergy = fpModelData->GetLowELimit(materialID, p);
181
182 G4double sigma = 0.;
183
184 if (ekin < fpModelData->GetHighELimit(materialID, p)) {
185 if (ekin < fKillBelowEnergy) {
186 return DBL_MAX;
187 }
188
189 auto tableData = fpModelData->GetData();
190
191 if ((*tableData)[materialID][p] == nullptr) {
192 G4Exception("G4DNACPA100ElasticModel::CrossSectionPerVolume", "em00236", FatalException,
193 "No model is registered");
194 }
195 sigma = (*tableData)[materialID][p]->FindValue(ekin);
196 }
197
198 if (verboseLevel > 2) {
199 auto MolDensity =
201 G4cout << "__________________________________" << G4endl;
202 G4cout << "°°° G4DNACPA100ElasticModel - XS INFO START" << G4endl;
203 G4cout << "°°° Kinetic energy(eV)=" << ekin / eV << " particle : " << particleName << G4endl;
204 G4cout << "°°° lowLim (eV) = " << GetLowELimit(materialID, p) / eV
205 << " highLim (eV) : " << GetHighELimit(materialID, p) / eV << G4endl;
206 G4cout << "°°° Materials = " << (*G4Material::GetMaterialTable())[materialID]->GetName()
207 << G4endl;
208 G4cout << "°°° Cross section per molecule (cm^2)=" << sigma / cm / cm << G4endl;
209 G4cout << "°°° Cross section per Phosphate molecule (cm^-1)=" << sigma * MolDensity / (1. / cm)
210 << G4endl;
211 G4cout << "°°° G4DNACPA100ElasticModel - XS INFO END" << G4endl;
212 }
213
214 auto MolDensity =
216 return sigma * MolDensity;
217}
@ 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()
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.
#define DBL_MAX
Definition templates.hh:62

◆ GetElasticLevel()

G4double G4DNACPA100ElasticModel::GetElasticLevel ( const std::size_t & l)
inline

Definition at line 84 of file G4DNACPA100ElasticModel.hh.

85 {
86 return fLevels[l];
87 }

Referenced by SampleSecondaries().

◆ Initialise()

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

Initialise Each model must implement an Initialize method.

Parameters
particle
cuts

Implements G4VDNAModel.

Definition at line 63 of file G4DNACPA100ElasticModel.cc.

65{
66 if (isInitialised) {
67 return;
68 }
69 if (verboseLevel > 3) {
70 G4cout << "Calling G4DNACPA100ExcitationModel::Initialise()" << G4endl;
71 }
72
74 if (p != fpParticle) {
75 std::ostringstream oss;
76 oss << " Model is not applied for this particle " << p->GetParticleName();
77 G4Exception("G4DNACPA100ElasticModel::G4DNACPA100ElasticModel", "CPA001", FatalException,
78 oss.str().c_str());
79 }
80
81 const char* path = G4FindDataDir("G4LEDATA");
82
83 if (path == nullptr) {
84 G4Exception("G4DNACPA100ElasticModel::Initialise", "em0006", FatalException,
85 "G4LEDATA environment variable not set.");
86 return;
87 }
88
89 std::size_t index;
90 if (fpG4_WATER != nullptr) {
91 index = fpG4_WATER->GetIndex();
92 fLevels[index] = 1.214e-4;
93 AddCrossSectionData(index, p, "dna/sigma_elastic_e_cpa100",
94 "dna/sigmadiff_cumulated_elastic_e_cpa100", 1e-20 * m * m);
95 SetLowELimit(index, p, 11. * eV);
96 SetHighELimit(index, p, 255955. * eV);
97 }
98 if (fpGuanine != nullptr) {
99 index = fpGuanine->GetIndex();
100 fLevels[index] = 1.4504480e-05;
101 AddCrossSectionData(index, p, "dna/sigma_elastic_e_cpa100_guanine",
102 "dna/sigmadiff_cumulated_elastic_e_cpa100_guanine", 1 * cm * cm);
103 SetLowELimit(index, p, 11 * eV);
104 SetHighELimit(index, p, 1 * MeV);
105 }
106 if (fpDeoxyribose != nullptr) {
107 index = fpDeoxyribose->GetIndex();
108 fLevels[index] = 1.6343100e-05;
109 AddCrossSectionData(index, p, "dna/sigma_elastic_e_cpa100_deoxyribose",
110 "dna/sigmadiff_cumulated_elastic_e_cpa100_deoxyribose", 1 * cm * cm);
111 SetLowELimit(index, p, 11 * eV);
112 SetHighELimit(index, p, 1 * MeV);
113 }
114 if (fpCytosine != nullptr) {
115 index = fpCytosine->GetIndex();
116 fLevels[index] = 1.9729660e-05;
117 AddCrossSectionData(index, p, "dna/sigma_elastic_e_cpa100_cytosine",
118 "dna/sigmadiff_cumulated_elastic_e_cpa100_cytosine", 1 * cm * cm);
119 SetLowELimit(index, p, 11 * eV);
120 SetHighELimit(index, p, 1 * MeV);
121 }
122 if (fpThymine != nullptr) {
123 index = fpThymine->GetIndex();
124 fLevels[index] = 1.7381300e-05;
125 AddCrossSectionData(index, p, "dna/sigma_elastic_e_cpa100_thymine",
126 "dna/sigmadiff_cumulated_elastic_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 fLevels[index] = 1.6221800e-05;
133 AddCrossSectionData(index, p, "dna/sigma_elastic_e_cpa100_adenine",
134 "dna/sigmadiff_cumulated_elastic_e_cpa100_adenine", 1 * cm * cm);
135 SetLowELimit(index, p, 11 * eV);
136 SetHighELimit(index, p, 1 * MeV);
137 }
138 if (fpPhosphate != nullptr) {
139 index = fpPhosphate->GetIndex();
140 fLevels[index] = 2.2369600e-05;
141 AddCrossSectionData(index, p, "dna/sigma_elastic_e_cpa100_phosphoric_acid",
142 "dna/sigmadiff_cumulated_elastic_e_cpa100_phosphoric_acid", 1 * cm * cm);
143 SetLowELimit(index, p, 11 * eV);
144 SetHighELimit(index, p, 1 * MeV);
145 }
146
147 // Load data
150 fpModelData = this;
151 }
152 else {
153 auto dataModel = dynamic_cast<G4DNACPA100ElasticModel*>(
155 if (dataModel == nullptr) {
156 G4cout << "G4DNACPA100ElasticModel::CrossSectionPerVolume:: not good modelData" << G4endl;
157 G4Exception("G4DNACPA100ElasticModel::CrossSectionPerVolume", "em004", FatalException,
158 "no modelData is registered");
159 }
160 else {
161 fpModelData = dataModel;
162 }
163 }
164
165 fParticleChangeForGamma = GetParticleChangeForGamma();
166 isInitialised = true;
167}
const char * G4FindDataDir(const char *)
static G4DNAMaterialManager * Instance()
void SetMasterDataModel(const DNAModelType &t, G4VEmModel *m)
G4VEmModel * GetModel(const DNAModelType &t)
std::size_t GetIndex() const
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=()

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

◆ SampleSecondaries()

void G4DNACPA100ElasticModel::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 221 of file G4DNACPA100ElasticModel.cc.

225{
226 G4double electronEnergy0 = aDynamicElectron->GetKineticEnergy();
227 auto materialID = couple->GetMaterial()->GetIndex();
228 auto p = aDynamicElectron->GetParticleDefinition();
229
230 if (p != fpParticle) {
231 G4Exception("G4DNACPA100ElasticModel::SampleSecondaries", "em00436", FatalException,
232 "This particle is not applied for this model");
233 }
234 if (electronEnergy0 < fKillBelowEnergy) {
235 return;
236 }
237 G4double cosTheta = fpModelData->RandomizeCosTheta(electronEnergy0, materialID);
238 G4double phi = 2. * CLHEP::pi * G4UniformRand();
239
240 const G4ThreeVector& zVers = aDynamicElectron->GetMomentumDirection();
241
242 G4double CT1, ST1, CF1, SF1, CT2, ST2, CF2, SF2;
243 G4double sinTheta = std::sqrt(1 - cosTheta * cosTheta);
244
245 CT1 = zVers.z();
246 ST1 = std::sqrt(1. - CT1 * CT1);
247
248 if (ST1 != 0)
249 CF1 = zVers.x() / ST1;
250 else
251 CF1 = std::cos(2. * CLHEP::pi * G4UniformRand());
252 if (ST1 != 0)
253 SF1 = zVers.y() / ST1;
254 else
255 SF1 = std::sqrt(1. - CF1 * CF1);
256
257 G4double A3, A4, A5, A2, A1;
258
259 A3 = sinTheta * std::cos(phi);
260 A4 = A3 * CT1 + ST1 * cosTheta;
261 A5 = sinTheta * std::sin(phi);
262 A2 = A4 * SF1 + A5 * CF1;
263 A1 = A4 * CF1 - A5 * SF1;
264
265 CT2 = CT1 * cosTheta - ST1 * A3;
266 ST2 = std::sqrt(1. - CT2 * CT2);
267
268 if (ST2 == 0) ST2 = 1E-6;
269 CF2 = A1 / ST2;
270 SF2 = A2 / ST2;
271 G4ThreeVector zPrimeVers(ST2 * CF2, ST2 * SF2, CT2);
272
273 fParticleChangeForGamma->ProposeMomentumDirection(zPrimeVers.unit());
274
275 auto EnergyDeposit = fpModelData->GetElasticLevel(materialID) * (1. - cosTheta) * electronEnergy0;
276 fParticleChangeForGamma->ProposeLocalEnergyDeposit(EnergyDeposit);
277 if (statCode) {
278 fParticleChangeForGamma->SetProposedKineticEnergy(electronEnergy0);
279 }
280 else {
281 auto newEnergy = electronEnergy0 - EnergyDeposit;
282 fParticleChangeForGamma->SetProposedKineticEnergy(newEnergy);
283 }
284}
#define G4UniformRand()
Definition Randomize.hh:52
double z() const
double x() const
double y() const
G4double GetElasticLevel(const std::size_t &l)
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(const G4ThreeVector &Pfinal)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)

◆ SelectStationary()

void G4DNACPA100ElasticModel::SelectStationary ( G4bool input)
inline

Definition at line 134 of file G4DNACPA100ElasticModel.hh.

135{
136 statCode = input;
137}

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