Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/examples/extended/exoticphysics/dmparticle/src/G4LDMBremModel.cc

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  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 //
 27 // -------------------------------------------------------------------
 28 //
 29 // 21.03.17 V. Grichine based on G4hBremsstrahlungModel
 30 //
 31 // Class Description:
 32 //
 33 // Implementation of energy loss for LDMPhoton emission by hadrons
 34 //
 35 
 36 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 37 
 38 #include "G4LDMBremModel.hh"
 39 
 40 #include "TestParameters.hh"
 41 
 42 #include "G4LDMPhoton.hh"
 43 #include "G4Log.hh"
 44 #include "G4ParticleChangeForLoss.hh"
 45 #include "G4PhysicalConstants.hh"
 46 #include "G4SystemOfUnits.hh"
 47 
 48 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 49 
 50 using namespace std;
 51 
 52 G4LDMBremModel::G4LDMBremModel(const G4ParticleDefinition* p, const G4String& nam)
 53   : G4MuBremsstrahlungModel(p, nam)
 54 {
 55   fEpsilon = TestParameters::GetPointer()->GetAlphaFactor();
 56   theLDMPhoton = G4LDMPhoton::LDMPhoton();
 57   fLDMPhotonMass = theLDMPhoton->GetPDGMass();
 58   minThreshold = 1.2 * fLDMPhotonMass;
 59 }
 60 
 61 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 62 
 63 G4LDMBremModel::~G4LDMBremModel() {}
 64 
 65 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 66 
 67 G4double G4LDMBremModel::ComputeDEDXPerVolume(const G4Material*, const G4ParticleDefinition*,
 68                                               G4double, G4double)
 69 {
 70   return 0.0;
 71 }
 72 
 73 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 74 
 75 G4double G4LDMBremModel::ComputeDMicroscopicCrossSection(G4double tkin, G4double Z,
 76                                                          G4double gammaEnergy)
 77 //  differential cross section
 78 {
 79   G4double dxsection = 0.;
 80 
 81   if (gammaEnergy > tkin || tkin < minThreshold) return dxsection;
 82   /*
 83   G4cout << "G4LDMBremModel m= " << mass
 84          << "  " << particle->GetParticleName()
 85          << "  Egamma(GeV)= " << gammaEnergy/GeV
 86          << "  Ekin(GeV)= " << tkin/GeV << G4endl;
 87   */
 88   G4double E = tkin + mass;
 89   G4double v = gammaEnergy / E;
 90   G4double delta = 0.5 * mass * mass * v / (E - gammaEnergy);
 91   G4double rab0 = delta * sqrte;
 92 
 93   G4int iz = std::max(1, std::min(G4lrint(Z), 99));
 94 
 95   G4double z13 = 1.0 / nist->GetZ13(iz);
 96   G4double dn = mass * nist->GetA27(iz) / (70. * MeV);
 97 
 98   G4double b = btf;
 99   if (1 == iz) b = bh;
100 
101   // nucleus contribution logarithm
102   G4double rab1 = b * z13;
103   G4double fn =
104     G4Log(rab1 / (dn * (electron_mass_c2 + rab0 * rab1)) * (mass + delta * (dn * sqrte - 2.)));
105   if (fn < 0.) fn = 0.;
106 
107   G4double x = 1.0 - v;
108 
109   if (particle->GetPDGSpin() != 0) {
110     x += 0.75 * v * v;
111   }
112 
113   dxsection = coeff * x * Z * Z * fn / gammaEnergy;
114   return dxsection;
115 }
116 
117 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
118 
119 G4double G4LDMBremModel::ComputeCrossSectionPerAtom(const G4ParticleDefinition*,
120                                                     G4double kineticEnergy, G4double Z, G4double,
121                                                     G4double cutEnergy, G4double maxEnergy)
122 {
123   G4double cross = 0.0;
124 
125   if (kineticEnergy <= lowestKinEnergy) return cross;
126 
127   G4double tmax = std::min(maxEnergy, kineticEnergy);
128   G4double cut = std::min(cutEnergy, kineticEnergy);
129 
130   cut = std::max(cut, minThreshold);
131   if (cut >= tmax) return cross;
132 
133   cross = ComputeMicroscopicCrossSection(kineticEnergy, Z, cut);
134 
135   if (tmax < kineticEnergy) {
136     cross -= ComputeMicroscopicCrossSection(kineticEnergy, Z, tmax);
137   }
138   cross *= fEpsilon * fEpsilon;
139 
140   return cross;
141 }
142 
143 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
144 
145 void G4LDMBremModel::SampleSecondaries(std::vector<G4DynamicParticle*>* vdp,
146                                        const G4MaterialCutsCouple* couple,
147                                        const G4DynamicParticle* dp, G4double minEnergy,
148                                        G4double maxEnergy)
149 {
150   G4double kineticEnergy = dp->GetKineticEnergy();
151   // check against insufficient energy
152   G4double tmax = std::min(kineticEnergy, maxEnergy);
153   G4double tmin = std::min(kineticEnergy, minEnergy);
154   tmin = std::max(tmin, minThreshold);
155   if (tmin >= tmax) return;
156 
157   // ===== sampling of energy transfer ======
158 
159   G4ParticleMomentum partDirection = dp->GetMomentumDirection();
160 
161   // select randomly one element constituing the material
162   const G4Element* anElement = SelectRandomAtom(couple, particle, kineticEnergy);
163   G4double Z = anElement->GetZ();
164 
165   G4double totalEnergy = kineticEnergy + mass;
166   G4double totalMomentum = sqrt(kineticEnergy * (kineticEnergy + 2.0 * mass));
167 
168   G4double func1 = tmin * ComputeDMicroscopicCrossSection(kineticEnergy, Z, tmin);
169 
170   G4double lnepksi, epksi;
171   G4double func2;
172 
173   G4double xmin = G4Log(tmin / MeV);
174   G4double xmax = G4Log(tmax / tmin);
175 
176   do {
177     lnepksi = xmin + G4UniformRand() * xmax;
178     epksi = MeV * G4Exp(lnepksi);
179     func2 = epksi * ComputeDMicroscopicCrossSection(kineticEnergy, Z, epksi);
180 
181     // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
182   } while (func2 < func1 * G4UniformRand());
183 
184   G4double gEnergy = std::max(epksi, fLDMPhotonMass);
185   G4double gMomentum = std::sqrt((epksi - fLDMPhotonMass) * (epksi + fLDMPhotonMass));
186 
187   // ===== sample angle =====
188 
189   G4double gam = totalEnergy / mass;
190   G4double rmax = gam * std::min(1.0, totalEnergy / gEnergy - 1.0);
191   G4double rmax2 = rmax * rmax;
192   G4double x = G4UniformRand() * rmax2 / (1.0 + rmax2);
193 
194   G4double theta = std::sqrt(x / (1.0 - x)) / gam;
195   G4double sint = std::sin(theta);
196   G4double phi = twopi * G4UniformRand();
197   G4double dirx = sint * cos(phi), diry = sint * sin(phi), dirz = cos(theta);
198 
199   G4ThreeVector gDirection(dirx, diry, dirz);
200   gDirection.rotateUz(partDirection);
201 
202   partDirection *= totalMomentum;
203   partDirection -= gMomentum * gDirection;
204   partDirection = partDirection.unit();
205 
206   // primary change
207 
208   kineticEnergy -= gEnergy;
209 
210   fParticleChange->SetProposedKineticEnergy(kineticEnergy);
211   fParticleChange->SetProposedMomentumDirection(partDirection);
212 
213   // save secondary
214   G4DynamicParticle* aLDMPhoton = new G4DynamicParticle(theLDMPhoton, gDirection, gEnergy);
215   vdp->push_back(aLDMPhoton);
216 }
217 
218 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
219