Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4TDNAOneStepThermalizationModel< MODEL > Class Template Reference

#include <G4DNAOneStepThermalizationModel.hh>

Inheritance diagram for G4TDNAOneStepThermalizationModel< MODEL >:
Collaboration 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 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 *, 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)
 
G4int SelectIsotopeNumber (const G4Element *)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, 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 *)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
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)
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetLPMFlag (G4bool val)
 
void SetDeexcitationFlag (G4bool val)
 
void SetForceBuildTable (G4bool val)
 
void SetMasterThread (G4bool val)
 
G4bool IsMaster () 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)
 

Protected Attributes

const std::vector< G4double > * fpWaterDensity
 
G4ParticleChangeForGammafParticleChangeForGamma
 
G4bool fIsInitialised
 
G4int fVerboseLevel
 
G4NavigatorfNavigator
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 

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 *)
 
- Static Protected Attributes inherited from G4VEmModel
static const G4double inveplus = 1.0/CLHEP::eplus
 

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 100 of file G4DNAOneStepThermalizationModel.hh.

Member Typedef Documentation

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

Definition at line 103 of file G4DNAOneStepThermalizationModel.hh.

Constructor & Destructor Documentation

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

Definition at line 56 of file G4DNAOneStepThermalizationModel.hpp.

Here is the call graph for this function:

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

Definition at line 72 of file G4DNAOneStepThermalizationModel.hpp.

73 {
74  if(fNavigator)
75  {
76  // if(fNavigator->GetNavigatorState())
77  // delete fNavigator->GetNavigatorState();
78  delete fNavigator;
79  }
80 }

Member Function Documentation

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 132 of file G4DNAOneStepThermalizationModel.hpp.

137 {
138 #ifdef MODEL_VERBOSE
139  if(fVerboseLevel > 1)
140  G4cout << "Calling CrossSectionPerVolume() of G4DNAOneStepThermalizationModel"
141  << G4endl;
142 #endif
143 
144  if(ekin > HighEnergyLimit()){
145  return 0.0;
146  }
147 
148  G4double waterDensity = (*fpWaterDensity)[material->GetIndex()];
149 
150  if(waterDensity!= 0.0){
151  return DBL_MAX;
152  }
153  return 0.;
154 }
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:633
size_t GetIndex() const
Definition: G4Material.hh:262
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
#define DBL_MAX
Definition: templates.hh:83

Here is the call graph for this function:

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

Definition at line 167 of file G4DNAOneStepThermalizationModel.hpp.

168 {
169  return MODEL::GetPenetration(k, displacement);
170 }
template<typename MODEL >
double G4TDNAOneStepThermalizationModel< MODEL >::GetRmean ( double  energy)

Definition at line 158 of file G4DNAOneStepThermalizationModel.hpp.

158  {
159  return MODEL::GetRmean(k);
160 }
template<typename MODEL >
void G4TDNAOneStepThermalizationModel< MODEL >::Initialise ( const G4ParticleDefinition particleDefinition,
const G4DataVector  
)
virtual

Implements G4VEmModel.

Definition at line 85 of file G4DNAOneStepThermalizationModel.hpp.

87 {
88 #ifdef MODEL_VERBOSE
89  if(fVerboseLevel)
90  G4cout << "Calling G4DNAOneStepThermalizationModel::Initialise()"
91  << G4endl;
92 #endif
93  if (particleDefinition->GetParticleName() != "e-")
94  {
96  errMsg << "G4DNAOneStepThermalizationModel can only be applied "
97  "to electrons";
98  G4Exception("G4DNAOneStepThermalizationModel::CrossSectionPerVolume",
99  "G4DNAOneStepThermalizationModel001",
100  FatalErrorInArgument,errMsg);
101  return;
102  }
103 
104  if(!fIsInitialised)
105  {
106  fIsInitialised = true;
108  }
109 
110  G4Navigator* navigator =
112  GetNavigatorForTracking();
113 
114  fNavigator = new G4Navigator();
115 
116  if(navigator){ // add these checks for testing mode
117  auto world=navigator->GetWorldVolume();
118  if(world){
119  fNavigator->SetWorldVolume(world);
120  //fNavigator->NewNavigatorState();
121  }
122  }
123 
126  GetNumMolPerVolTableFor(G4Material::GetMaterial("G4_WATER"));
127 }
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
const std::vector< G4double > * fpWaterDensity
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:602
const G4String & GetParticleName() const
G4GLOB_DLL std::ostream G4cout
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static G4TransportationManager * GetTransportationManager()
static G4DNAMolecularMaterial * Instance()
void SetWorldVolume(G4VPhysicalVolume *pWorld)
#define G4endl
Definition: G4ios.hh:61
G4VPhysicalVolume * GetWorldVolume() const
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:132

