source: trunk/source/processes/electromagnetic/adjoint/src/G4AdjointPhotoElectricModel.cc @ 1317

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

update geant4.9.3 tag

File size: 9.8 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: G4AdjointPhotoElectricModel.cc,v 1.5 2009/12/16 17:50:05 gunter Exp $
27// GEANT4 tag $Name: geant4-09-03 $
28//
29#include "G4AdjointPhotoElectricModel.hh"
30#include "G4AdjointCSManager.hh"
31
32
33#include "G4Integrator.hh"
34#include "G4TrackStatus.hh"
35#include "G4ParticleChange.hh"
36#include "G4AdjointElectron.hh"
37#include  "G4Gamma.hh"
38#include "G4AdjointGamma.hh"
39
40
41////////////////////////////////////////////////////////////////////////////////
42//
43G4AdjointPhotoElectricModel::G4AdjointPhotoElectricModel():
44 G4VEmAdjointModel("AdjointPEEffect")
45
46{ SetUseMatrix(false);
47  SetApplyCutInRange(false);
48  current_eEnergy =0.;
49  totAdjointCS=0.;
50  theAdjEquivOfDirectPrimPartDef =G4AdjointGamma::AdjointGamma();
51  theAdjEquivOfDirectSecondPartDef=G4AdjointElectron::AdjointElectron();
52  theDirectPrimaryPartDef=G4Gamma::Gamma();
53  second_part_of_same_type=false;
54  theDirectPEEffectModel = new G4PEEffectModel();
55 
56}
57////////////////////////////////////////////////////////////////////////////////
58//
59G4AdjointPhotoElectricModel::~G4AdjointPhotoElectricModel()
60{;}
61
62////////////////////////////////////////////////////////////////////////////////
63//
64void G4AdjointPhotoElectricModel::SampleSecondaries(const G4Track& aTrack,
65                                G4bool IsScatProjToProjCase,
66                                G4ParticleChange* fParticleChange)
67{ if (IsScatProjToProjCase) return ;
68
69  //Compute the totAdjointCS vectors if not already done for the current couple and electron energy
70  //-----------------------------------------------------------------------------------------------
71  const G4MaterialCutsCouple* aCouple = aTrack.GetMaterialCutsCouple();
72  const G4DynamicParticle* aDynPart =  aTrack.GetDynamicParticle() ;
73  G4double electronEnergy = aDynPart->GetKineticEnergy();
74  G4ThreeVector electronDirection= aDynPart->GetMomentumDirection() ;
75  pre_step_AdjointCS = totAdjointCS; //The last computed CS was  at pre step point
76  G4double adjCS;
77  adjCS = AdjointCrossSection(aCouple, electronEnergy,IsScatProjToProjCase);
78  post_step_AdjointCS = totAdjointCS; 
79                               
80
81
82
83  //Sample element
84  //-------------
85   const G4ElementVector* theElementVector = currentMaterial->GetElementVector();
86   size_t nelm =  currentMaterial->GetNumberOfElements();
87   G4double rand_CS= G4UniformRand()*xsec[nelm-1];
88   for (index_element=0; index_element<nelm-1; index_element++){
89        if (rand_CS<xsec[index_element]) break;
90   }
91       
92   //Sample shell and binding energy
93   //-------------
94   G4int nShells = (*theElementVector)[index_element]->GetNbOfAtomicShells();
95   rand_CS= shell_prob[index_element][nShells-1]*G4UniformRand();
96   G4int i  = 0; 
97   for (i=0; i<nShells-1; i++){
98        if (rand_CS<shell_prob[index_element][i]) break;
99   }
100   G4double gammaEnergy= electronEnergy+(*theElementVector)[index_element]->GetAtomicShell(i);
101       
102  //Sample cos theta
103  //Copy of the G4PEEffectModel cos theta sampling method ElecCosThetaDistribution.   
104  //This method cannot be used directly from G4PEEffectModel because it is a friend method. I should ask Vladimir to change that 
105  //------------------------------------------------------------------------------------------------   
106  //G4double cos_theta = theDirectPEEffectModel->ElecCosThetaDistribution(electronEnergy);
107       
108   G4double  cos_theta = 1.;
109   G4double gamma   = 1. + electronEnergy/electron_mass_c2;
110   if (gamma <= 5.) {
111        G4double beta  = std::sqrt(gamma*gamma-1.)/gamma;
112        G4double b     = 0.5*gamma*(gamma-1.)*(gamma-2);
113
114        G4double rndm,term,greject,grejsup;
115        if (gamma < 2.) grejsup = gamma*gamma*(1.+b-beta*b);
116        else            grejsup = gamma*gamma*(1.+b+beta*b);
117
118        do {    rndm = 1.-2*G4UniformRand();
119                cos_theta = (rndm+beta)/(rndm*beta+1.);
120                term = 1.-beta*cos_theta;
121                greject = (1.-cos_theta*cos_theta)*(1.+b*term)/(term*term);
122        } while(greject < G4UniformRand()*grejsup);
123  }
124       
125  // direction of the adjoint gamma electron
126  //---------------------------------------
127 
128     
129  G4double sin_theta = std::sqrt(1.-cos_theta*cos_theta);
130  G4double Phi     = twopi * G4UniformRand();
131  G4double dirx = sin_theta*std::cos(Phi),diry = sin_theta*std::sin(Phi),dirz = cos_theta;
132  G4ThreeVector adjoint_gammaDirection(dirx,diry,dirz);
133  adjoint_gammaDirection.rotateUz(electronDirection);
134 
135 
136 
137  //Weight correction
138 //-----------------------                                         
139  CorrectPostStepWeight(fParticleChange, aTrack.GetWeight(), electronEnergy,gammaEnergy,IsScatProjToProjCase); 
140 
141 
142 
143  //Create secondary and modify fParticleChange
144  //--------------------------------------------
145  G4DynamicParticle* anAdjointGamma = new G4DynamicParticle (
146                       G4AdjointGamma::AdjointGamma(),adjoint_gammaDirection, gammaEnergy);
147 
148 
149 
150 
151 
152  fParticleChange->ProposeTrackStatus(fStopAndKill);
153  fParticleChange->AddSecondary(anAdjointGamma);   
154       
155
156 
157
158}
159
160////////////////////////////////////////////////////////////////////////////////
161//
162void G4AdjointPhotoElectricModel::CorrectPostStepWeight(G4ParticleChange* fParticleChange, 
163                                                            G4double old_weight, 
164                                                            G4double adjointPrimKinEnergy, 
165                                                            G4double projectileKinEnergy ,
166                                                            G4bool  ) 
167{
168 G4double new_weight=old_weight;
169
170 G4double w_corr =G4AdjointCSManager::GetAdjointCSManager()->GetPostStepWeightCorrection()/factorCSBiasing;
171 w_corr*=post_step_AdjointCS/pre_step_AdjointCS; 
172 new_weight*=w_corr; 
173 new_weight*=projectileKinEnergy/adjointPrimKinEnergy;
174 fParticleChange->SetParentWeightByProcess(false);
175 fParticleChange->SetSecondaryWeightByProcess(false);
176 fParticleChange->ProposeParentWeight(new_weight);     
177}       
178
179////////////////////////////////////////////////////////////////////////////////
180//                     
181
182G4double G4AdjointPhotoElectricModel::AdjointCrossSection(const G4MaterialCutsCouple* aCouple,
183                                G4double electronEnergy,
184                                G4bool IsScatProjToProjCase)
185{ 
186 
187
188  if (IsScatProjToProjCase)  return 0.;
189
190       
191  if (aCouple !=currentCouple || current_eEnergy !=electronEnergy) {
192        totAdjointCS = 0.;
193        DefineCurrentMaterialAndElectronEnergy(aCouple, electronEnergy);
194        const G4ElementVector* theElementVector = currentMaterial->GetElementVector();
195        const double* theAtomNumDensityVector = currentMaterial->GetVecNbOfAtomsPerVolume();
196        size_t nelm =  currentMaterial->GetNumberOfElements();
197        for (index_element=0;index_element<nelm;index_element++){
198               
199                totAdjointCS +=AdjointCrossSectionPerAtom((*theElementVector)[index_element],electronEnergy)*theAtomNumDensityVector[index_element];
200                xsec[index_element] = totAdjointCS;
201        } 
202
203        totBiasedAdjointCS=std::min(totAdjointCS,0.01);
204//      totBiasedAdjointCS=totAdjointCS;
205        factorCSBiasing = totBiasedAdjointCS/totAdjointCS;
206        lastCS=totBiasedAdjointCS;
207       
208       
209  }
210  return totBiasedAdjointCS;
211
212 
213}
214////////////////////////////////////////////////////////////////////////////////
215//                     
216
217G4double G4AdjointPhotoElectricModel::AdjointCrossSectionPerAtom(const G4Element* anElement,G4double electronEnergy)
218{ 
219  G4int nShells = anElement->GetNbOfAtomicShells();
220  G4double Z= anElement->GetZ();
221  G4int i  = 0; 
222  G4double B0=anElement->GetAtomicShell(0);
223  G4double gammaEnergy = electronEnergy+B0;
224  G4double CS= theDirectPEEffectModel->ComputeCrossSectionPerAtom(G4Gamma::Gamma(),gammaEnergy,Z,0.,0.,0.);
225  G4double adjointCS =0.;
226  if (CS >0) adjointCS += CS/gammaEnergy; 
227  shell_prob[index_element][0] = adjointCS;                                         
228  for (i=1;i<nShells;i++){
229        //G4cout<<i<<G4endl;
230        G4double Bi_= anElement->GetAtomicShell(i-1);
231        G4double Bi = anElement->GetAtomicShell(i);
232        //G4cout<<Bi_<<'\t'<<Bi<<G4endl;
233        if (electronEnergy <Bi_-Bi) {
234                gammaEnergy = electronEnergy+Bi;
235               
236                CS=theDirectPEEffectModel->ComputeCrossSectionPerAtom(G4Gamma::Gamma(),gammaEnergy,Z,0.,0.,0.);
237                if (CS>0) adjointCS +=CS/gammaEnergy;
238        }
239        shell_prob[index_element][i] = adjointCS;       
240 
241  }
242  adjointCS*=electronEnergy;
243  return adjointCS;
244 
245}                               
246////////////////////////////////////////////////////////////////////////////////
247//                     
248
249void G4AdjointPhotoElectricModel::DefineCurrentMaterialAndElectronEnergy(const G4MaterialCutsCouple* couple, G4double anEnergy)
250{ currentCouple   = const_cast<G4MaterialCutsCouple*> (couple);
251  currentMaterial = const_cast<G4Material*> (couple->GetMaterial());
252  currentCoupleIndex = couple->GetIndex();
253  currentMaterialIndex = currentMaterial->GetIndex();
254  current_eEnergy = anEnergy;   
255}
Note: See TracBrowser for help on using the repository browser.