Geant4  10.03.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4DNAMeltonAttachmentModel Class Reference

#include <G4DNAMeltonAttachmentModel.hh>

Inheritance diagram for G4DNAMeltonAttachmentModel:
Collaboration diagram for G4DNAMeltonAttachmentModel:

Public Member Functions

 G4DNAMeltonAttachmentModel (const G4ParticleDefinition *p=0, const G4String &nam="DNAMeltonAttachmentModel")
 
virtual ~G4DNAMeltonAttachmentModel ()
 
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 SetDissociationFlag (G4bool)
 
G4bool GetDissociationFlag ()
 
void SelectStationary (G4bool input)
 
- 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 SetFluctuationFlag (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

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

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

Definition at line 44 of file G4DNAMeltonAttachmentModel.hh.

Constructor & Destructor Documentation

G4DNAMeltonAttachmentModel::G4DNAMeltonAttachmentModel ( const G4ParticleDefinition p = 0,
const G4String nam = "DNAMeltonAttachmentModel" 
)

Definition at line 44 of file G4DNAMeltonAttachmentModel.cc.

45  :
46  G4VEmModel(nam), isInitialised(false)
47 {
48  fpWaterDensity = 0;
49 
52 
53  verboseLevel = 0;
54  // Verbosity scale:
55  // 0 = nothing
56  // 1 = warning for energy non-conservation
57  // 2 = details of energy budget
58  // 3 = calculation of cross sections, file openings, sampling of atoms
59  // 4 = entering in methods
60 
61 #ifdef MELTON_VERBOSE
62  if (verboseLevel > 0)
63  {
64  G4cout << "Melton Attachment model is constructed "
65  << G4endl
66  << "Energy range: "
67  << LowEnergyLimit() / eV << " eV - "
68  << HighEnergyLimit() / eV << " eV"
69  << G4endl;
70  }
71 #endif
72 
74  fDissociationFlag = true;
75  fData = 0;
76 
77  // Selection of stationary mode
78 
79  statCode = false;
80 }
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:643
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:636
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:68
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:732
G4ParticleChangeForGamma * fParticleChangeForGamma
G4GLOB_DLL std::ostream G4cout
static constexpr double eV
Definition: G4SIunits.hh:215
#define G4endl
Definition: G4ios.hh:61
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:739

Here is the call graph for this function:

G4DNAMeltonAttachmentModel::~G4DNAMeltonAttachmentModel ( )
virtual

Definition at line 84 of file G4DNAMeltonAttachmentModel.cc.

85 {
86  if(fData) delete fData;
87 }

Member Function Documentation

G4double G4DNAMeltonAttachmentModel::CrossSectionPerVolume ( const G4Material material,
const G4ParticleDefinition p,
G4double  ekin,
G4double  emin,
G4double  emax 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 183 of file G4DNAMeltonAttachmentModel.cc.

188 {
189 #ifdef MELTON_VERBOSE
190  if (verboseLevel > 3)
191  G4cout
192  << "Calling CrossSectionPerVolume() of G4DNAMeltonAttachmentModel"
193  << G4endl;
194 #endif
195 
196  // Calculate total cross section for model
197 
198  G4double sigma = 0.;
199 
200  G4double waterDensity = (*fpWaterDensity)[material->GetIndex()];
201 
202  if(waterDensity != 0.0)
203  {
204  if (ekin >= LowEnergyLimit() && ekin < HighEnergyLimit())
205  // necessaire ?
206  {
207  sigma = fData->FindValue(ekin);
208  }
209 
210 #ifdef MELTON_VERBOSE
211  if (verboseLevel > 2)
212  {
213  G4cout << "__________________________________" << G4endl;
214  G4cout << "=== G4DNAMeltonAttachmentModel - XS INFO START" << G4endl;
215  G4cout << "--- Kinetic energy(eV)=" << ekin/eV
216  << " particle : " << particleDefinition->GetParticleName()
217  << G4endl;
218  G4cout << "--- Cross section per water molecule (cm^2)="
219  << sigma/cm/cm << G4endl;
220  G4cout << "--- Cross section per water molecule (cm^-1)="
221  << sigma*waterDensity/(1./cm) << G4endl;
222  G4cout << "--- G4DNAMeltonAttachmentModel - XS INFO END" << G4endl;
223  }
224 #endif
225  } // if water
226 
227  return sigma*waterDensity;
228 }
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:643
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:636
size_t GetIndex() const
Definition: G4Material.hh:262
G4GLOB_DLL std::ostream G4cout
static constexpr double cm
Definition: G4SIunits.hh:119
static constexpr double eV
Definition: G4SIunits.hh:215
virtual G4double FindValue(G4double e, G4int componentId=0) const
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

G4bool G4DNAMeltonAttachmentModel::GetDissociationFlag ( )
inline

Definition at line 97 of file G4DNAMeltonAttachmentModel.hh.

98 {
99  return fDissociationFlag;
100 }
void G4DNAMeltonAttachmentModel::Initialise ( const G4ParticleDefinition particle,
const G4DataVector  
)
virtual

Implements G4VEmModel.

Definition at line 91 of file G4DNAMeltonAttachmentModel.cc.

93 {
94 #ifdef MELTON_VERBOSE
95  if (verboseLevel > 3)
96  G4cout
97  << "Calling G4DNAMeltonAttachmentModel::Initialise()" << G4endl;
98 #endif
99 
100  // ONLY ELECTRON
101 
102  if(particle->GetParticleName() != "e-")
103  {
104  G4Exception("G4DNAMeltonAttachmentModel::Initialise",
105  "em0002",
107  "Model not applicable to particle type.");
108  }
109 
110  // Energy limits
111 
112  if (LowEnergyLimit() < 4.*eV)
113  {
114  G4ExceptionDescription errMsg;
115  errMsg << "G4DNAMeltonAttachmentModel: low energy limit increased from " <<
116  LowEnergyLimit()/eV << " eV to " << 4. << " eV" << G4endl;
117 
118  G4Exception("G4DNAMeltonAttachmentModel::Initialise",
119  "Melton_LowerEBoundary",
120  JustWarning,
121  errMsg);
122 
124  }
125 
126  if (HighEnergyLimit() > 13.*eV)
127  {
128  G4ExceptionDescription errMsg;
129  errMsg << "G4DNAMeltonAttachmentModel: high energy limit decreased from " <<
130  HighEnergyLimit()/eV << " eV to " << 13. << " eV" << G4endl;
131 
132  G4Exception("G4DNAMeltonAttachmentModel::Initialise",
133  "Melton_HigherEBoundary",
134  JustWarning,
135  errMsg);
136 
137  SetHighEnergyLimit(13.*eV);
138  }
139 
140  // Reading of data files
141 
142  G4double scaleFactor = 1e-18*cm2;
143 
144  // For total cross section
145  G4String fileElectron("dna/sigma_attachment_e_melton");
146 
148  eV, scaleFactor);
149  fData->LoadData(fileElectron);
150 
151 
152 #ifdef MELTON_VERBOSE
153  if( verboseLevel >0)
154  {
155  if (verboseLevel > 2)
156  {
157  G4cout << "Loaded cross section data for Melton Attachment model" << G4endl;
158  }
159 
160  G4cout << "Melton Attachment model is initialized " << G4endl
161  << "Energy range: "
162  << LowEnergyLimit() / eV << " eV - "
163  << HighEnergyLimit() / eV << " eV"
164  << G4endl;
165  }
166 #endif
167 
168  // Initialize water density pointer
169  fpWaterDensity = G4DNAMolecularMaterial::Instance()->
170  GetNumMolPerVolTableFor(G4Material::GetMaterial("G4_WATER"));
171 
172  if (isInitialised)
173  {
174  return;
175  }
177  isInitialised = true;
178 }
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:643
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
static constexpr double cm2
Definition: G4SIunits.hh:120
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:636
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:602
virtual G4bool LoadData(const G4String &argFileName)
const G4String & GetParticleName() const
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:732
G4ParticleChangeForGamma * fParticleChangeForGamma
G4GLOB_DLL std::ostream G4cout
static constexpr double eV
Definition: G4SIunits.hh:215
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static G4DNAMolecularMaterial * Instance()
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:739
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:132

Here is the call graph for this function:

void G4DNAMeltonAttachmentModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  ,
const G4MaterialCutsCouple ,
const G4DynamicParticle aDynamicElectron,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Definition at line 234 of file G4DNAMeltonAttachmentModel.cc.

239 {
240 
241 #ifdef MELTON_VERBOSE
242  if (verboseLevel > 3)
243  G4cout
244  << "Calling SampleSecondaries() of G4DNAMeltonAttachmentModel" << G4endl;
245 #endif
246 
247  // Electron is killed
248 
249  G4double electronEnergy0 = aDynamicElectron->GetKineticEnergy();
250 
251  if (!statCode)
252  {
256  }
257 
258  else
259  {
262  }
263 
264  if(fDissociationFlag)
265  {
267  CreateWaterMolecule(eDissociativeAttachment,
268  -1,
270  }
271  return;
272 }
G4double GetKineticEnergy() const
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4ParticleChangeForGamma * fParticleChangeForGamma
G4GLOB_DLL std::ostream G4cout
static G4DNAChemistryManager * Instance()
const G4Track * GetCurrentTrack() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)

Here is the call graph for this function:

void G4DNAMeltonAttachmentModel::SelectStationary ( G4bool  input)
inline

Definition at line 104 of file G4DNAMeltonAttachmentModel.hh.

105 {
106  statCode = input;
107 }
void G4DNAMeltonAttachmentModel::SetDissociationFlag ( G4bool  flag)
inline

Definition at line 92 of file G4DNAMeltonAttachmentModel.hh.

93 {
94  fDissociationFlag = flag;
95 }

Member Data Documentation

G4ParticleChangeForGamma* G4DNAMeltonAttachmentModel::fParticleChangeForGamma
protected

Definition at line 72 of file G4DNAMeltonAttachmentModel.hh.


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