Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4DNACPA100ExcitationModel.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 // CPA100 excitation model class for electrons
27 //
28 // Based on the work of M. Terrissol and M. C. Bordage
29 //
30 // Users are requested to cite the following papers:
31 // - M. Terrissol, A. Baudre, Radiat. Prot. Dosim. 31 (1990) 175-177
32 // - M.C. Bordage, J. Bordes, S. Edel, M. Terrissol, X. Franceries,
33 // M. Bardies, N. Lampe, S. Incerti, Phys. Med. 32 (2016) 1833-1840
34 //
35 // Authors of this class:
36 // M.C. Bordage, M. Terrissol, S. Edel, J. Bordes, S. Incerti
37 //
38 // 15.01.2014: creation
39 //
40 
42 #include "G4SystemOfUnits.hh"
43 #include "G4PhysicalConstants.hh"
44 #include "G4DNAChemistryManager.hh"
46 
47 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
48 
49 using namespace std;
50 
51 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
52 
54  const G4String& nam)
55 :G4VEmModel(nam),isInitialised(false)
56 {
58 
59  verboseLevel= 0;
60  // Verbosity scale:
61  // 0 = nothing
62  // 1 = warning for energy non-conservation
63  // 2 = details of energy budget
64  // 3 = calculation of cross sections, file openings, sampling of atoms
65  // 4 = entering in methods
66 
67  if( verboseLevel>0 )
68  {
69  G4cout << "CPA100 excitation model is constructed " << G4endl;
70  }
72 
73  // Selection of stationary mode
74 
75  statCode = false;
76 }
77 
78 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
79 
81 {
82  // Cross section
83 
84  std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
85  for (pos = tableData.begin(); pos != tableData.end(); ++pos)
86  {
87  G4DNACrossSectionDataSet* table = pos->second;
88  delete table;
89  }
90 
91 }
92 
93 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
94 
96  const G4DataVector& /*cuts*/)
97 {
98 
99  if (verboseLevel > 3)
100  G4cout << "Calling G4DNACPA100ExcitationModel::Initialise()" << G4endl;
101 
102  G4String fileElectron("dna/sigma_excitation_e_cpa100");
103 
105 
107 
108  G4double scaleFactor = 1.e-20 *m*m;
109 
110  // *** ELECTRON
111 
112  electron = electronDef->GetParticleName();
113 
114  tableFile[electron] = fileElectron;
115 
116  lowEnergyLimit[electron] = 11 * eV; // Default low enetgy limit
117  //lowEnergyLimit[electron] = 10.481 * eV; // From sigexi2 file by MCB
118  highEnergyLimit[electron] = 255955 * eV; // idem
119 
120  // Cross section
121 
122  G4DNACrossSectionDataSet* tableE
123  = new G4DNACrossSectionDataSet(new G4LogLogInterpolation, eV, scaleFactor );
124 
125  /*
126  G4DNACrossSectionDataSet* tableE =
127  new G4DNACrossSectionDataSet(new G4DNACPA100LogLogInterpolation, eV, scaleFactor );
128  */
129 
130  tableE->LoadData(fileElectron);
131 
132  tableData[electron] = tableE;
133 
134  //
135 
136  if (particle==electronDef)
137  {
140  }
141 
142 
143 // if( verboseLevel>0 )
144  {
145  G4cout << "CPA100 excitation model is initialized " << G4endl
146  << "Energy range: "
147  << LowEnergyLimit() / eV << " eV - "
148  << HighEnergyLimit() / keV << " keV for "
149  << particle->GetParticleName()
150  << G4endl;
151  }
152 
153  // Initialize water density pointer
156 
157  if (isInitialised) { return; }
159  isInitialised = true;
160 }
161 
162 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
163 
165  const G4ParticleDefinition* particleDefinition,
166  G4double ekin,
167  G4double,
168  G4double)
169 {
170 
171  if (verboseLevel > 3)
172  G4cout << "Calling CrossSectionPerVolume() of G4DNACPA100ExcitationModel" << G4endl;
173 
174  if (particleDefinition != G4Electron::ElectronDefinition()) return 0;
175 
176  // Calculate total cross section for model
177 
178  G4double lowLim = 0;
179  G4double highLim = 0;
180  G4double sigma=0;
181 
182  G4double waterDensity = (*fpMolWaterDensity)[material->GetIndex()];
183 
184  if(waterDensity!= 0.0)
185  {
186  const G4String& particleName = particleDefinition->GetParticleName();
187 
188  std::map< G4String,G4double,std::less<G4String> >::iterator pos1;
189  pos1 = lowEnergyLimit.find(particleName);
190  if (pos1 != lowEnergyLimit.end())
191  {
192  lowLim = pos1->second;
193  }
194 
195  std::map< G4String,G4double,std::less<G4String> >::iterator pos2;
196  pos2 = highEnergyLimit.find(particleName);
197  if (pos2 != highEnergyLimit.end())
198  {
199  highLim = pos2->second;
200  }
201 
202  if (ekin > lowLim && ekin < highLim)
203  {
204  std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
205  pos = tableData.find(particleName);
206 
207  if (pos != tableData.end())
208  {
209  G4DNACrossSectionDataSet* table = pos->second;
210  if (table != 0)
211  {
212  sigma = table->FindValue(ekin);
213  }
214  }
215  else
216  {
217  G4Exception("G4DNACPA100ExcitationModel::CrossSectionPerVolume","em0002",
218  FatalException,"Model not applicable to particle type.");
219  }
220  }
221 
222  if (verboseLevel > 2)
223  {
224  G4cout << "__________________________________" << G4endl;
225  G4cout << "G4DNACPA100ExcitationModel - XS INFO START" << G4endl;
226  G4cout << "Kinetic energy(eV)=" << ekin/eV << " particle : " << particleName << G4endl;
227  G4cout << "Cross section per water molecule (cm^2)=" << sigma/cm/cm << G4endl;
228  G4cout << "Cross section per water molecule (cm^-1)=" << sigma*waterDensity/(1./cm) << G4endl;
229  // G4cout << " - Cross section per water molecule (cm^-1)=" << sigma*material->GetAtomicNumDensityVector()[1]/(1./cm) << G4endl;
230  G4cout << "G4DNACPA100ExcitationModel - XS INFO END" << G4endl;
231  }
232 
233  } // if (waterMaterial)
234 
235  return sigma*waterDensity;
236 
237 }
238 
239 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
240 
241 void G4DNACPA100ExcitationModel::SampleSecondaries(std::vector<G4DynamicParticle*>* ,
242  const G4MaterialCutsCouple*,
243  const G4DynamicParticle* aDynamicParticle,
244  G4double,
245  G4double)
246 {
247 
248  if (verboseLevel > 3)
249  G4cout << "Calling SampleSecondaries() of G4DNACPA100ExcitationModel" << G4endl;
250 
251  G4double k = aDynamicParticle->GetKineticEnergy();
252 
253  const G4String& particleName = aDynamicParticle->GetDefinition()->GetParticleName();
254 
255  G4int level = RandomSelect(k,particleName);
256  G4double excitationEnergy = waterStructure.ExcitationEnergy(level);
257  G4double newEnergy = k - excitationEnergy;
258 
259  if (newEnergy > 0)
260  {
261  // fParticleChangeForGamma->ProposeMomentumDirection(aDynamicParticle->GetMomentumDirection());
262 
263  // We take into account direction change as described page 87 (II.92) in thesis by S. Edel
264 
265  G4double cosTheta =
266 
267  (excitationEnergy/k) / (1. + (k/(2*electron_mass_c2))*(1.-excitationEnergy/k) );
268 
269  cosTheta = std::sqrt(1.-cosTheta);
270 
271  G4double phi = 2. * pi * G4UniformRand();
272 
273  G4ThreeVector zVers = aDynamicParticle->GetMomentumDirection();
274 
275  //G4ThreeVector xVers = zVers.orthogonal();
276  //G4ThreeVector yVers = zVers.cross(xVers);
277  //G4double xDir = std::sqrt(1. - cosTheta*cosTheta);
278  //G4double yDir = xDir;
279  //xDir *= std::cos(phi);
280  //yDir *= std::sin(phi);
281  // G4ThreeVector zPrimeVers((xDir*xVers + yDir*yVers + cosTheta*zVers));
282 
283  // Computation of scattering angles (from Subroutine DIRAN in CPA100)
284 
285  G4double CT1, ST1, CF1, SF1, CT2, ST2, CF2, SF2;
286  G4double sinTheta = std::sqrt (1-cosTheta*cosTheta);
287 
288  CT1=0;
289  ST1=0;
290  CF1=0;
291  SF1=0;
292  CT2=0;
293  ST2=0;
294  CF2=0;
295  SF2=0;
296 
297  CT1 = zVers.z();
298  ST1=std::sqrt(1.-CT1*CT1);
299 
300  if (ST1!=0) CF1 = zVers.x()/ST1; else CF1 = std::cos(2. * pi * G4UniformRand());
301  if (ST1!=0) SF1 = zVers.y()/ST1; else SF1 = std::sqrt(1.-CF1*CF1);
302 
303  G4double A3, A4, A5, A2, A1;
304  A3=0;
305  A4=0;
306  A5=0;
307  A2=0;
308  A1=0;
309 
310  A3 = sinTheta*std::cos(phi);
311  A4 = A3*CT1 + ST1*cosTheta;
312  A5 = sinTheta * std::sin(phi);
313  A2 = A4 * SF1 + A5 * CF1;
314  A1 = A4 * CF1 - A5 * SF1;
315 
316  CT2 = CT1*cosTheta - ST1*A3;
317  ST2 = std::sqrt(1.-CT2*CT2);
318 
319  if (ST2==0) ST2=1E-6;
320  CF2 = A1/ST2;
321  SF2 = A2/ST2;
322 
323  /*
324  G4cout << "CT1=" << CT1 << G4endl;
325  G4cout << "ST1=" << ST1 << G4endl;
326  G4cout << "CF1=" << CF1 << G4endl;
327  G4cout << "SF1=" << SF1 << G4endl;
328  G4cout << "cosTheta=" << cosTheta << G4endl;
329  G4cout << "sinTheta=" << sinTheta << G4endl;
330  G4cout << "cosPhi=" << std::cos(phi) << G4endl;
331  G4cout << "sinPhi=" << std::sin(phi) << G4endl;
332  G4cout << "CT2=" << CT2 << G4endl;
333  G4cout << "ST2=" << ST2 << G4endl;
334  G4cout << "CF2=" << CF2 << G4endl;
335  G4cout << "SF2=" << SF2 << G4endl;
336  */
337 
338  G4ThreeVector zPrimeVers(ST2*CF2,ST2*SF2,CT2);
339 
340  //
341 
343 
344  //
345 
348 
350  }
351 
352  // Chemistry
353 
354  const G4Track * theIncomingTrack = fParticleChangeForGamma->GetCurrentTrack();
356  level,
357  theIncomingTrack);
358 }
359 
360 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
361 
363 {
364  G4int level = 0;
365 
366  std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
367  pos = tableData.find(particle);
368 
369  if (pos != tableData.end())
370  {
371  G4DNACrossSectionDataSet* table = pos->second;
372 
373  if (table != 0)
374  {
375  G4double* valuesBuffer = new G4double[table->NumberOfComponents()];
376  const size_t n(table->NumberOfComponents());
377  size_t i(n);
378  G4double value = 0.;
379 
380  //Verification
381  /*
382  G4double tmp=10.481*eV;
383  G4cout << table->GetComponent(0)->FindValue(tmp)/(1e-20*m*m) << G4endl;
384  G4cout << table->GetComponent(1)->FindValue(tmp)/(1e-20*m*m) << G4endl;
385  G4cout << table->GetComponent(2)->FindValue(tmp)/(1e-20*m*m) << G4endl;
386  G4cout << table->GetComponent(3)->FindValue(tmp)/(1e-20*m*m) << G4endl;
387  G4cout << table->GetComponent(4)->FindValue(tmp)/(1e-20*m*m) << G4endl;
388  G4cout <<
389  table->GetComponent(0)->FindValue(tmp)/(1e-20*m*m) +
390  table->GetComponent(1)->FindValue(tmp)/(1e-20*m*m) +
391  table->GetComponent(2)->FindValue(tmp)/(1e-20*m*m) +
392  table->GetComponent(3)->FindValue(tmp)/(1e-20*m*m) +
393  table->GetComponent(4)->FindValue(tmp)/(1e-20*m*m)
394  << G4endl;
395  abort();
396  */
397  //
398  //Dump
399  //
400  /*
401  G4double minEnergy = 10.481 * eV;
402  G4double maxEnergy = 255955. * eV;
403  G4int nEnergySteps = 1000;
404  G4double energy(minEnergy);
405  G4double stpEnergy(std::pow(maxEnergy/energy, 1./static_cast<G4double>(nEnergySteps-1)));
406  G4int step(nEnergySteps);
407  system ("rm -rf excitation-cap100.out");
408  FILE* myFile=fopen("excitation-cpa100.out","a");
409  while (step>0)
410  {
411  step--;
412  fprintf (myFile,"%16.9le %16.9le %16.9le %16.9le %16.9le %16.9le %16.9le \n",
413  energy/eV,
414  table->GetComponent(0)->FindValue(energy)/(1e-20*m*m),
415  table->GetComponent(1)->FindValue(energy)/(1e-20*m*m),
416  table->GetComponent(2)->FindValue(energy)/(1e-20*m*m),
417  table->GetComponent(3)->FindValue(energy)/(1e-20*m*m),
418  table->GetComponent(4)->FindValue(energy)/(1e-20*m*m),
419  table->GetComponent(0)->FindValue(energy)/(1e-20*m*m)+
420  table->GetComponent(1)->FindValue(energy)/(1e-20*m*m)+
421  table->GetComponent(2)->FindValue(energy)/(1e-20*m*m)+
422  table->GetComponent(3)->FindValue(energy)/(1e-20*m*m)+
423  table->GetComponent(4)->FindValue(energy)/(1e-20*m*m)
424  );
425  energy*=stpEnergy;
426  }
427  fclose (myFile);
428  abort();
429  */
430  //
431  // end of dump
432  //
433 
434  while (i>0)
435  {
436  i--;
437  valuesBuffer[i] = table->GetComponent(i)->FindValue(k);
438  value += valuesBuffer[i];
439  }
440 
441  value *= G4UniformRand();
442 
443  i = n;
444 
445  while (i > 0)
446  {
447  i--;
448 
449  if (valuesBuffer[i] > value)
450  {
451  delete[] valuesBuffer;
452  return i;
453  }
454  value -= valuesBuffer[i];
455  }
456 
457  if (valuesBuffer) delete[] valuesBuffer;
458 
459  }
460  }
461  else
462  {
463  G4Exception("G4DNACPA100ExcitationModel::RandomSelect","em0002",
464  FatalException,"Model not applicable to particle type.");
465  }
466  return level;
467 }
468 
469 
virtual G4bool LoadData(const G4String &argFileName)
static const G4double pos
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:707
virtual size_t NumberOfComponents(void) const
G4DNACPA100WaterExcitationStructure waterStructure
G4DNACPA100ExcitationModel(const G4ParticleDefinition *p=0, const G4String &nam="DNACPA100ExcitationModel")
static constexpr double keV
Definition: G4SIunits.hh:216
size_t GetIndex() const
Definition: G4Material.hh:262
const std::vector< G4double > * fpMolWaterDensity
#define G4endl
Definition: G4ios.hh:61
const G4ThreeVector & GetMomentumDirection() const
const G4String & GetParticleName() const
G4ParticleChangeForGamma * fParticleChangeForGamma
void SetProposedKineticEnergy(G4double proposedKinEnergy)
double z() const
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:608
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:616
const std::vector< double > * GetNumMolPerVolTableFor(const G4Material *) const
Retrieve a table of molecular densities (number of molecules per unit volume) in the G4 unit system f...
std::map< G4String, G4double, std::less< G4String > > highEnergyLimit
static constexpr double m
Definition: G4SIunits.hh:129
double G4double
Definition: G4Types.hh:76
G4ParticleDefinition * GetDefinition() const
const XML_Char int const XML_Char * value
Definition: expat.h:331
static constexpr double electron_mass_c2
static G4DNAMolecularMaterial * Instance()
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
#define G4UniformRand()
Definition: Randomize.hh:53
void CreateWaterMolecule(ElectronicModification, G4int, const G4Track *)
static constexpr double eV
Definition: G4SIunits.hh:215
const G4Track * GetCurrentTrack() const
G4int RandomSelect(G4double energy, const G4String &particle)
std::map< G4String, G4double, std::less< G4String > > lowEnergyLimit
virtual G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:89
Hep3Vector unit() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.hh:65
int G4int
Definition: G4Types.hh:78
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
double x() const
virtual const G4VEMDataSet * GetComponent(G4int componentId) const
Char_t n[5]
static constexpr double pi
Definition: G4SIunits.hh:75
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:131
double y() const
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &=*(new G4DataVector()))
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:609
static G4DNAChemistryManager * Instance()
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
virtual G4double FindValue(G4double x, G4int componentId=0) const =0