source: trunk/source/processes/electromagnetic/lowenergy/src/G4DNAMillerGreenExcitationModel.cc @ 1315

Last change on this file since 1315 was 1315, checked in by garnier, 14 years ago

update geant4-09-04-beta-cand-01 interfaces-V09-03-09 vis-V09-03-08

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: G4DNAMillerGreenExcitationModel.cc,v 1.9 2010/06/08 21:50:00 sincerti Exp $
27// GEANT4 tag $Name: geant4-09-04-beta-cand-01 $
28//
29
30#include "G4DNAMillerGreenExcitationModel.hh"
31
32//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
33
34using namespace std;
35
36//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
37
38G4DNAMillerGreenExcitationModel::G4DNAMillerGreenExcitationModel(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  if( verboseLevel>0 ) 
52  { 
53    G4cout << "Miller & Green excitation model is constructed " << G4endl;
54  }
55}
56
57//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
58
59G4DNAMillerGreenExcitationModel::~G4DNAMillerGreenExcitationModel()
60{}
61
62//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
63
64void G4DNAMillerGreenExcitationModel::Initialise(const G4ParticleDefinition* particle,
65                                       const G4DataVector& /*cuts*/)
66{
67
68  if (verboseLevel > 3)
69    G4cout << "Calling G4DNAMillerGreenExcitationModel::Initialise()" << G4endl;
70
71  // Energy limits
72 
73  G4DNAGenericIonsManager *instance;
74  instance = G4DNAGenericIonsManager::Instance();
75  G4ParticleDefinition* protonDef = G4Proton::ProtonDefinition();
76  G4ParticleDefinition* alphaPlusPlusDef = instance->GetIon("alpha++");
77  G4ParticleDefinition* alphaPlusDef = instance->GetIon("alpha+");
78  G4ParticleDefinition* heliumDef = instance->GetIon("helium");
79
80  G4String proton;
81  G4String alphaPlusPlus;
82  G4String alphaPlus;
83  G4String helium;
84
85  if (protonDef != 0)
86  {
87    proton = protonDef->GetParticleName();
88    lowEnergyLimit[proton] = 10. * eV;
89    highEnergyLimit[proton] = 500. * keV;
90   
91    kineticEnergyCorrection[0] = 1.;
92    slaterEffectiveCharge[0][0] = 0.;
93    slaterEffectiveCharge[1][0] = 0.;
94    slaterEffectiveCharge[2][0] = 0.;
95    sCoefficient[0][0] = 0.;
96    sCoefficient[1][0] = 0.;
97    sCoefficient[2][0] = 0.;
98  }
99  else
100  {
101    G4Exception("G4DNAMillerGreenExcitationModel::Initialise: proton is not defined");
102  }
103
104  if (alphaPlusPlusDef != 0)
105  {
106    alphaPlusPlus = alphaPlusPlusDef->GetParticleName();
107    lowEnergyLimit[alphaPlusPlus] = 1. * keV;
108    highEnergyLimit[alphaPlusPlus] = 10. * MeV;
109
110    kineticEnergyCorrection[1] = 0.9382723/3.727417;
111    slaterEffectiveCharge[0][1]=0.;
112    slaterEffectiveCharge[1][1]=0.;
113    slaterEffectiveCharge[2][1]=0.;
114    sCoefficient[0][1]=0.;
115    sCoefficient[1][1]=0.;
116    sCoefficient[2][1]=0.;
117  }
118  else
119  {
120      G4Exception("G4DNAMillerGreenExcitationModel::Initialise: alphaPlusPlus is not defined");
121  }
122
123  if (alphaPlusDef != 0)
124  {
125    alphaPlus = alphaPlusDef->GetParticleName();
126    lowEnergyLimit[alphaPlus] = 1. * keV;
127    highEnergyLimit[alphaPlus] = 10. * MeV;
128
129    kineticEnergyCorrection[2] = 0.9382723/3.727417;
130    slaterEffectiveCharge[0][2]=2.0;
131    slaterEffectiveCharge[1][2]=1.15;
132    slaterEffectiveCharge[2][2]=1.15;
133    sCoefficient[0][2]=0.7;
134    sCoefficient[1][2]=0.15;
135    sCoefficient[2][2]=0.15;
136  }
137  else
138  {
139    G4Exception("G4DNAMillerGreenExcitationModel::Initialise: alphaPlus is not defined");
140  }
141
142  if (heliumDef != 0)
143  {
144    helium = heliumDef->GetParticleName();
145    lowEnergyLimit[helium] = 1. * keV;
146    highEnergyLimit[helium] = 10. * MeV;
147   
148    kineticEnergyCorrection[3] = 0.9382723/3.727417;
149    slaterEffectiveCharge[0][3]=1.7;
150    slaterEffectiveCharge[1][3]=1.15;
151    slaterEffectiveCharge[2][3]=1.15;
152    sCoefficient[0][3]=0.5;
153    sCoefficient[1][3]=0.25;
154    sCoefficient[2][3]=0.25;
155  }
156  else
157  {
158    G4Exception("G4DNAMillerGreenExcitationModel::Initialise: helium is not defined");
159  }
160
161  if (particle==protonDef) 
162  {
163    SetLowEnergyLimit(lowEnergyLimit[proton]);
164    SetHighEnergyLimit(highEnergyLimit[proton]);
165  }
166
167  if (particle==alphaPlusPlusDef) 
168  {
169    SetLowEnergyLimit(lowEnergyLimit[alphaPlusPlus]);
170    SetHighEnergyLimit(highEnergyLimit[alphaPlusPlus]);
171  }
172
173  if (particle==alphaPlusDef) 
174  {
175    SetLowEnergyLimit(lowEnergyLimit[alphaPlus]);
176    SetHighEnergyLimit(highEnergyLimit[alphaPlus]);
177  }
178
179  if (particle==heliumDef) 
180  {
181    SetLowEnergyLimit(lowEnergyLimit[helium]);
182    SetHighEnergyLimit(highEnergyLimit[helium]);
183  }
184
185  //
186 
187  nLevels = waterExcitation.NumberOfLevels();
188
189  //
190  if( verboseLevel>0 ) 
191  { 
192    G4cout << "Miller & Green excitation model is initialized " << G4endl
193           << "Energy range: "
194           << LowEnergyLimit() / eV << " eV - "
195           << HighEnergyLimit() / keV << " keV for "
196           << particle->GetParticleName()
197           << G4endl;
198  }
199 
200  if(!isInitialised) 
201  {
202    isInitialised = true;
203 
204    if(pParticleChange)
205      fParticleChangeForGamma = reinterpret_cast<G4ParticleChangeForGamma*>(pParticleChange);
206    else
207      fParticleChangeForGamma = new G4ParticleChangeForGamma();
208  }   
209
210  // InitialiseElementSelectors(particle,cuts);
211 
212}
213
214//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
215
216G4double G4DNAMillerGreenExcitationModel::CrossSectionPerVolume(const G4Material* material,
217                                           const G4ParticleDefinition* particleDefinition,
218                                           G4double k,
219                                           G4double,
220                                           G4double)
221{
222  if (verboseLevel > 3)
223    G4cout << "Calling CrossSectionPerVolume() of G4DNAMillerGreenExcitationModel" << G4endl;
224
225 // Calculate total cross section for model
226
227  G4DNAGenericIonsManager *instance;
228  instance = G4DNAGenericIonsManager::Instance();
229
230  if (
231      particleDefinition != G4Proton::ProtonDefinition()
232      &&
233      particleDefinition != instance->GetIon("alpha++")
234      &&
235      particleDefinition != instance->GetIon("alpha+")
236      &&
237      particleDefinition != instance->GetIon("helium")
238     )
239           
240    return 0;
241
242  G4double lowLim = 0;
243  G4double highLim = 0;
244  G4double crossSection = 0.;
245
246  if (material->GetName() == "G4_WATER")
247  {
248    const G4String& particleName = particleDefinition->GetParticleName();
249
250    std::map< G4String,G4double,std::less<G4String> >::iterator pos1;
251    pos1 = lowEnergyLimit.find(particleName);
252
253    if (pos1 != lowEnergyLimit.end())
254    {
255      lowLim = pos1->second;
256    }
257
258    std::map< G4String,G4double,std::less<G4String> >::iterator pos2;
259    pos2 = highEnergyLimit.find(particleName);
260
261    if (pos2 != highEnergyLimit.end())
262    {
263      highLim = pos2->second;
264    }
265
266    if (k >= lowLim && k < highLim)
267    {
268      crossSection = Sum(k,particleDefinition);
269
270      G4DNAGenericIonsManager *instance;
271      instance = G4DNAGenericIonsManager::Instance();
272
273      // add ONE or TWO electron-water excitation for alpha+ and helium
274 
275      if ( particleDefinition == instance->GetIon("alpha+") 
276           ||
277           particleDefinition == instance->GetIon("helium")
278         ) 
279      {
280          G4DNAEmfietzoglouExcitationModel * excitationXS = new G4DNAEmfietzoglouExcitationModel();
281
282          G4double sigmaExcitation=0;
283          G4double tmp =0.;
284         
285          if (k*0.511/3728 > 8.23*eV && k*0.511/3728 < 10*MeV ) sigmaExcitation = 
286            excitationXS->CrossSectionPerVolume(material,G4Electron::ElectronDefinition(),k*0.511/3728,tmp,tmp)
287            /material->GetAtomicNumDensityVector()[1];
288       
289          if ( particleDefinition == instance->GetIon("alpha+") ) 
290            crossSection = crossSection +  sigmaExcitation ;
291         
292          if ( particleDefinition == instance->GetIon("helium") ) 
293            crossSection = crossSection + 2*sigmaExcitation ;
294         
295          delete excitationXS;
296      }     
297
298    }
299
300    if (verboseLevel > 3)
301    {
302      G4cout << "---> Kinetic energy(eV)=" << k/eV << G4endl;
303      G4cout << " - Cross section per water molecule (cm^2)=" << crossSection/cm/cm << G4endl;
304      G4cout << " - Cross section per water molecule (cm^-1)=" << crossSection*material->GetAtomicNumDensityVector()[1]/(1./cm) << G4endl;
305    } 
306
307  } 
308
309 return crossSection*material->GetAtomicNumDensityVector()[1];             
310
311}
312
313//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
314
315void G4DNAMillerGreenExcitationModel::SampleSecondaries(std::vector<G4DynamicParticle*>* /*fvect*/,
316                                              const G4MaterialCutsCouple* /*couple*/,
317                                              const G4DynamicParticle* aDynamicParticle,
318                                              G4double,
319                                              G4double)
320{
321
322  if (verboseLevel > 3)
323    G4cout << "Calling SampleSecondaries() of G4DNAMillerGreenExcitationModel" << G4endl;
324
325  G4double particleEnergy0 = aDynamicParticle->GetKineticEnergy();
326 
327  G4int level = RandomSelect(particleEnergy0,aDynamicParticle->GetDefinition());
328
329  G4double excitationEnergy = waterExcitation.ExcitationEnergy(level);
330  G4double newEnergy = particleEnergy0 - excitationEnergy;
331 
332  if (newEnergy>0)
333  {
334      fParticleChangeForGamma->ProposeMomentumDirection(aDynamicParticle->GetMomentumDirection());
335      fParticleChangeForGamma->SetProposedKineticEnergy(newEnergy);
336      fParticleChangeForGamma->ProposeLocalEnergyDeposit(excitationEnergy);
337  }
338
339}
340
341//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
342
343G4double G4DNAMillerGreenExcitationModel::PartialCrossSection(G4double k, G4int excitationLevel, 
344                                                                  const G4ParticleDefinition* particleDefinition)
345{
346  //                               ( ( z * aj ) ^ omegaj ) * ( t - ej ) ^ nu
347  // sigma(t) = zEff^2 * sigma0 * --------------------------------------------
348  //                               jj ^ ( omegaj + nu ) + t ^ ( omegaj + nu )
349  //
350  // where t is the kinetic energy corrected by Helium mass over proton mass for Helium ions
351  //
352  // zEff is:
353  //  1 for protons
354  //  2 for alpha++
355  //  and  2 - c1 S_1s - c2 S_2s - c3 S_2p for alpha+ and He
356  //
357  // Dingfelder et al., RPC 59, 255-275, 2000 from Miller and Green (1973)
358  // Formula (34) and Table 2
359 
360  const G4double sigma0(1.E+8 * barn);
361  const G4double nu(1.);
362  const G4double aj[]={876.*eV, 2084.* eV, 1373.*eV, 692.*eV, 900.*eV};
363  const G4double jj[]={19820.*eV, 23490.*eV, 27770.*eV, 30830.*eV, 33080.*eV};
364  const G4double omegaj[]={0.85, 0.88, 0.88, 0.78, 0.78};
365 
366  G4int particleTypeIndex = 0;
367  G4DNAGenericIonsManager* instance;
368  instance = G4DNAGenericIonsManager::Instance();
369
370  if (particleDefinition == G4Proton::ProtonDefinition()) particleTypeIndex=0;
371  if (particleDefinition == instance->GetIon("alpha++")) particleTypeIndex=1;
372  if (particleDefinition == instance->GetIon("alpha+")) particleTypeIndex=2;
373  if (particleDefinition == instance->GetIon("helium")) particleTypeIndex=3;
374 
375  G4double tCorrected;
376  tCorrected = k * kineticEnergyCorrection[particleTypeIndex];
377
378  // SI - added protection
379  if (tCorrected < waterExcitation.ExcitationEnergy(excitationLevel)) return 0;
380  //
381 
382  G4int z = 10;
383
384  G4double numerator;
385  numerator = std::pow(z * aj[excitationLevel], omegaj[excitationLevel]) * 
386    std::pow(tCorrected - waterExcitation.ExcitationEnergy(excitationLevel), nu);
387
388  G4double power;
389  power = omegaj[excitationLevel] + nu;
390
391  G4double denominator;
392  denominator = std::pow(jj[excitationLevel], power) + std::pow(tCorrected, power);
393
394  G4double zEff = particleDefinition->GetPDGCharge() / eplus + particleDefinition->GetLeptonNumber();
395
396  zEff -= ( sCoefficient[0][particleTypeIndex] * S_1s(k, waterExcitation.ExcitationEnergy(excitationLevel), slaterEffectiveCharge[0][particleTypeIndex], 1.) +
397            sCoefficient[1][particleTypeIndex] * S_2s(k, waterExcitation.ExcitationEnergy(excitationLevel), slaterEffectiveCharge[1][particleTypeIndex], 2.) +
398            sCoefficient[2][particleTypeIndex] * S_2p(k, waterExcitation.ExcitationEnergy(excitationLevel), slaterEffectiveCharge[2][particleTypeIndex], 2.) );
399
400  G4double cross = sigma0 * zEff * zEff * numerator / denominator;
401
402  return cross;
403}
404
405//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
406
407G4int G4DNAMillerGreenExcitationModel::RandomSelect(G4double k,const G4ParticleDefinition* particle)
408{
409  G4int i = nLevels;
410  G4double value = 0.;
411  std::deque<double> values;
412 
413  G4DNAGenericIonsManager *instance;
414  instance = G4DNAGenericIonsManager::Instance();
415
416  if ( particle == instance->GetIon("alpha++") ||
417       particle == G4Proton::ProtonDefinition()  )
418  { 
419     while (i > 0)
420     {
421        i--;
422        G4double partial = PartialCrossSection(k,i,particle);
423        values.push_front(partial);
424        value += partial;
425     }
426
427     value *= G4UniformRand();
428   
429     i = nLevels;
430
431     while (i > 0)
432     {
433        i--;
434        if (values[i] > value) return i;
435        value -= values[i];
436     }
437  } 
438
439  // add ONE or TWO electron-water excitation for alpha+ and helium
440   
441  if ( particle == instance->GetIon("alpha+") 
442       ||
443       particle == instance->GetIon("helium")
444     ) 
445  {
446    while (i>0)
447    {
448          i--;
449         
450          G4DNAEmfietzoglouExcitationModel * excitationXS = new G4DNAEmfietzoglouExcitationModel();
451         
452          G4double sigmaExcitation=0;
453
454          if (k*0.511/3728 > 8.23*eV && k*0.511/3728 < 10*MeV ) sigmaExcitation = excitationXS->PartialCrossSection(k*0.511/3728,i);
455 
456          G4double partial = PartialCrossSection(k,i,particle);
457          if (particle == instance->GetIon("alpha+")) partial = PartialCrossSection(k,i,particle) + sigmaExcitation;
458          if (particle == instance->GetIon("helium")) partial = PartialCrossSection(k,i,particle) + 2*sigmaExcitation;
459          values.push_front(partial);
460          value += partial;
461          delete excitationXS;
462    }
463 
464    value*=G4UniformRand();
465 
466    i=5;
467    while (i>0)
468    {
469          i--;
470   
471          if (values[i]>value) return i;
472 
473          value-=values[i];
474    }
475  }     
476
477  return 0;
478}
479
480//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
481
482G4double G4DNAMillerGreenExcitationModel::Sum(G4double k, const G4ParticleDefinition* particle)
483{
484  G4double totalCrossSection = 0.;
485
486  for (G4int i=0; i<nLevels; i++)
487  {
488    totalCrossSection += PartialCrossSection(k,i,particle);
489  }
490  return totalCrossSection;
491}
492
493//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
494
495G4double G4DNAMillerGreenExcitationModel::S_1s(G4double t,
496                                                          G4double energyTransferred,
497                                                          G4double slaterEffectiveCharge,
498                                                          G4double shellNumber)
499{
500  // 1 - e^(-2r) * ( 1 + 2 r + 2 r^2)
501  // Dingfelder, in Chattanooga 2005 proceedings, formula (7)
502 
503  G4double r = R(t, energyTransferred, slaterEffectiveCharge, shellNumber);
504  G4double value = 1. - std::exp(-2 * r) * ( ( 2. * r + 2. ) * r + 1. );
505 
506  return value;
507}
508
509
510//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
511
512G4double G4DNAMillerGreenExcitationModel::S_2s(G4double t,
513                                                          G4double energyTransferred,
514                                                          G4double slaterEffectiveCharge,
515                                                          G4double shellNumber)
516{
517  // 1 - e^(-2 r) * ( 1 + 2 r + 2 r^2 + 2 r^4)
518  // Dingfelder, in Chattanooga 2005 proceedings, formula (8)
519
520  G4double r = R(t, energyTransferred, slaterEffectiveCharge, shellNumber);
521  G4double value =  1. - std::exp(-2 * r) * (((2. * r * r + 2.) * r + 2.) * r + 1.);
522
523  return value;
524 
525}
526
527//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
528
529G4double G4DNAMillerGreenExcitationModel::S_2p(G4double t,
530                                                          G4double energyTransferred,
531                                                          G4double slaterEffectiveCharge,
532                                                          G4double shellNumber)
533{
534  // 1 - e^(-2 r) * ( 1 + 2 r + 2 r^2 + 4/3 r^3 + 2/3 r^4)
535  // Dingfelder, in Chattanooga 2005 proceedings, formula (9)
536
537  G4double r = R(t, energyTransferred, slaterEffectiveCharge, shellNumber);
538  G4double value =  1. - std::exp(-2 * r) * (((( 2./3. * r + 4./3.) * r + 2.) * r + 2.) * r  + 1.);
539
540  return value;
541}
542
543//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
544
545G4double G4DNAMillerGreenExcitationModel::R(G4double t,
546                                                       G4double energyTransferred,
547                                                       G4double slaterEffectiveCharge,
548                                                       G4double shellNumber) 
549{
550  // tElectron = m_electron / m_alpha * t
551  // Dingfelder, in Chattanooga 2005 proceedings, p 4
552
553  G4double tElectron = 0.511/3728. * t;
554  G4double value = 2. * tElectron * slaterEffectiveCharge / (energyTransferred * shellNumber);
555 
556  return value;
557}
558
Note: See TracBrowser for help on using the repository browser.