Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4DNAMeltonAttachmentModel.cc
이 파일의 문서화 페이지로 가기
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 98733 2016-08-09 10:51:58Z 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 //#define MELTON_VERBOSE // prevent checking conditions at run time
41 
42 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
43 
45  const G4String& nam) :
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 }
81 
82 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
83 
85 {
86  if(fData) delete fData;
87 }
88 
89 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
90 
92  const G4DataVector& /*cuts*/)
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
170  GetNumMolPerVolTableFor(G4Material::GetMaterial("G4_WATER"));
171 
172  if (isInitialised)
173  {
174  return;
175  }
177  isInitialised = true;
178 }
179 
180 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
181 
182 G4double
184  const G4ParticleDefinition*,
185  G4double ekin,
186  G4double,
187  G4double)
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 }
229 
230 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
231 
232 void
234 SampleSecondaries(std::vector<G4DynamicParticle*>* /*fvect*/,
235  const G4MaterialCutsCouple* /*couple*/,
236  const G4DynamicParticle* aDynamicElectron,
237  G4double,
238  G4double)
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 
265  {
267  CreateWaterMolecule(eDissociativeAttachment,
268  -1,
270  }
271  return;
272 }
virtual G4bool LoadData(const G4String &argFileName)
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:45
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:707
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
size_t GetIndex() const
Definition: G4Material.hh:262
#define G4endl
Definition: G4ios.hh:61
const G4String & GetParticleName() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:608
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:616
double G4double
Definition: G4Types.hh:76
static constexpr double cm2
Definition: G4SIunits.hh:120
static G4DNAMolecularMaterial * Instance()
G4DNAMeltonAttachmentModel(const G4ParticleDefinition *p=0, const G4String &nam="DNAMeltonAttachmentModel")
virtual G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
static constexpr double eV
Definition: G4SIunits.hh:215
const G4Track * GetCurrentTrack() const
G4ParticleChangeForGamma * fParticleChangeForGamma
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.hh:65
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:714
G4double GetKineticEnergy() const
virtual G4double FindValue(G4double e, G4int componentId=0) const
static constexpr double cm
Definition: G4SIunits.hh:119
G4GLOB_DLL std::ostream G4cout
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
G4DNACrossSectionDataSet * fData
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:131
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:609
void ProposeTrackStatus(G4TrackStatus status)
static G4DNAChemistryManager * Instance()
const std::vector< G4double > * fpWaterDensity