Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4LDMBremModel.cc
이 파일의 문서화 페이지로 가기
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: G4hLDMBremModel.cc 74020 2013-09-19 13:38:38Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // 21.03.17 V. Grichine based on G4hBremsstrahlungModel
31 //
32 // Class Description:
33 //
34 // Implementation of energy loss for LDMPhoton emission by hadrons
35 //
36 
37 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
38 
39 #include "G4LDMBremModel.hh"
40 #include "G4PhysicalConstants.hh"
41 #include "G4SystemOfUnits.hh"
42 #include "G4Log.hh"
43 #include "G4LDMPhoton.hh"
45 #include "TestParameters.hh"
46 
47 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
48 
49 using namespace std;
50 
52  const G4String& nam)
53  : G4MuBremsstrahlungModel(p, nam)
54 {
59 }
60 
61 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
62 
64 {}
65 
66 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
67 
69  const G4ParticleDefinition*,
71 {
72  return 0.0;
73 }
74 
75 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
76 
78  G4double tkin,
79  G4double Z,
80  G4double gammaEnergy)
81 // differential cross section
82 {
83  G4double dxsection = 0.;
84 
85  if( gammaEnergy > tkin || tkin < minThreshold) return dxsection;
86  /*
87  G4cout << "G4LDMBremModel m= " << mass
88  << " " << particle->GetParticleName()
89  << " Egamma(GeV)= " << gammaEnergy/GeV
90  << " Ekin(GeV)= " << tkin/GeV << G4endl;
91  */
92  G4double E = tkin + mass ;
93  G4double v = gammaEnergy/E ;
94  G4double delta = 0.5*mass*mass*v/(E-gammaEnergy) ;
95  G4double rab0=delta*sqrte ;
96 
97  G4int iz = std::max(1,std::min(G4lrint(Z),99));
98 
99  G4double z13 = 1.0/nist->GetZ13(iz);
100  G4double dn = mass*nist->GetA27(iz)/(70.*MeV);
101 
102  G4double b = btf;
103  if(1 == iz) b = bh;
104 
105  // nucleus contribution logarithm
106  G4double rab1=b*z13;
107  G4double fn=G4Log(rab1/(dn*(electron_mass_c2+rab0*rab1))*
108  (mass+delta*(dn*sqrte-2.))) ;
109  if(fn <0.) fn = 0. ;
110 
111  G4double x = 1.0 - v;
112 
113  if(particle->GetPDGSpin() != 0) { x += 0.75*v*v; }
114 
115  dxsection = coeff*x*Z*Z*fn/gammaEnergy;
116  return dxsection;
117 }
118 
119 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
120 
122  const G4ParticleDefinition*,
125  G4double cutEnergy,
126  G4double maxEnergy)
127 {
128  G4double cross = 0.0;
129 
130  if (kineticEnergy <= lowestKinEnergy) return cross;
131 
132  G4double tmax = std::min(maxEnergy, kineticEnergy);
133  G4double cut = std::min(cutEnergy, kineticEnergy);
134 
135  cut = std::max(cut, minThreshold);
136  if (cut >= tmax) return cross;
137 
138  cross = ComputeMicroscopicCrossSection (kineticEnergy, Z, cut);
139 
140  if(tmax < kineticEnergy)
141  {
142  cross -= ComputeMicroscopicCrossSection(kineticEnergy, Z, tmax);
143  }
144  cross *= fEpsilon*fEpsilon;
145 
146  return cross;
147 }
148 
149 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
150 
152  std::vector<G4DynamicParticle*>* vdp,
153  const G4MaterialCutsCouple* couple,
154  const G4DynamicParticle* dp,
155  G4double minEnergy,
156  G4double maxEnergy)
157 {
159  // check against insufficient energy
160  G4double tmax = std::min(kineticEnergy, maxEnergy);
161  G4double tmin = std::min(kineticEnergy, minEnergy);
162  tmin = std::max(tmin, minThreshold);
163  if(tmin >= tmax) return;
164 
165  // ===== sampling of energy transfer ======
166 
167  G4ParticleMomentum partDirection = dp->GetMomentumDirection();
168 
169  // select randomly one element constituing the material
170  const G4Element* anElement = SelectRandomAtom(couple,particle,kineticEnergy);
171  G4double Z = anElement->GetZ();
172 
173  G4double totalEnergy = kineticEnergy + mass;
174  G4double totalMomentum = sqrt(kineticEnergy*(kineticEnergy + 2.0*mass));
175 
176  G4double func1 = tmin*
177  ComputeDMicroscopicCrossSection(kineticEnergy,Z,tmin);
178 
179  G4double lnepksi, epksi;
180  G4double func2;
181 
182  G4double xmin = G4Log(tmin/MeV);
183  G4double xmax = G4Log(tmax/tmin);
184 
185  do
186  {
187  lnepksi = xmin + G4UniformRand()*xmax;
188  epksi = MeV*G4Exp(lnepksi);
189  func2 = epksi*ComputeDMicroscopicCrossSection(kineticEnergy,Z,epksi);
190 
191  // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
192  } while(func2 < func1*G4UniformRand());
193 
194  G4double gEnergy = std::max(epksi, fLDMPhotonMass);
195  G4double gMomentum =
196  std::sqrt((epksi - fLDMPhotonMass)*(epksi + fLDMPhotonMass));
197 
198  // ===== sample angle =====
199 
200  G4double gam = totalEnergy/mass;
201  G4double rmax = gam*std::min(1.0, totalEnergy/gEnergy - 1.0);
202  G4double rmax2= rmax*rmax;
203  G4double x = G4UniformRand()*rmax2/(1.0 + rmax2);
204 
205  G4double theta = std::sqrt(x/(1.0 - x))/gam;
206  G4double sint = std::sin(theta);
207  G4double phi = twopi * G4UniformRand() ;
208  G4double dirx = sint*cos(phi), diry = sint*sin(phi), dirz = cos(theta) ;
209 
210  G4ThreeVector gDirection(dirx, diry, dirz);
211  gDirection.rotateUz(partDirection);
212 
213  partDirection *= totalMomentum;
214  partDirection -= gMomentum*gDirection;
215  partDirection = partDirection.unit();
216 
217  // primary change
218 
219  kineticEnergy -= gEnergy;
220 
223 
224  // save secondary
225  G4DynamicParticle* aLDMPhoton =
226  new G4DynamicParticle(theLDMPhoton,gDirection,gEnergy);
227  vdp->push_back(aLDMPhoton);
228 }
229 
230 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
Float_t x
Definition: compare.C:6
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
T max(const T t1, const T t2)
brief Return the largest of the two arguments
static constexpr double MeV
Definition: G4SIunits.hh:214
G4LDMBremModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="ldmBrem")
G4double fLDMPhotonMass
const G4ParticleDefinition * theLDMPhoton
const G4ThreeVector & GetMomentumDirection() const
const char * p
Definition: xmltok.h:285
int func2(int i)
Definition: XFunc.cc:51
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:38
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:552
void SetProposedMomentumDirection(const G4ThreeVector &dir)
G4double GetPDGMass() const
G4double G4Log(G4double x)
Definition: G4Log.hh:230
static G4LDMPhoton * LDMPhoton()
Definition: G4LDMPhoton.cc:125
Float_t Z
virtual ~G4LDMBremModel()
G4double GetPDGSpin() const
double G4double
Definition: G4Types.hh:76
const G4ParticleDefinition * particle
static constexpr double electron_mass_c2
#define G4UniformRand()
Definition: Randomize.hh:53
static constexpr double twopi
Definition: G4SIunits.hh:76
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
int func1(int i)
Definition: XFunc.cc:40
Hep3Vector unit() const
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy) override
int G4lrint(double ad)
Definition: templates.hh:151
int G4int
Definition: G4Types.hh:78
G4double GetA27(G4int Z) const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
G4double GetKineticEnergy() const
G4double GetZ13(G4double Z) const
G4double GetZ() const
Definition: G4Element.hh:131
virtual G4double ComputeDMicroscopicCrossSection(G4double tkin, G4double Z, G4double gammaEnergy) override
G4ParticleChangeForLoss * fParticleChange
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
G4double ComputeMicroscopicCrossSection(G4double tkin, G4double Z, G4double cut)