Geant4  10.02
G4DNAMeltonAttachmentModel.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 // $Id: G4DNAMeltonAttachmentModel.cc 85244 2014-10-27 08:24:13Z gcosmo $
27 //
28 
29 // Created by Z. Francis
30 
32 #include "G4SystemOfUnits.hh"
33 #include "G4DNAChemistryManager.hh"
35 
36 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
37 
38 using namespace std;
39 
40 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
41 
43  const G4String& nam) :
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 }
73 
74 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
75 
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 }
91 
92 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
93 
95  const G4DataVector& /*cuts*/)
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 }
162 
163 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
164 
165 G4double
167  const G4ParticleDefinition* particleDefinition,
168  G4double ekin,
169  G4double,
170  G4double)
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 }
230 
231 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
232 
233 void G4DNAMeltonAttachmentModel::SampleSecondaries(std::vector<G4DynamicParticle*>* /*fvect*/,
234  const G4MaterialCutsCouple* /*couple*/,
235  const G4DynamicParticle* aDynamicElectron,
236  G4double,
237  G4double)
238 {
239 
240  if (verboseLevel > 3) G4cout
241  << "Calling SampleSecondaries() of G4DNAMeltonAttachmentModel" << G4endl;
242 
243  // Electron is killed
244 
245  G4double electronEnergy0 = aDynamicElectron->GetKineticEnergy();
249 
251  {
254  }
255  return;
256  }
static const double cm
Definition: G4SIunits.hh:118
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:89
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:641
G4double GetKineticEnergy() const
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:634
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:604
size_t GetIndex() const
Definition: G4Material.hh:262
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
virtual G4bool LoadData(const G4String &argFileName)
const G4String & GetParticleName() const
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:725
virtual G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
G4ParticleChangeForGamma * fParticleChangeForGamma
const std::vector< G4double > * fpWaterDensity
G4GLOB_DLL std::ostream G4cout
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
G4DNAMeltonAttachmentModel(const G4ParticleDefinition *p=0, const G4String &nam="DNAMeltonAttachmentModel")
virtual G4double FindValue(G4double e, G4int componentId=0) const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static G4DNAChemistryManager * Instance()
static const double eV
Definition: G4SIunits.hh:212
static G4DNAMolecularMaterial * Instance()
void CreateWaterMolecule(ElectronicModification, G4int, const G4Track *)
Method used by DNA physics model to create a water molecule.
const G4Track * GetCurrentTrack() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:732
static const G4double pos
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:134