source: trunk/source/processes/electromagnetic/lowenergy/src/G4DNADingfelderChargeIncreaseModel.cc @ 968

Last change on this file since 968 was 968, checked in by garnier, 15 years ago

fichier ajoutes

File size: 17.4 KB
Line 
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: G4DNADingfelderChargeIncreaseModel.cc,v 1.3 2009/02/16 11:00:11 sincerti Exp $
27// GEANT4 tag $Name: geant4-09-02-ref-02 $
28//
29
30#include "G4DNADingfelderChargeIncreaseModel.hh"
31
32//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
33
34using namespace std;
35
36//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
37
38G4DNADingfelderChargeIncreaseModel::G4DNADingfelderChargeIncreaseModel(const G4ParticleDefinition*,
39                                             const G4String& nam)
40:G4VEmModel(nam),isInitialised(false)
41{
42
43  verboseLevel= 0;
44  // Verbosity scale:
45  // 0 = nothing
46  // 1 = warning for energy non-conservation
47  // 2 = details of energy budget
48  // 3 = calculation of cross sections, file openings, sampling of atoms
49  // 4 = entering in methods
50 
51  G4cout << "Dingfelder charge increase model is constructed " << G4endl;
52
53}
54
55//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
56
57G4DNADingfelderChargeIncreaseModel::~G4DNADingfelderChargeIncreaseModel()
58{}
59
60//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
61
62void G4DNADingfelderChargeIncreaseModel::Initialise(const G4ParticleDefinition* particle,
63                                       const G4DataVector& /*cuts*/)
64{
65
66  if (verboseLevel > 3)
67    G4cout << "Calling G4DNADingfelderChargeIncreaseModel::Initialise()" << G4endl;
68
69  // Energy limits
70 
71  G4DNAGenericIonsManager *instance;
72  instance = G4DNAGenericIonsManager::Instance();
73  G4ParticleDefinition* hydrogenDef = instance->GetIon("hydrogen");
74  G4ParticleDefinition* alphaPlusDef = instance->GetIon("alpha+");
75  G4ParticleDefinition* heliumDef = instance->GetIon("helium");
76
77  G4String hydrogen;
78  G4String alphaPlus;
79  G4String helium;
80
81  if (hydrogenDef != 0)
82  {
83    hydrogen = hydrogenDef->GetParticleName();
84    lowEnergyLimit[hydrogen] = 1. * keV;
85    highEnergyLimit[hydrogen] = 10. * MeV;
86  }
87  else
88  {
89    G4Exception("G4DNADingfelderChargeIncreaseModel::Initialise: hydrogen is not defined");
90  }
91
92  if (alphaPlusDef != 0)
93  {
94    alphaPlus = alphaPlusDef->GetParticleName();
95    lowEnergyLimit[alphaPlus] = 1. * keV;
96    highEnergyLimit[alphaPlus] = 10. * MeV;
97  }
98  else
99  {
100    G4Exception("G4DNADingfelderChargeIncreaseModel::Initialise: alphaPlus is not defined");
101  }
102
103  if (heliumDef != 0)
104  {
105    helium = heliumDef->GetParticleName();
106    lowEnergyLimit[helium] = 1. * keV;
107    highEnergyLimit[helium] = 10. * MeV;
108  }
109  else
110  {
111    G4Exception("G4DNADingfelderChargeIncreaseModel::Initialise: helium is not defined");
112  }
113 
114  if (particle==hydrogenDef) 
115  {
116    SetLowEnergyLimit(lowEnergyLimit[hydrogen]);
117    SetHighEnergyLimit(highEnergyLimit[hydrogen]);
118  }
119
120  if (particle==alphaPlusDef) 
121  {
122    SetLowEnergyLimit(lowEnergyLimit[alphaPlus]);
123    SetHighEnergyLimit(highEnergyLimit[alphaPlus]);
124  }
125
126  if (particle==heliumDef) 
127  {
128    SetLowEnergyLimit(lowEnergyLimit[helium]);
129    SetHighEnergyLimit(highEnergyLimit[helium]);
130  }
131
132  // Final state
133 
134  //ALPHA+
135
136  f0[0][0]=1.;
137  a0[0][0]=2.25;
138  a1[0][0]=-0.75;
139  b0[0][0]=-32.10;
140  c0[0][0]=0.600;
141  d0[0][0]=2.40;
142  x0[0][0]=4.60;
143
144  x1[0][0]=-1.;
145  b1[0][0]=-1.;
146   
147  numberOfPartialCrossSections[0]=1;
148
149  //HELIUM
150
151  f0[0][1]=1.;
152  a0[0][1]=2.25;
153  a1[0][1]=-0.75;
154  b0[0][1]=-30.93;
155  c0[0][1]=0.590;
156  d0[0][1]=2.35;
157  x0[0][1]=4.29;
158
159  f0[1][1]=1.;
160  a0[1][1]=2.25;
161  a1[1][1]=-0.75;
162  b0[1][1]=-32.61;
163  c0[1][1]=0.435;
164  d0[1][1]=2.70;
165  x0[1][1]=4.45;
166
167  x1[0][1]=-1.;
168  b1[0][1]=-1.;
169
170  x1[1][1]=-1.;
171  b1[1][1]=-1.;
172
173  numberOfPartialCrossSections[1]=2;
174 
175  //
176 
177  G4cout << "Dingfelder charge increase model is initialized " << G4endl
178         << "Energy range: "
179         << LowEnergyLimit() / keV << " keV - "
180         << HighEnergyLimit() / MeV << " MeV for "
181         << particle->GetParticleName()
182         << G4endl;
183
184  if(!isInitialised) 
185  {
186    isInitialised = true;
187 
188    if(pParticleChange)
189      fParticleChangeForGamma = reinterpret_cast<G4ParticleChangeForGamma*>(pParticleChange);
190    else
191      fParticleChangeForGamma = new G4ParticleChangeForGamma();
192  }   
193
194  // InitialiseElementSelectors(particle,cuts);
195 
196  // Test if water material
197
198  flagMaterialIsWater= false;
199  densityWater = 0;
200
201  const G4ProductionCutsTable* theCoupleTable = G4ProductionCutsTable::GetProductionCutsTable();
202
203  if(theCoupleTable) 
204  {
205    G4int numOfCouples = theCoupleTable->GetTableSize();
206 
207    if(numOfCouples>0) 
208    {
209          for (G4int i=0; i<numOfCouples; i++) 
210          {
211            const G4MaterialCutsCouple* couple = theCoupleTable->GetMaterialCutsCouple(i);
212            const G4Material* material = couple->GetMaterial();
213
214            size_t j = material->GetNumberOfElements();
215            while (j>0)
216            {
217               j--;
218               const G4Element* element(material->GetElement(j));
219               if (element->GetZ() == 8.)
220               {
221                  G4double density = material->GetAtomicNumDensityVector()[j];
222                  if (density > 0.) 
223                  { 
224                    flagMaterialIsWater = true; 
225                    densityWater = density; 
226                   
227                    if (verboseLevel > 3) 
228                    G4cout << "Water material is found with density(cm^-3)=" << density/(cm*cm*cm) << G4endl;
229                  }
230               }
231            }
232 
233          }
234   } // if(numOfCouples>0)
235
236  } // if (theCoupleTable)
237
238}
239
240//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
241
242G4double G4DNADingfelderChargeIncreaseModel::CrossSectionPerVolume(const G4Material*,
243                                           const G4ParticleDefinition* particleDefinition,
244                                           G4double k,
245                                           G4double,
246                                           G4double)
247{
248  if (verboseLevel > 3)
249    G4cout << "Calling CrossSectionPerVolume() of G4DNADingfelderChargeIncreaseModel" << G4endl;
250
251 // Calculate total cross section for model
252
253  G4DNAGenericIonsManager *instance;
254  instance = G4DNAGenericIonsManager::Instance();
255
256  if (
257      particleDefinition != instance->GetIon("hydrogen")
258      &&
259      particleDefinition != instance->GetIon("alpha+")
260      &&
261      particleDefinition != instance->GetIon("helium")
262      )
263                   
264    return 0;
265
266  G4double lowLim = 0;
267  G4double highLim = 0;
268  G4double totalCrossSection = 0.;
269 
270  if (flagMaterialIsWater)
271  {
272    const G4String& particleName = particleDefinition->GetParticleName();
273
274    std::map< G4String,G4double,std::less<G4String> >::iterator pos1;
275    pos1 = lowEnergyLimit.find(particleName);
276
277    if (pos1 != lowEnergyLimit.end())
278    {
279      lowLim = pos1->second;
280    }
281
282    std::map< G4String,G4double,std::less<G4String> >::iterator pos2;
283    pos2 = highEnergyLimit.find(particleName);
284
285    if (pos2 != highEnergyLimit.end())
286    {
287      highLim = pos2->second;
288    }
289
290    if (k >= lowLim && k < highLim)
291    {
292      //HYDROGEN
293      if (particleDefinition == instance->GetIon("hydrogen"))
294      {
295          const  G4double aa = 2.835;
296          const  G4double bb = 0.310;
297          const  G4double cc = 2.100;
298          const  G4double dd = 0.760;
299          const  G4double fac = 1.0e-18;
300          const  G4double rr = 13.606 * eV; 
301         
302          G4double t = k / (proton_mass_c2/electron_mass_c2); 
303          G4double x = t / rr; 
304          G4double temp = 4.0 * pi * Bohr_radius/nm * Bohr_radius/nm * fac;
305          G4double sigmal =  temp * cc * (std::pow(x,dd));
306          G4double sigmah = temp * (aa * std::log(1.0 + x) + bb) / x;
307          totalCrossSection = 1.0/(1.0/sigmal + 1.0/sigmah) *m*m;
308      }
309      else
310      {
311          totalCrossSection = Sum(k,particleDefinition);
312      }
313    }
314 
315    if (verboseLevel > 3)
316    {
317      G4cout << "---> Kinetic energy(eV)=" << k/eV << G4endl;
318      G4cout << " - Cross section per water molecule (cm^2)=" << totalCrossSection/cm/cm << G4endl;
319      G4cout << " - Cross section per water molecule (cm^-1)=" << totalCrossSection*densityWater/(1./cm) << G4endl;
320    } 
321
322  } // if (flagMaterialIsWater)
323
324  return totalCrossSection*densityWater;                   
325
326}
327
328//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
329
330void G4DNADingfelderChargeIncreaseModel::SampleSecondaries(std::vector<G4DynamicParticle*>* fvect,
331                                              const G4MaterialCutsCouple* /*couple*/,
332                                              const G4DynamicParticle* aDynamicParticle,
333                                              G4double,
334                                              G4double)
335{
336  if (verboseLevel > 3)
337    G4cout << "Calling SampleSecondaries() of G4DNADingfelderChargeIncreaseModel" << G4endl;
338
339  fParticleChangeForGamma->ProposeTrackStatus(fStopAndKill);
340  fParticleChangeForGamma->ProposeLocalEnergyDeposit(0.);
341
342  G4ParticleDefinition* definition = aDynamicParticle->GetDefinition();
343 
344  G4double inK = aDynamicParticle->GetKineticEnergy();
345
346  G4int finalStateIndex = RandomSelect(inK,definition);
347
348  G4int n = NumberOfFinalStates(definition,finalStateIndex);
349 
350  G4double outK = inK - IncomingParticleBindingEnergyConstant(definition,finalStateIndex);
351
352  G4DNAGenericIonsManager* instance;
353  instance = G4DNAGenericIonsManager::Instance();
354
355  G4double electronK;
356  if (definition == instance->GetIon("hydrogen")) electronK = inK*electron_mass_c2/proton_mass_c2;
357  else electronK = inK*electron_mass_c2/(3728*MeV);
358 
359  if (outK<0)
360  {
361    G4String message;
362    message="G4DNADingfelderChargeIncreaseModel::SampleSecondaries: final kinetic energy is below 0! Process ";
363  }
364 
365  G4DynamicParticle* dp = new G4DynamicParticle
366 
367  (OutgoingParticleDefinition(definition,finalStateIndex), 
368   aDynamicParticle->GetMomentumDirection(), 
369   outK) ;
370
371  fvect->push_back(dp);
372
373  n = n - 1;
374 
375  while (n>0)
376  {
377    n--;
378    fvect->push_back(new G4DynamicParticle
379    (G4Electron::Electron(), aDynamicParticle->GetMomentumDirection(), electronK) );
380  }
381       
382}
383
384//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
385
386G4int G4DNADingfelderChargeIncreaseModel::NumberOfFinalStates(G4ParticleDefinition* particleDefinition, 
387                                                      G4int finalStateIndex )
388 
389{
390  G4DNAGenericIonsManager* instance;
391  instance = G4DNAGenericIonsManager::Instance();
392  if (particleDefinition == instance->GetIon("hydrogen")) return 2;
393  if (particleDefinition == instance->GetIon("alpha+")) return 2;
394 
395  if (particleDefinition == instance->GetIon("helium")) 
396  {    if (finalStateIndex==0) return 2;
397       return 3;
398  }
399  return 0;
400}
401
402//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
403
404G4ParticleDefinition* G4DNADingfelderChargeIncreaseModel::OutgoingParticleDefinition (G4ParticleDefinition* particleDefinition, 
405                                                                              G4int finalStateIndex) 
406{
407  G4DNAGenericIonsManager * instance(G4DNAGenericIonsManager::Instance());
408  if (particleDefinition == instance->GetIon("hydrogen")) return G4Proton::Proton(); 
409  if (particleDefinition == instance->GetIon("alpha+")) return instance->GetIon("alpha++");
410
411  if (particleDefinition == instance->GetIon("alpha+")) return instance->GetIon("helium");
412  {
413    if (finalStateIndex==0) return instance->GetIon("alpha+"); 
414    return instance->GetIon("alpha++");
415  }
416  return 0;
417}
418
419//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
420
421G4double G4DNADingfelderChargeIncreaseModel::IncomingParticleBindingEnergyConstant(G4ParticleDefinition* particleDefinition, 
422                                                                           G4int finalStateIndex )
423{
424  G4DNAGenericIonsManager * instance(G4DNAGenericIonsManager::Instance());
425  if(particleDefinition == instance->GetIon("hydrogen")) return 13.6*eV;
426 
427  if(particleDefinition == instance->GetIon("alpha+"))
428  {
429      // Binding energy for    He+ -> He++ + e-    54.509 eV
430      // Binding energy for    He  -> He+  + e-    24.587 eV
431      return 54.509*eV;
432  }
433   
434  if(particleDefinition == instance->GetIon("helium"))
435  {
436      // Binding energy for    He+ -> He++ + e-    54.509 eV
437      // Binding energy for    He  -> He+  + e-    24.587 eV
438
439      if (finalStateIndex==0) return 24.587*eV;
440      return (54.509 + 24.587)*eV;
441  } 
442
443  return 0;
444}
445
446//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
447
448//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
449
450G4double G4DNADingfelderChargeIncreaseModel::PartialCrossSection(G4double k, G4int index, 
451                                                           const G4ParticleDefinition* particleDefinition)
452{
453
454  G4int particleTypeIndex = 0;
455  G4DNAGenericIonsManager *instance;
456  instance = G4DNAGenericIonsManager::Instance();
457
458  if (particleDefinition == instance->GetIon("alpha+")) particleTypeIndex=0;
459  if (particleDefinition == instance->GetIon("helium")) particleTypeIndex=1;
460
461  //
462  // sigma(T) = f0 10 ^ y(log10(T/eV))
463  //
464  //         /  a0 x + b0                    x < x0
465  //         |
466  // y(x) = <   a0 x + b0 - c0 (x - x0)^d0   x0 <= x < x1
467  //         |
468  //         \  a1 x + b1                    x >= x1
469  //
470  //
471  // f0, a0, a1, b0, b1, c0, d0, x0, x1 are parameters that change for protons and helium (0, +, ++)
472  //
473  // f0 has been added to the code in order to manage partial (shell-dependent) cross sections (if no shell dependence is present. f0=1. Sum of f0 over the considered shells should give 1)
474  //
475  // From Rad. Phys. and Chem. 59 (2000) 255-275, M. Dingfelder et al.
476  // Inelastic-collision cross sections of liquid water for interactions of energetic proton
477  //
478 
479  if (x1[index][particleTypeIndex]<x0[index][particleTypeIndex])
480  {
481      //
482      // if x1 < x0 means that x1 and b1 will be calculated with the following formula (this piece of code is run on all alphas and not on protons)
483      //
484      // x1 = x0 + ((a0 - a1)/(c0 * d0)) ^ (1 / (d0 - 1))
485      //
486      // b1 = (a0 - a1) * x1 + b0 - c0 * (x1 - x0) ^ d0
487      //
488 
489      x1[index][particleTypeIndex]=x0[index][particleTypeIndex] + std::pow((a0[index][particleTypeIndex] - a1[index][particleTypeIndex]) / (c0[index][particleTypeIndex] * d0[index][particleTypeIndex]), 1. / (d0[index][particleTypeIndex] - 1.));
490      b1[index][particleTypeIndex]=(a0[index][particleTypeIndex] - a1[index][particleTypeIndex]) * x1[index][particleTypeIndex] + b0[index][particleTypeIndex] - c0[index][particleTypeIndex] * std::pow(x1[index][particleTypeIndex] - x0[index][particleTypeIndex], d0[index][particleTypeIndex]);
491  }
492
493  G4double x(std::log10(k/eV));
494  G4double y;
495 
496  if (x<x0[index][particleTypeIndex])
497    y=a0[index][particleTypeIndex] * x + b0[index][particleTypeIndex];
498  else if (x<x1[index][particleTypeIndex])
499    y=a0[index][particleTypeIndex] * x + b0[index][particleTypeIndex] - c0[index][particleTypeIndex] * std::pow(x - x0[index][particleTypeIndex], d0[index][particleTypeIndex]);
500  else
501    y=a1[index][particleTypeIndex] * x + b1[index][particleTypeIndex];
502
503  return f0[index][particleTypeIndex] * std::pow(10., y)*m*m;
504 
505}
506
507//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
508
509G4int G4DNADingfelderChargeIncreaseModel::RandomSelect(G4double k, 
510                                                        const G4ParticleDefinition* particleDefinition)
511{
512
513  G4int particleTypeIndex = 0;
514  G4DNAGenericIonsManager *instance;
515  instance = G4DNAGenericIonsManager::Instance();
516
517  if (particleDefinition == instance->GetIon("hydrogen")) return 0;
518  if (particleDefinition == instance->GetIon("alpha+")) particleTypeIndex=0;
519  if (particleDefinition == instance->GetIon("helium")) particleTypeIndex=1;
520
521  const G4int n = numberOfPartialCrossSections[particleTypeIndex];
522  G4double* values(new G4double[n]);
523  G4double value = 0;
524  G4int i = n;
525 
526  while (i>0)
527  {
528      i--;
529      values[i]=PartialCrossSection(k, i, particleDefinition);
530      value+=values[i];
531  }
532 
533  value*=G4UniformRand();
534 
535  i=n;
536  while (i>0)
537  {
538      i--;
539   
540      if (values[i]>value)
541        break;
542 
543      value-=values[i];
544  }
545 
546  delete[] values;
547 
548  return i;
549}
550
551//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
552
553G4double G4DNADingfelderChargeIncreaseModel::Sum(G4double k, const G4ParticleDefinition* particleDefinition)
554{
555  G4int particleTypeIndex = 0;
556  G4DNAGenericIonsManager *instance;
557  instance = G4DNAGenericIonsManager::Instance();
558
559  if (particleDefinition == instance->GetIon("alpha+")) particleTypeIndex=0;
560  if (particleDefinition == instance->GetIon("helium")) particleTypeIndex=1;
561
562  G4double totalCrossSection = 0.;
563
564  for (G4int i=0; i<numberOfPartialCrossSections[particleTypeIndex]; i++)
565  {
566    totalCrossSection += PartialCrossSection(k,i,particleDefinition);
567  }
568  return totalCrossSection;
569}
Note: See TracBrowser for help on using the repository browser.