Geant4  10.02
G4DNABornExcitationModel1.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: G4DNABornExcitationModel1.cc 90054 2015-05-11 18:47:32Z matkara $
27 //
28 
30 #include "G4SystemOfUnits.hh"
31 #include "G4DNAChemistryManager.hh"
33 #include <map>
34 
35 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
36 
37 using namespace std;
38 
39 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
40 
42  const G4String& nam) :
43  G4VEmModel(nam), isInitialised(false), fTableData(0)
44 {
46  fHighEnergy = 0;
47  fLowEnergy = 0;
49 
50  verboseLevel = 0;
51  // Verbosity scale:
52  // 0 = nothing
53  // 1 = warning for energy non-conservation
54  // 2 = details of energy budget
55  // 3 = calculation of cross sections, file openings, sampling of atoms
56  // 4 = entering in methods
57 
58  if (verboseLevel > 0)
59  {
60  G4cout << "Born excitation model is constructed " << G4endl;
61  }
63 
64  // Selection of stationary mode
65 
66  statCode = false;
67 }
68 
69 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
70 
72 {
73  // Cross section
74  if (fTableData)
75  delete fTableData;
76 }
77 
78 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
79 
81  const G4DataVector& /*cuts*/)
82 {
83 
84  if (verboseLevel > 3)
85  {
86  G4cout << "Calling G4DNABornExcitationModel1::Initialise()" << G4endl;
87  }
88 
89  if(fParticleDefinition != 0 && fParticleDefinition != particle)
90  {
91  G4Exception("G4DNABornExcitationModel1::Initialise","em0001",
92  FatalException,"Model already initialized for another particle type.");
93  }
94 
95  fParticleDefinition = particle;
96 
97  if(particle->GetParticleName() == "e-")
98  {
99  fTableFile = "dna/sigma_excitation_e_born";
100  fLowEnergy = 9*eV;
101  fHighEnergy = 1*MeV;
102  }
103  else if(particle->GetParticleName() == "proton")
104  {
105  fTableFile = "dna/sigma_excitation_p_born";
106  fLowEnergy = 500. * keV;
107  fHighEnergy = 100. * MeV;
108  }
109 
112 
113  G4double scaleFactor = (1.e-22 / 3.343) * m*m;
116 
117  if( verboseLevel>0 )
118  {
119  G4cout << "Born excitation model is initialized " << G4endl
120  << "Energy range: "
121  << LowEnergyLimit() / eV << " eV - "
122  << HighEnergyLimit() / keV << " keV for "
123  << particle->GetParticleName()
124  << G4endl;
125  }
126 
127  // Initialize water density pointer
129 
130  if (isInitialised)
131  { return;}
133  isInitialised = true;
134 }
135 
136 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
137 
139  const G4ParticleDefinition* particleDefinition,
140  G4double ekin,
141  G4double,
142  G4double)
143 {
144  if (verboseLevel > 3)
145  {
146  G4cout << "Calling CrossSectionPerVolume() of G4DNABornExcitationModel1"
147  << G4endl;
148  }
149 
150  if(particleDefinition != fParticleDefinition) return 0;
151 
152  // Calculate total cross section for model
153 
154  G4double sigma=0;
155 
156  G4double waterDensity = (*fpMolWaterDensity)[material->GetIndex()];
157 
158  if(waterDensity!= 0.0)
159  {
160  if (ekin >= fLowEnergy && ekin < fHighEnergy)
161  {
162  sigma = fTableData->FindValue(ekin);
163  }
164 
165  if (verboseLevel > 2)
166  {
167  G4cout << "__________________________________" << G4endl;
168  G4cout << "G4DNABornExcitationModel1 - XS INFO START" << G4endl;
169  G4cout << "Kinetic energy(eV)=" << ekin/eV << " particle : " << particleDefinition->GetParticleName() << G4endl;
170  G4cout << "Cross section per water molecule (cm^2)=" << sigma/cm/cm << G4endl;
171  G4cout << "Cross section per water molecule (cm^-1)=" << sigma*waterDensity/(1./cm) << G4endl;
172  G4cout << "G4DNABornExcitationModel1 - XS INFO END" << G4endl;
173  }
174  } // if (waterMaterial)
175 
176  return sigma*waterDensity;
177 }
178 
179 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
180 
181 void G4DNABornExcitationModel1::SampleSecondaries(std::vector<G4DynamicParticle*>* /*fvect*/,
182  const G4MaterialCutsCouple* /*couple*/,
183  const G4DynamicParticle* aDynamicParticle,
184  G4double,
185  G4double)
186 {
187 
188  if (verboseLevel > 3)
189  {
190  G4cout << "Calling SampleSecondaries() of G4DNABornExcitationModel1"
191  << G4endl;
192  }
193 
194  G4double k = aDynamicParticle->GetKineticEnergy();
195 
196  G4int level = RandomSelect(k);
197  G4double excitationEnergy = waterStructure.ExcitationEnergy(level);
198  G4double newEnergy = k - excitationEnergy;
199 
200  if (newEnergy > 0)
201  {
203 
206 
208  }
209 
210  const G4Track * theIncomingTrack = fParticleChangeForGamma->GetCurrentTrack();
212  level,
213  theIncomingTrack);
214 }
215 
216 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
217 
219  G4int level,
220  const G4ParticleDefinition* particle,
221  G4double kineticEnergy)
222 {
223  if (fParticleDefinition != particle)
224  {
225  G4Exception("G4DNABornExcitationModel1::GetPartialCrossSection",
226  "bornParticleType",
228  "Model initialized for another particle type.");
229  }
230 
231  return fTableData->GetComponent(level)->FindValue(kineticEnergy);
232 }
233 
234 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
235 
237 {
238  G4int level = 0;
239 
240  G4double* valuesBuffer = new G4double[fTableData->NumberOfComponents()];
241  const size_t n(fTableData->NumberOfComponents());
242  size_t i(n);
243  G4double value = 0.;
244 
245  while (i > 0)
246  {
247  i--;
248  valuesBuffer[i] = fTableData->GetComponent(i)->FindValue(k);
249  value += valuesBuffer[i];
250  }
251 
252  value *= G4UniformRand();
253  i = n;
254 
255  while (i > 0)
256  {
257  i--;
258 
259  if (valuesBuffer[i] > value)
260  {
261  delete[] valuesBuffer;
262  return i;
263  }
264  value -= valuesBuffer[i];
265  }
266 
267  if (valuesBuffer)
268  delete[] valuesBuffer;
269 
270  return level;
271 }
272 
static const double cm
Definition: G4SIunits.hh:118
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:641
static const double MeV
Definition: G4SIunits.hh:211
virtual G4double FindValue(G4double x, G4int componentId=0) const =0
G4double GetKineticEnergy() const
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:634
const std::vector< G4double > * fpMolWaterDensity
virtual const G4VEMDataSet * GetComponent(G4int componentId) const
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:604
size_t GetIndex() const
Definition: G4Material.hh:262
virtual G4bool LoadData(const G4String &argFileName)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &=*(new G4DataVector()))
int G4int
Definition: G4Types.hh:78
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
const G4String & GetParticleName() const
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:725
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
const std::vector< double > * GetNumMolPerVolTableFor(const G4Material *) const
virtual G4double GetPartialCrossSection(const G4Material *, G4int level, const G4ParticleDefinition *, G4double kineticEnergy)
const G4ThreeVector & GetMomentumDirection() const
const G4int n
virtual G4double FindValue(G4double e, G4int componentId=0) const
virtual size_t NumberOfComponents(void) const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4DNAWaterExcitationStructure waterStructure
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.
G4DNABornExcitationModel1(const G4ParticleDefinition *p=0, const G4String &nam="DNABornExcitationModel")
const G4Track * GetCurrentTrack() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
#define G4endl
Definition: G4ios.hh:61
static const double m
Definition: G4SIunits.hh:128
static const double keV
Definition: G4SIunits.hh:213
G4DNACrossSectionDataSet * fTableData
G4ParticleChangeForGamma * fParticleChangeForGamma
double G4double
Definition: G4Types.hh:76
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:732
const G4ParticleDefinition * fParticleDefinition
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:134
virtual G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)