Geant4  10.02.p03
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 ()
 
- 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, const G4ParticleDefinition *, G4double)
 
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 *> *)
 
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=0)
 
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

G4ParticleChangeForGamma * fParticleChangeForGamma
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChange * pParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 

Private Types

typedef std::map< G4String, G4String, std::less< G4String > > MapFile
 
typedef std::map< G4String, G4DNACrossSectionDataSet *, std::less< G4String > > MapData
 

Private Member Functions

G4DNAMeltonAttachmentModeloperator= (const G4DNAMeltonAttachmentModel &right)
 
 G4DNAMeltonAttachmentModel (const G4DNAMeltonAttachmentModel &)
 

Private Attributes

const std::vector< G4double > * fpWaterDensity
 
G4double lowEnergyLimit
 
G4double highEnergyLimit
 
G4double lowEnergyLimitOfModel
 
G4bool isInitialised
 
G4int verboseLevel
 
G4bool fDissociationFlag
 
MapFile tableFile
 
MapData tableData
 

Additional Inherited Members

- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLoss * GetParticleChangeForLoss ()
 
G4ParticleChangeForGamma * GetParticleChangeForGamma ()
 
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.

Member Typedef Documentation

◆ MapData

Definition at line 91 of file G4DNAMeltonAttachmentModel.hh.

◆ MapFile

typedef std::map<G4String,G4String,std::less<G4String> > G4DNAMeltonAttachmentModel::MapFile
private

Definition at line 88 of file G4DNAMeltonAttachmentModel.hh.

Constructor & Destructor Documentation

◆ G4DNAMeltonAttachmentModel() [1/2]

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

Definition at line 42 of file G4DNAMeltonAttachmentModel.cc.

43  :
44  G4VEmModel(nam), isInitialised(false)
45 {
46 // nistwater = G4NistManager::Instance()->FindOrBuildMaterial("G4_WATER");
47  fpWaterDensity = 0;
48 
49  lowEnergyLimit = 4 * eV;
51  highEnergyLimit = 13 * eV;
54 
55  verboseLevel = 0;
56  // Verbosity scale:
57  // 0 = nothing
58  // 1 = warning for energy non-conservation
59  // 2 = details of energy budget
60  // 3 = calculation of cross sections, file openings, sampling of atoms
61  // 4 = entering in methods
62 
63  if (verboseLevel > 0)
64  {
65  G4cout << "Melton Attachment model is constructed " << G4endl<< "Energy range: "
66  << lowEnergyLimit / eV << " eV - "
67  << highEnergyLimit / eV << " eV"
68  << G4endl;
69  }
71  fDissociationFlag = true;
72 }
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:69
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:725
G4ParticleChangeForGamma * fParticleChangeForGamma
const std::vector< G4double > * fpWaterDensity
G4GLOB_DLL std::ostream G4cout
static const double eV
Definition: G4SIunits.hh:212
#define G4endl
Definition: G4ios.hh:61
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:732
Here is the call graph for this function:

◆ ~G4DNAMeltonAttachmentModel()

G4DNAMeltonAttachmentModel::~G4DNAMeltonAttachmentModel ( )
virtual

Definition at line 76 of file G4DNAMeltonAttachmentModel.cc.

77 {
78  // For total cross section
79 
80  std::map<G4String, G4DNACrossSectionDataSet*, std::less<G4String> >::iterator pos;
81 
82  for (pos = tableData.begin(); pos != tableData.end(); ++pos)
83  {
84  G4DNACrossSectionDataSet* table = pos->second;
85  delete table;
86  }
87 
88  // For final state
89 
90 }
static const G4double pos

◆ G4DNAMeltonAttachmentModel() [2/2]

G4DNAMeltonAttachmentModel::G4DNAMeltonAttachmentModel ( const G4DNAMeltonAttachmentModel )
private

Member Function Documentation

◆ CrossSectionPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 166 of file G4DNAMeltonAttachmentModel.cc.