Here is the call graph for this function:

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 175 of file G4DNAOneStepThermalizationModel.hpp.

180 {
181 #ifdef MODEL_VERBOSE
182  if(fVerboseLevel)
183  G4cout << "Calling SampleSecondaries() of G4DNAOneStepThermalizationModel"
184  << G4endl;
185 #endif
186 
187  G4double k = particle->GetKineticEnergy();
188 
189  if (k <= HighEnergyLimit())
190  {
193 
195  {
196  G4ThreeVector displacement(0,0,0);
197  GetPenetration(k, displacement);
198 
199  //______________________________________________________________
200  const G4Track * theIncomingTrack =
202  G4ThreeVector finalPosition(theIncomingTrack->GetPosition()+displacement);
203 
204  fNavigator->SetWorldVolume(theIncomingTrack->GetTouchable()->
205  GetVolume(theIncomingTrack->GetTouchable()->
206  GetHistoryDepth()));
207 
208  double displacementMag = displacement.mag();
209  double safety = DBL_MAX;
210  G4ThreeVector direction = displacement/displacementMag;
211 
212  //--
213  // 6/09/16 - recupere de molecular dissocation
214  double mag_displacement = displacement.mag();
215  G4ThreeVector displacement_direction = displacement/mag_displacement;
216 
217  // double step = DBL_MAX;
218  // step = fNavigator->CheckNextStep(theIncomingTrack->GetPosition(),
219  // displacement_direction,
220  // mag_displacement,
221  // safety);
222  //
223  //
224  // if(safety < mag_displacement)
225  // {
227  // finalPosition = theIncomingTrack->GetPosition()
228  // + (displacement/displacementMag)*safety*0.80;
229  // }
230  //--
231 
232  fNavigator->ResetHierarchyAndLocate(theIncomingTrack->GetPosition(),
233  direction,
234  *((G4TouchableHistory*)
235  theIncomingTrack->GetTouchable()));
236 
237  fNavigator->ComputeStep(theIncomingTrack->GetPosition(),
238  displacement/displacementMag,
239  displacementMag,
240  safety);
241 
242  if(safety <= displacementMag)
243  {
244  finalPosition = theIncomingTrack->GetPosition()
245  + (displacement/displacementMag)*safety*0.80;
246  }
247 
249  &finalPosition);
250 
252  }
253  }
254 }
G4double GetKineticEnergy() const
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:633
void CreateSolvatedElectron(const G4Track *, G4ThreeVector *finalPosition=0)
const G4ThreeVector & GetPosition() const
virtual G4double ComputeStep(const G4ThreeVector &pGlobalPoint, const G4ThreeVector &pDirection, const G4double pCurrentProposedStepLength, G4double &pNewSafety)
Definition: G4Navigator.cc:747
virtual G4VPhysicalVolume * ResetHierarchyAndLocate(const G4ThreeVector &point, const G4ThreeVector &direction, const G4TouchableHistory &h)
Definition: G4Navigator.cc:95
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4GLOB_DLL std::ostream G4cout
void GetPenetration(G4double energy, G4ThreeVector &displacement)
static constexpr double eV
Definition: G4SIunits.hh:215
static G4DNAChemistryManager * Instance()
const G4VTouchable * GetTouchable() const
void SetWorldVolume(G4VPhysicalVolume *pWorld)
const G4Track * GetCurrentTrack() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
double mag() const
#define DBL_MAX
Definition: templates.hh:83

Here is the call graph for this function:

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

Definition at line 123 of file G4DNAOneStepThermalizationModel.hh.

Member Data Documentation

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

Definition at line 136 of file G4DNAOneStepThermalizationModel.hh.

template<typename MODEL = DNA::Penetration::Meesungnoen2002>
G4Navigator* G4TDNAOneStepThermalizationModel< MODEL >::fNavigator
protected

Definition at line 138 of file G4DNAOneStepThermalizationModel.hh.

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

Definition at line 135 of file G4DNAOneStepThermalizationModel.hh.

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

Definition at line 133 of file G4DNAOneStepThermalizationModel.hh.

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

Definition at line 137 of file G4DNAOneStepThermalizationModel.hh.


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