Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4TDNAOneStepThermalizationModel< MODEL > Class Template Reference

#include <G4DNAOneStepThermalizationModel.hh>

+ Inheritance diagram for G4TDNAOneStepThermalizationModel< MODEL >:

Public Types

using Model = MODEL
 

Public Member Functions

 G4TDNAOneStepThermalizationModel (const G4ParticleDefinition *p=nullptr, const G4String &nam="DNAOneStepThermalizationModel")
 
 ~G4TDNAOneStepThermalizationModel () override
 
void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
G4double CrossSectionPerVolume (const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax) override
 
void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
 
void SetVerbose (int flag)
 
void GetPenetration (G4double energy, G4ThreeVector &displacement)
 
double GetRmean (double energy)
 
- 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)
 
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)
 
void SetMasterThread (G4bool)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

Protected Attributes

const std::vector< G4double > * fpWaterDensity
 
G4ParticleChangeForGammafpParticleChangeForGamma
 
G4bool fIsInitialised {false}
 
G4int fVerboseLevel
 
std::unique_ptr< G4NavigatorfpNavigator
 
- 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

template<typename MODEL = DNA::Penetration::Meesungnoen2002>
class G4TDNAOneStepThermalizationModel< MODEL >

When an electron reaches the highest energy domain of G4DNAOneStepThermalizationModel, it is then automatically converted into a solvated electron and displace from its original position using a published thermalization statistic.

Definition at line 137 of file G4DNAOneStepThermalizationModel.hh.

Member Typedef Documentation

◆ Model

template<typename MODEL = DNA::Penetration::Meesungnoen2002>
using G4TDNAOneStepThermalizationModel< MODEL >::Model = MODEL

Definition at line 140 of file G4DNAOneStepThermalizationModel.hh.

Constructor & Destructor Documentation

◆ G4TDNAOneStepThermalizationModel()

template<typename MODEL>
G4TDNAOneStepThermalizationModel< MODEL >::G4TDNAOneStepThermalizationModel ( const G4ParticleDefinition * p = nullptr,
const G4String & nam = "DNAOneStepThermalizationModel" )

◆ ~G4TDNAOneStepThermalizationModel()

template<typename MODEL>
G4TDNAOneStepThermalizationModel< MODEL >::~G4TDNAOneStepThermalizationModel ( )
overridedefault

Member Function Documentation

◆ CrossSectionPerVolume()

template<typename MODEL>
G4double G4TDNAOneStepThermalizationModel< MODEL >::CrossSectionPerVolume ( const G4Material * material,
const G4ParticleDefinition * p,
G4double ekin,
G4double emin,
G4double emax )
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 124 of file G4DNAOneStepThermalizationModel.hpp.

130{
131#ifdef MODEL_VERBOSE
132 if(fVerboseLevel > 1)
133 G4cout << "Calling CrossSectionPerVolume() of G4DNAOneStepThermalizationModel"
134 << G4endl;
135#endif
136
137 if(ekin > HighEnergyLimit()){
138 return 0.0;
139 }
140
141 G4double waterDensity = (*fpWaterDensity)[material->GetIndex()];
142
143 if(waterDensity!= 0.0){
144 return DBL_MAX;
145 }
146 return 0.;
147}
G4double HighEnergyLimit() const

◆ GetPenetration()

template<typename MODEL>
void G4TDNAOneStepThermalizationModel< MODEL >::GetPenetration ( G4double energy,
G4ThreeVector & displacement )

Definition at line 159 of file G4DNAOneStepThermalizationModel.hpp.

161{
163}

◆ GetRmean()

template<typename MODEL>
double G4TDNAOneStepThermalizationModel< MODEL >::GetRmean ( double energy)

Definition at line 151 of file G4DNAOneStepThermalizationModel.hpp.

151 {
152 return MODEL::GetRmean(k);
153}

◆ Initialise()

template<typename MODEL>
void G4TDNAOneStepThermalizationModel< MODEL >::Initialise ( const G4ParticleDefinition * particleDefinition,
const G4DataVector &  )
overridevirtual

Implements G4VEmModel.

Definition at line 77 of file G4DNAOneStepThermalizationModel.hpp.