171 {
172  if (verboseLevel > 3) G4cout
173  << "Calling CrossSectionPerVolume() of G4DNAMeltonAttachmentModel"
174  << G4endl;
175 
176  // Calculate total cross section for model
177 
178  G4double sigma=0;
179 
180  G4double waterDensity = (*fpWaterDensity)[material->GetIndex()];
181 
182  if(waterDensity!= 0.0)
183  // if (material == nistwater || material->GetBaseMaterial() == nistwater)
184  {
185  const G4String& particleName = particleDefinition->GetParticleName();
186 
187  if (ekin >= lowEnergyLimit && ekin < highEnergyLimit)
188  {
189 
190  std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
191  pos = tableData.find(particleName);
192 
193  if (pos != tableData.end())
194  {
195  G4DNACrossSectionDataSet* table = pos->second;
196  if (table != 0)
197  {
198  sigma = table->FindValue(ekin);
199  }
200  }
201  else
202  {
203  G4Exception("G4DNAMeltonAttachmentModel::ComputeCrossSectionPerVolume",
204  "em0002",
205  FatalException,"Model not applicable to particle type.");
206  }
207  }
208 
209  if (verboseLevel > 2)
210  {
211  G4cout << "__________________________________" << G4endl;
212  G4cout << "=== G4DNAMeltonAttachmentModel - XS INFO START" << G4endl;
213  G4cout << "--- Kinetic energy(eV)=" << ekin/eV
214  << " particle : " << particleDefinition->GetParticleName() << G4endl;
215  G4cout << "--- Cross section per water molecule (cm^2)="
216  << sigma/cm/cm << G4endl;
217  G4cout << "--- Cross section per water molecule (cm^-1)="
218  << sigma*waterDensity/(1./cm) << G4endl;
219  // G4cout << "--- Cross section per water molecule (cm^-1)="
220  // << sigma*material->GetAtomicNumDensityVector()[1]/(1./cm)
221  // << G4endl;
222  G4cout << "--- G4DNAMeltonAttachmentModel - XS INFO END" << G4endl;
223  }
224 
225  } // if water
226 
227  return sigma*waterDensity;
228 // return sigma*material->GetAtomicNumDensityVector()[1];
229 }
static const double cm
Definition: G4SIunits.hh:118
size_t GetIndex() const
Definition: G4Material.hh:262
virtual G4double FindValue(G4double e, G4int componentId=0) const
G4GLOB_DLL std::ostream G4cout
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static const double eV
Definition: G4SIunits.hh:212
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
static const G4double pos
Here is the call graph for this function:

◆ GetDissociationFlag()

G4bool G4DNAMeltonAttachmentModel::GetDissociationFlag ( )
inline

Definition at line 106 of file G4DNAMeltonAttachmentModel.hh.

107 {
108  return fDissociationFlag;
109 }

◆ Initialise()

void G4DNAMeltonAttachmentModel::Initialise ( const G4ParticleDefinition ,
const G4DataVector  
)
virtual

Implements G4VEmModel.

Definition at line 94 of file G4DNAMeltonAttachmentModel.cc.

96 {
97 
98  if (verboseLevel > 3) G4cout
99  << "Calling G4DNAMeltonAttachmentModel::Initialise()" << G4endl;
100 
101  // Energy limits
102 
104  {
105  G4cout << "G4DNAMeltonAttachmentModel: low energy limit increased from " <<
106  LowEnergyLimit()/eV << " eV to " << lowEnergyLimit/eV << " eV" << G4endl;
108  }
109 
111  {
112  G4cout << "G4DNAMeltonAttachmentModel: high energy limit decreased from " <<
113  HighEnergyLimit()/eV << " eV to " << highEnergyLimit/eV << " eV" << G4endl;
115  }
116 
117  // Reading of data files
118 
119  G4double scaleFactor = 1e-18*cm*cm;
120 
121  G4String fileElectron("dna/sigma_attachment_e_melton");
122 
125 
126  // ELECTRON
127 
128  // For total cross section
129 
130  electron = electronDef->GetParticleName();
131 
132  tableFile[electron] = fileElectron;
133 
134  G4DNACrossSectionDataSet* tableE =
135  new G4DNACrossSectionDataSet(new G4LogLogInterpolation, eV,scaleFactor );
136  tableE->LoadData(fileElectron);
137  tableData[electron] = tableE;
138 
139  //
140 
141  if (verboseLevel > 2)
142  G4cout << "Loaded cross section data for Melton Attachment model" << G4endl;
143 
144  if( verboseLevel>0 )
145  {
146  G4cout << "Melton Attachment model is initialized " << G4endl
147  << "Energy range: "
148  << LowEnergyLimit() / eV << " eV - "
149  << HighEnergyLimit() / eV << " eV"
150  << G4endl;
151  }
152  // Initialize water density pointer
154  GetNumMolPerVolTableFor(G4Material::GetMaterial("G4_WATER"));
155 
156  if (isInitialised)
157  { return;}
159  isInitialised = true;
160 
161 }
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:641
static const double cm
Definition: G4SIunits.hh:118
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:89
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:604
virtual G4bool LoadData(const G4String &argFileName)
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:725
G4ParticleChangeForGamma * fParticleChangeForGamma
const std::vector< G4double > * fpWaterDensity
const G4String & GetParticleName() const
G4GLOB_DLL std::ostream G4cout
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:634
static const double eV
Definition: G4SIunits.hh:212
static G4DNAMolecularMaterial * Instance()
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:732
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:134
Here is the call graph for this function:

◆ operator=()

G4DNAMeltonAttachmentModel& G4DNAMeltonAttachmentModel::operator= ( const G4DNAMeltonAttachmentModel right)
private

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 233 of file G4DNAMeltonAttachmentModel.cc.

238 {
239 
240  if (verboseLevel > 3) G4cout
241  << "Calling SampleSecondaries() of G4DNAMeltonAttachmentModel" << G4endl;
242 
243  // Electron is killed
244 
245  G4double electronEnergy0 = aDynamicElectron->GetKineticEnergy();
246  fParticleChangeForGamma->SetProposedKineticEnergy(0.);
247  fParticleChangeForGamma->ProposeTrackStatus(fStopAndKill);
248  fParticleChangeForGamma->ProposeLocalEnergyDeposit(electronEnergy0);
249 
251  {
253  fParticleChangeForGamma->GetCurrentTrack());
254  }
255  return;
256  }
G4ParticleChangeForGamma * fParticleChangeForGamma
G4double GetKineticEnergy() const
G4GLOB_DLL std::ostream G4cout
static G4DNAChemistryManager * Instance()
void CreateWaterMolecule(ElectronicModification, G4int, const G4Track *)
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:

◆ SetDissociationFlag()

void G4DNAMeltonAttachmentModel::SetDissociationFlag ( G4bool  flag)
inline

Definition at line 101 of file G4DNAMeltonAttachmentModel.hh.

102 {
103  fDissociationFlag = flag;
104 }

Member Data Documentation

◆ fDissociationFlag

G4bool G4DNAMeltonAttachmentModel::fDissociationFlag
private

Definition at line 84 of file G4DNAMeltonAttachmentModel.hh.

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4DNAMeltonAttachmentModel::fParticleChangeForGamma
protected

Definition at line 73 of file G4DNAMeltonAttachmentModel.hh.

◆ fpWaterDensity

const std::vector<G4double>* G4DNAMeltonAttachmentModel::fpWaterDensity
private

Definition at line 77 of file G4DNAMeltonAttachmentModel.hh.

◆ highEnergyLimit

G4double G4DNAMeltonAttachmentModel::highEnergyLimit
private

Definition at line 80 of file G4DNAMeltonAttachmentModel.hh.

◆ isInitialised

G4bool G4DNAMeltonAttachmentModel::isInitialised
private

Definition at line 82 of file G4DNAMeltonAttachmentModel.hh.

◆ lowEnergyLimit

G4double G4DNAMeltonAttachmentModel::lowEnergyLimit
private

Definition at line 79 of file G4DNAMeltonAttachmentModel.hh.

◆ lowEnergyLimitOfModel

G4double G4DNAMeltonAttachmentModel::lowEnergyLimitOfModel
private

Definition at line 81 of file G4DNAMeltonAttachmentModel.hh.

◆ tableData

MapData G4DNAMeltonAttachmentModel::tableData
private

Definition at line 92 of file G4DNAMeltonAttachmentModel.hh.

◆ tableFile

MapFile G4DNAMeltonAttachmentModel::tableFile
private

Definition at line 89 of file G4DNAMeltonAttachmentModel.hh.

◆ verboseLevel

G4int G4DNAMeltonAttachmentModel::verboseLevel
private

Definition at line 83 of file G4DNAMeltonAttachmentModel.hh.


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