source: trunk/source/processes/electromagnetic/utils/src/G4VAtomDeexcitation.cc @ 1350

Last change on this file since 1350 was 1350, checked in by garnier, 13 years ago

update to last version 4.9.4

  • Property svn:executable set to *
File size: 7.7 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: G4VAtomDeexcitation.cc,v 1.3 2010/11/04 12:55:09 vnivanch Exp $
27// GEANT4 tag $Name: emutils-V09-03-23 $
28//
29// -------------------------------------------------------------------
30//
31// GEANT4 Class class file
32//
33//
34// File name:     G4VAtomDeexcitation
35//
36// Author:        Alfonso Mantero & Vladimir Ivanchenko
37//
38// Creation date: 21.04.2010
39//
40// Modifications:
41//
42// Class Description:
43//
44// Abstract interface to energy loss models
45
46// -------------------------------------------------------------------
47//
48
49#include "G4VAtomDeexcitation.hh"
50#include "G4ParticleDefinition.hh"
51#include "G4DynamicParticle.hh"
52#include "G4Step.hh"
53#include "G4Region.hh"
54#include "G4RegionStore.hh"
55#include "G4MaterialCutsCouple.hh"
56#include "G4MaterialCutsCouple.hh"
57#include "G4Material.hh"
58#include "G4Element.hh"
59#include "G4ElementVector.hh"
60#include "Randomize.hh"
61#include "G4VParticleChange.hh"
62
63G4VAtomDeexcitation::G4VAtomDeexcitation(const G4String& modname, 
64                                         const G4String& pname) 
65  : verbose(1), name(modname), namePIXE(pname)
66{
67  activeZ.resize(93, false);
68  vdyn.reserve(5);
69  secVect.reserve(5);
70  theCoupleTable = 0;
71}
72
73G4VAtomDeexcitation::~G4VAtomDeexcitation()
74{}
75
76void G4VAtomDeexcitation::InitialiseAtomicDeexcitation()
77{
78  // Define list of couples
79  theCoupleTable = G4ProductionCutsTable::GetProductionCutsTable();
80  size_t numOfCouples = theCoupleTable->GetTableSize();
81  activeDeexcitationMedia.resize(numOfCouples, false);
82
83  // Define list of regions
84  G4RegionStore* regionStore = G4RegionStore::GetInstance();
85  size_t nRegions = regionStore->size();
86
87  // There is no active regions
88  if(0 == nRegions) { return; }
89
90  if(0 < verbose) {
91    G4cout << "### ================ Deexcitation model " << name
92           << " is activated for regions:" << G4endl; 
93  }
94
95  // Identify active media
96  for(size_t j=0; j<nRegions; ++j) {
97    const G4Region* reg = regionStore->GetRegion(activeRegions[j], false);
98    const G4ProductionCuts* rpcuts = reg->GetProductionCuts();
99    if(0 < verbose) {
100      G4cout << "###                   " << activeRegions[j] << G4endl; 
101    }
102 
103    for(size_t i=0; i<numOfCouples; ++i) {
104      if( !activeDeexcitationMedia[i] ) {
105
106        const G4MaterialCutsCouple* couple =
107          theCoupleTable->GetMaterialCutsCouple(i);
108        if (couple->GetProductionCuts() == rpcuts) {
109          activeDeexcitationMedia[i] = true;
110          const G4Material* mat = couple->GetMaterial();
111          const G4ElementVector* theElementVector = 
112            mat->GetElementVector();
113          G4int Z = (G4int)((*theElementVector)[i])->GetZ();
114          activeZ[Z] = true;
115        }
116      }
117    }
118  }
119
120  if(0 < verbose) {
121    G4cout << "### ================ PIXE model " << namePIXE
122           << G4endl; 
123  }
124
125  // Initialise derived class
126  InitialiseForNewRun();
127}
128
129void 
130G4VAtomDeexcitation::SetDeexcitationActiveRegion(const G4String& rname)
131{
132  G4String s = rname;
133  if(s == "" || s == "world" || s == "World" || s == "WORLD") {
134    s = "DefaultRegionForTheWorld";
135  }
136  size_t n = activeRegions.size();
137  if(n > 0) {
138    for(size_t i=0; i<n; ++i) { if(s == activeRegions[i]) { return; } }
139  }
140  activeRegions.push_back(s);
141}
142
143void 
144G4VAtomDeexcitation::AlongStepDeexcitation(G4VParticleChange* pParticleChange,
145                                           const G4Step& step, 
146                                           G4double& eLoss,
147                                           G4int coupleIndex)
148{
149  if(!CheckDeexcitationActiveRegion(coupleIndex) || eLoss == 0.0) { return; }
150
151  // step parameters
152  const G4StepPoint* preStep = step.GetPreStepPoint();
153  G4ThreeVector prePos = preStep->GetPosition();
154  G4ThreeVector delta = step.GetPostStepPoint()->GetPosition() - prePos;
155  G4double preTime = preStep->GetGlobalTime();
156  G4double dt = step.GetPostStepPoint()->GetGlobalTime() - preTime;
157  G4double truelength = step.GetStepLength();
158
159  // particle parameters
160  const G4Track* track = step.GetTrack();
161  const G4ParticleDefinition* part = track->GetDefinition();
162  G4double ekin = preStep->GetKineticEnergy() - 0.5*eLoss;
163
164  // media parameters
165  G4double gCut = (*theCoupleTable->GetEnergyCutsVector(coupleIndex))[0];
166  G4double eCut = (*theCoupleTable->GetEnergyCutsVector(coupleIndex))[1];
167  const G4Material* material = preStep->GetMaterial();
168  const G4ElementVector* theElementVector = material->GetElementVector();
169  const G4double* theAtomNumDensityVector = material->GetVecNbOfAtomsPerVolume();
170  G4int nelm = material->GetNumberOfElements();
171
172  // loop over deexcitations
173  for(G4int i=0; i<nelm; ++i) {
174    G4int Z = G4int((*theElementVector)[i]->GetZ());
175    if(Z > 5) {
176      G4double x = truelength*theAtomNumDensityVector[i];
177      if(x > 0.0) {
178        for(G4int ii=0; ii<9; ++ii) {
179          G4AtomicShellEnumerator as = G4AtomicShellEnumerator(ii);
180          const G4AtomicShell* shell = GetAtomicShell(Z, as);
181          if(gCut < shell->BindingEnergy()) {
182            G4double mfp = 
183              GetShellIonisationCrossSectionPerAtom(part, Z, as, ekin); 
184
185            // mfp is mean free path in units of step size
186            if(mfp > 0.0) {
187              mfp = 1.0/(x*mfp);
188              G4double stot = 0.0;
189
190              // sample ionisation points
191              do {
192                stot -= mfp*std::log(G4UniformRand());
193                if( stot > 1.0) { break; }
194
195                // sample deexcitation
196                vdyn.clear();
197                GenerateParticles(&vdyn, shell, Z, gCut, eCut); 
198                G4int nsec = vdyn.size();
199                if(nsec > 0) {
200                  secVect.clear();
201                  G4ThreeVector r = prePos  + stot*delta;
202                  G4double time   = preTime + stot*dt;
203                  for(G4int j=0; j<nsec; ++j) {
204                    G4DynamicParticle* dp = vdyn[j];
205                    G4double e = dp->GetKineticEnergy();
206
207                    // save new secondary if there is enough energy
208                    if(e <= eLoss) {
209                      G4Track* t = new G4Track(dp, time, r);
210                      secVect.push_back(t); 
211                      eLoss -= e;
212                    } else {
213                      delete dp;
214                    }                   
215                  }
216                }
217              } while ( stot < 1.0 && eLoss > 0.0);
218            }
219          }
220        }
221      } 
222    }
223  }
224  G4int nsec = secVect.size(); 
225  if(nsec > 0) {
226    G4int secondariesBefore = pParticleChange->GetNumberOfSecondaries();
227    pParticleChange->SetNumberOfSecondaries(nsec+secondariesBefore);
228    for(G4int j=0; j<nsec; ++j) {
229      pParticleChange->AddSecondary(secVect[j]);
230    }
231  }
232}
Note: See TracBrowser for help on using the repository browser.