80{
81#ifdef MODEL_VERBOSE
83 G4cout << "Calling G4DNAOneStepThermalizationModel::Initialise()"
84 << G4endl;
85#endif
86 if (particleDefinition->GetParticleName() != "e-")
87 {
89 errMsg << "G4DNAOneStepThermalizationModel can only be applied "
90 "to electrons";
91 G4Exception("G4DNAOneStepThermalizationModel::CrossSectionPerVolume",
92 "G4DNAOneStepThermalizationModel001",
94 return;
95 }
96
98 {
99 fIsInitialised = true;
101 }
102
106
108
109 if(navigator != nullptr){ // add these checks for testing mode
110 auto world=navigator->GetWorldVolume();
111 if(world != nullptr){
112 fpNavigator->SetWorldVolume(world);
113 //fNavigator->NewNavigatorState();
114 }
115 }
116
120}
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
static G4DNAMolecularMaterial * Instance()
static G4TransportationManager * GetTransportationManager()
G4ParticleChangeForGamma * GetParticleChangeForGamma()

◆ SampleSecondaries()

template<typename MODEL>
void G4TDNAOneStepThermalizationModel< MODEL >::SampleSecondaries ( std::vector< G4DynamicParticle * > * ,
const G4MaterialCutsCouple * ,
const G4DynamicParticle * particle,
G4double tmin,
G4double maxEnergy )
overridevirtual

Implements G4VEmModel.

Definition at line 167 of file G4DNAOneStepThermalizationModel.hpp.

173{
174#ifdef MODEL_VERBOSE
175 if(fVerboseLevel)
176 G4cout << "Calling SampleSecondaries() of G4DNAOneStepThermalizationModel"
177 << G4endl;
178#endif
179
180 G4double k = particle->GetKineticEnergy();
181
182 if (k <= HighEnergyLimit())
183 {
184 fpParticleChangeForGamma->ProposeTrackStatus(fStopAndKill);
185 fpParticleChangeForGamma->ProposeLocalEnergyDeposit(k);
186
188 {
191
192 //______________________________________________________________
193 const G4Track * theIncomingTrack =
194 fpParticleChangeForGamma->GetCurrentTrack();
196
197 fpNavigator->SetWorldVolume(theIncomingTrack->GetTouchable()->
198 GetVolume(theIncomingTrack->GetTouchable()->
199 GetHistoryDepth()));
200
201 double displacementMag = displacement.mag();
202 double safety = DBL_MAX;
204
205 //--
206 // 6/09/16 - recupere de molecular dissocation
207 double mag_displacement = displacement.mag();
209
210 // double step = DBL_MAX;
211 // step = fNavigator->CheckNextStep(theIncomingTrack->GetPosition(),
212 // displacement_direction,
213 // mag_displacement,
214 // safety);
215 //
216 //
217 // if(safety < mag_displacement)
218 // {
219 //// mag_displacement = prNewSafety;
220 // finalPosition = theIncomingTrack->GetPosition()
221 // + (displacement/displacementMag)*safety*0.80;
222 // }
223 //--
224
225 fpNavigator->ResetHierarchyAndLocate(theIncomingTrack->GetPosition(),
226 direction,
228 theIncomingTrack->GetTouchable()));
229
230 fpNavigator->ComputeStep(theIncomingTrack->GetPosition(),
233 safety);
234
236 {
237 finalPosition = theIncomingTrack->GetPosition()
239 }
240
243
244 fpParticleChangeForGamma->SetProposedKineticEnergy(25.e-3*eV);
245 }
246 }
247}
void CreateSolvatedElectron(const G4Track *, G4ThreeVector *pFinalPosition=nullptr)
static G4DNAChemistryManager * Instance()
void GetPenetration(G4double energy, G4ThreeVector &displacement)

◆ SetVerbose()

template<typename MODEL = DNA::Penetration::Meesungnoen2002>
void G4TDNAOneStepThermalizationModel< MODEL >::SetVerbose ( int flag)
inline

Definition at line 160 of file G4DNAOneStepThermalizationModel.hh.

160 {
162 }

Member Data Documentation

◆ fIsInitialised

template<typename MODEL = DNA::Penetration::Meesungnoen2002>
G4bool G4TDNAOneStepThermalizationModel< MODEL >::fIsInitialised {false}
protected

Definition at line 173 of file G4DNAOneStepThermalizationModel.hh.

173{false};

Referenced by Initialise().

◆ fpNavigator

template<typename MODEL = DNA::Penetration::Meesungnoen2002>
std::unique_ptr<G4Navigator> G4TDNAOneStepThermalizationModel< MODEL >::fpNavigator
protected

Definition at line 175 of file G4DNAOneStepThermalizationModel.hh.

Referenced by Initialise().

◆ fpParticleChangeForGamma

template<typename MODEL = DNA::Penetration::Meesungnoen2002>
G4ParticleChangeForGamma* G4TDNAOneStepThermalizationModel< MODEL >::fpParticleChangeForGamma
protected

◆ fpWaterDensity

template<typename MODEL = DNA::Penetration::Meesungnoen2002>
const std::vector<G4double>* G4TDNAOneStepThermalizationModel< MODEL >::fpWaterDensity
protected

◆ fVerboseLevel

template<typename MODEL = DNA::Penetration::Meesungnoen2002>
G4int G4TDNAOneStepThermalizationModel< MODEL >::fVerboseLevel
protected

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