Geant4 10.7.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

typedef MODEL Model
 

Public Member Functions

 G4TDNAOneStepThermalizationModel (const G4ParticleDefinition *p=0, const G4String &nam="DNAOneStepThermalizationModel")
 
virtual ~G4TDNAOneStepThermalizationModel ()
 
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)
 
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 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 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 CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=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 *, G4double &eloss, G4double &niel, G4double length)
 
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 ModelDescription (std::ostream &outFile) const
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
std::vector< G4EmElementSelector * > * GetElementSelectors ()
 
void SetElementSelectors (std::vector< G4EmElementSelector * > *)
 
virtual 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)
 
G4int SelectRandomAtomNumber (const G4Material *)
 
G4int SelectIsotopeNumber (const G4Element *)
 
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 LPMFlag () 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 SetLPMFlag (G4bool val)
 
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 *)
 
const G4ElementGetCurrentElement () const
 
const G4IsotopeGetCurrentIsotope () const
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

Protected Attributes

const std::vector< G4double > * fpWaterDensity
 
G4ParticleChangeForGammafpParticleChangeForGamma
 
G4bool fIsInitialised
 
G4int fVerboseLevel
 
std::unique_ptr< G4NavigatorfpNavigator
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const G4MaterialpBaseMaterial
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 
G4bool lossFlucFlag
 
G4double inveplus
 
G4double pFactor
 

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>
typedef MODEL G4TDNAOneStepThermalizationModel< MODEL >::Model

Definition at line 140 of file G4DNAOneStepThermalizationModel.hh.

Constructor & Destructor Documentation

◆ G4TDNAOneStepThermalizationModel()

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

◆ ~G4TDNAOneStepThermalizationModel()

template<typename MODEL >
G4TDNAOneStepThermalizationModel< MODEL >::~G4TDNAOneStepThermalizationModel
virtual

Definition at line 70 of file G4DNAOneStepThermalizationModel.hpp.

71{
72 // if(fpNavigator && fpNavigator->GetNavigatorState())
73 // delete fpNavigator->GetNavigatorState();
74}

Member Function Documentation

◆ CrossSectionPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 125 of file G4DNAOneStepThermalizationModel.hpp.

131{
132#ifdef MODEL_VERBOSE
133 if(fVerboseLevel > 1)
134 G4cout << "Calling CrossSectionPerVolume() of G4DNAOneStepThermalizationModel"
135 << G4endl;
136#endif
137
138 if(ekin > HighEnergyLimit()){
139 return 0.0;
140 }
141
142 G4double waterDensity = (*fpWaterDensity)[material->GetIndex()];
143
144 if(waterDensity!= 0.0){
145 return DBL_MAX;
146 }
147 return 0.;
148}
double G4double
Definition: G4Types.hh:83
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
size_t GetIndex() const
Definition: G4Material.hh:258
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:645
#define DBL_MAX
Definition: templates.hh:62

◆ GetPenetration()

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

Definition at line 160 of file G4DNAOneStepThermalizationModel.hpp.

162{
163 return MODEL::GetPenetration(k, displacement);
164}

◆ GetRmean()

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

Definition at line 152 of file G4DNAOneStepThermalizationModel.hpp.

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

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 78 of file G4DNAOneStepThermalizationModel.hpp.

81{
82#ifdef MODEL_VERBOSE
84 G4cout << "Calling G4DNAOneStepThermalizationModel::Initialise()"
85 << G4endl;
86#endif
87 if (particleDefinition->GetParticleName() != "e-")
88 {
90 errMsg << "G4DNAOneStepThermalizationModel can only be applied "
91 "to electrons";
92 G4Exception("G4DNAOneStepThermalizationModel::CrossSectionPerVolume",
93 "G4DNAOneStepThermalizationModel001",
95 return;
96 }
97
99 {
100 fIsInitialised = true;
102 }
103
104 G4Navigator* navigator =
106 GetNavigatorForTracking();
107
108 fpNavigator.reset(new G4Navigator());
109
110 if(navigator){ // add these checks for testing mode
111 auto world=navigator->GetWorldVolume();
112 if(world){
113 fpNavigator->SetWorldVolume(world);
114 //fNavigator->NewNavigatorState();
115 }
116 }
117
120 GetNumMolPerVolTableFor(G4Material::GetMaterial("G4_WATER"));
121}
@ FatalErrorInArgument
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
static G4DNAMolecularMaterial * Instance()
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:651
G4VPhysicalVolume * GetWorldVolume() const
const G4String & GetParticleName() const
static G4TransportationManager * GetTransportationManager()
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:133

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 168 of file G4DNAOneStepThermalizationModel.hpp.

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

◆ SetVerbose()

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

Definition at line 160 of file G4DNAOneStepThermalizationModel.hh.

160 {
161 fVerboseLevel = flag;
162 }

Member Data Documentation

◆ fIsInitialised

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

Definition at line 173 of file G4DNAOneStepThermalizationModel.hh.

◆ fpNavigator

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

Definition at line 175 of file G4DNAOneStepThermalizationModel.hh.

◆ 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: