source: trunk/source/processes/electromagnetic/adjoint/src/G4AdjointeIonisationModel.cc @ 1197

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

nvx fichiers dans CVS

File size: 8.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: G4AdjointeIonisationModel.cc,v 1.2 2009/11/20 10:31:20 ldesorgh Exp $
27// GEANT4 tag $Name: geant4-09-03-cand-01 $
28//
29#include "G4AdjointeIonisationModel.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//
43G4AdjointeIonisationModel::G4AdjointeIonisationModel():
44 G4VEmAdjointModel("Inv_eIon_model")
45
46{
47 
48  UseMatrix =true;
49  UseMatrixPerElement = true;
50  ApplyCutInRange = true;
51  UseOnlyOneMatrixForAllElements = true;
52  CS_biasing_factor =1.;
53  WithRapidSampling = false; 
54 
55  theAdjEquivOfDirectPrimPartDef =G4AdjointElectron::AdjointElectron();
56  theAdjEquivOfDirectSecondPartDef=G4AdjointElectron::AdjointElectron();
57  theDirectPrimaryPartDef=G4Electron::Electron();
58  second_part_of_same_type=true;
59}
60////////////////////////////////////////////////////////////////////////////////
61//
62G4AdjointeIonisationModel::~G4AdjointeIonisationModel()
63{;}
64////////////////////////////////////////////////////////////////////////////////
65//
66void G4AdjointeIonisationModel::SampleSecondaries(const G4Track& aTrack,
67                                G4bool IsScatProjToProjCase,
68                                G4ParticleChange* fParticleChange)
69{ 
70
71
72 const G4DynamicParticle* theAdjointPrimary =aTrack.GetDynamicParticle();
73 
74 //Elastic inverse scattering
75 //---------------------------------------------------------
76 G4double adjointPrimKinEnergy = theAdjointPrimary->GetKineticEnergy();
77 G4double adjointPrimP =theAdjointPrimary->GetTotalMomentum();
78
79 if (adjointPrimKinEnergy>HighEnergyLimit*0.999){
80        return;
81 }
82 
83 //Sample secondary energy
84 //-----------------------
85 G4double projectileKinEnergy;
86 if (!WithRapidSampling ) { //used by default
87        projectileKinEnergy = SampleAdjSecEnergyFromCSMatrix(adjointPrimKinEnergy, IsScatProjToProjCase);
88       
89        CorrectPostStepWeight(fParticleChange, 
90                              aTrack.GetWeight(), 
91                              adjointPrimKinEnergy,
92                              projectileKinEnergy,
93                              IsScatProjToProjCase); //Caution !!!this weight correction should be always applied
94 }
95  else { //only  for test at the moment
96       
97        G4double Emin,Emax;
98        if (IsScatProjToProjCase) {
99                Emin=GetSecondAdjEnergyMinForScatProjToProjCase(adjointPrimKinEnergy,currentTcutForDirectSecond);
100                Emax=GetSecondAdjEnergyMaxForScatProjToProjCase(adjointPrimKinEnergy);
101        }
102        else {
103                Emin=GetSecondAdjEnergyMinForProdToProjCase(adjointPrimKinEnergy);
104                Emax=GetSecondAdjEnergyMaxForProdToProjCase(adjointPrimKinEnergy);
105        }
106        projectileKinEnergy = Emin*std::pow(Emax/Emin,G4UniformRand());
107       
108       
109       
110        lastCS=lastAdjointCSForScatProjToProjCase;
111        if ( !IsScatProjToProjCase) lastCS=lastAdjointCSForProdToProjCase;
112       
113        G4double new_weight=aTrack.GetWeight();
114        G4double used_diffCS=lastCS*std::log(Emax/Emin)/projectileKinEnergy;
115        G4double needed_diffCS=adjointPrimKinEnergy/projectileKinEnergy;
116        if (!IsScatProjToProjCase) needed_diffCS *=DiffCrossSectionPerVolumePrimToSecond(currentMaterial,projectileKinEnergy,adjointPrimKinEnergy);
117        else needed_diffCS *=DiffCrossSectionPerVolumePrimToScatPrim(currentMaterial,projectileKinEnergy,adjointPrimKinEnergy);
118        new_weight*=needed_diffCS/used_diffCS;
119        fParticleChange->SetParentWeightByProcess(false);
120        fParticleChange->SetSecondaryWeightByProcess(false);
121        fParticleChange->ProposeParentWeight(new_weight);
122       
123 
124 }     
125                                         
126 
127                 
128 //Kinematic:
129 //we consider a two body elastic scattering for the forward processes where the projectile knock on an e- at rest  and gives
130 // him part of its  energy
131 //----------------------------------------------------------------------------------------
132 
133 G4double projectileM0 = theAdjEquivOfDirectPrimPartDef->GetPDGMass();
134 G4double projectileTotalEnergy = projectileM0+projectileKinEnergy;
135 G4double projectileP2 = projectileTotalEnergy*projectileTotalEnergy - projectileM0*projectileM0;       
136 
137 
138 
139 //Companion
140 //-----------
141 G4double companionM0 = theAdjEquivOfDirectPrimPartDef->GetPDGMass();
142 if (IsScatProjToProjCase) {
143        companionM0=theAdjEquivOfDirectSecondPartDef->GetPDGMass();
144 } 
145 G4double companionTotalEnergy =companionM0+ projectileKinEnergy-adjointPrimKinEnergy;
146 G4double companionP2 = companionTotalEnergy*companionTotalEnergy - companionM0*companionM0;   
147 
148 
149 //Projectile momentum
150 //--------------------
151 G4double  P_parallel = (adjointPrimP*adjointPrimP +  projectileP2 - companionP2)/(2.*adjointPrimP);
152 G4double  P_perp = std::sqrt( projectileP2 -  P_parallel*P_parallel);
153 G4ThreeVector dir_parallel=theAdjointPrimary->GetMomentumDirection();
154 G4double phi =G4UniformRand()*2.*3.1415926;
155 G4ThreeVector projectileMomentum = G4ThreeVector(P_perp*std::cos(phi),P_perp*std::sin(phi),P_parallel);
156 projectileMomentum.rotateUz(dir_parallel);
157 
158 
159 
160 if (!IsScatProjToProjCase ){ //kill the primary and add a secondary
161        fParticleChange->ProposeTrackStatus(fStopAndKill);
162        fParticleChange->AddSecondary(new G4DynamicParticle(theAdjEquivOfDirectPrimPartDef,projectileMomentum));
163        //G4cout<<"projectileMomentum "<<projectileMomentum<<G4endl;
164 }
165 else {
166        fParticleChange->ProposeEnergy(projectileKinEnergy);
167        fParticleChange->ProposeMomentumDirection(projectileMomentum.unit());
168 }     
169       
170
171 
172
173}
174////////////////////////////////////////////////////////////////////////////////
175//
176//The implementation here is correct for energy loss process, for the photoelectric and compton scattering the method should be redefine 
177G4double G4AdjointeIonisationModel::DiffCrossSectionPerAtomPrimToSecond(
178                                      G4double kinEnergyProj, 
179                                      G4double kinEnergyProd,
180                                      G4double Z, 
181                                      G4double )
182{
183 G4double dSigmadEprod=0;
184 G4double Emax_proj = GetSecondAdjEnergyMaxForProdToProjCase(kinEnergyProd);
185 G4double Emin_proj = GetSecondAdjEnergyMinForProdToProjCase(kinEnergyProd);
186 
187 
188 if (kinEnergyProj>Emin_proj && kinEnergyProj<=Emax_proj){ //the produced particle should have a kinetic energy smaller than the projectile
189        dSigmadEprod=Z*DiffCrossSectionMoller(kinEnergyProj,kinEnergyProd);
190 }
191 return dSigmadEprod;   
192 
193 
194 
195}
196
197//////////////////////////////////////////////////////////////////////////////
198//
199G4double G4AdjointeIonisationModel::DiffCrossSectionMoller(G4double kinEnergyProj,G4double kinEnergyProd){
200  G4double electron_mass_c2=0.51099906*MeV;
201  G4double energy = kinEnergyProj + electron_mass_c2;
202  G4double x   = kinEnergyProd/kinEnergyProj;
203  G4double gam    = energy/electron_mass_c2;
204  G4double gamma2 = gam*gam;
205  G4double beta2  = 1.0 - 1.0/gamma2;
206 
207  G4double g = (2.0*gam - 1.0)/gamma2;
208  G4double y = 1.0 - x;
209  G4double fac=twopi_mc2_rcl2/electron_mass_c2;
210  G4double dCS = fac*( 1.-g + ((1.0 - g*x)/(x*x)) + ((1.0 - g*y)/(y*y)))/(beta2*(gam-1));
211  return dCS/kinEnergyProj;
212 
213 
214
215} 
216
Note: See TracBrowser for help on using the repository browser.