CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
G4SeltzerBergerModel95.cc
Go to the documentation of this file.
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 // GEANT4 Class file
29 //
30 //
31 // File name: G4SeltzerBergerModel95
32 //
33 // Author: Andreas Schaelicke
34 //
35 // Creation date: 12.08.2008
36 //
37 // Modifications:
38 //
39 // 13.11.08 add SetLPMflag and SetLPMconstant methods
40 // 13.11.08 change default LPMconstant value
41 // 13.10.10 add angular distributon interface (VI)
42 //
43 // Main References:
44 // Y.-S.Tsai, Rev. Mod. Phys. 46 (1974) 815; Rev. Mod. Phys. 49 (1977) 421.
45 // S.Klein, Rev. Mod. Phys. 71 (1999) 1501.
46 // T.Stanev et.al., Phys. Rev. D25 (1982) 1291.
47 // M.L.Ter-Mikaelian, High-energy Electromagnetic Processes in Condensed Media, Wiley, 1972.
48 //
49 // -------------------------------------------------------------------
50 //
51 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
52 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
53 
54 #include "G4SeltzerBergerModel95.hh"
55 #include "G4Electron.hh"
56 #include "G4Positron.hh"
57 #include "G4Gamma.hh"
58 #include "Randomize.hh"
59 #include "G4Material.hh"
60 #include "G4Element.hh"
61 #include "G4ElementVector.hh"
62 #include "G4ProductionCutsTable.hh"
63 #include "G4ParticleChangeForLoss.hh"
64 #include "G4LossTableManager.hh"
65 #include "G4ModifiedTsai.hh"
66 
67 #include "G4Physics2DVector95.hh"
68 
69 #include "G4ios.hh"
70 #include <fstream>
71 #include <iomanip>
72 
73 
74 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
75 
76 using namespace std;
77 
78 G4SeltzerBergerModel95::G4SeltzerBergerModel95(const G4ParticleDefinition* p,
79  const G4String& name)
80  : G4eBremsstrahlungRelModel95(p,name)
81 {
82  SetLowEnergyLimit(0.0);
83  SetLPMFlag(false);
84  dataSB.resize(101,0);
85 }
86 
87 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
88 
89 G4SeltzerBergerModel95::~G4SeltzerBergerModel95()
90 {
91  for(size_t i=0; i<101; ++i) { delete dataSB[i]; }
92 }
93 
94 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
95 
96 void G4SeltzerBergerModel95::Initialise(const G4ParticleDefinition* p,
97  const G4DataVector& cuts)
98 {
99  // check environment variable
100  // Build the complete string identifying the file with the data set
101  char* path = getenv("G4LEDATA");
102 
103  // Access to elements
104  const G4ElementTable* theElmTable = G4Element::GetElementTable();
105  size_t numOfElm = G4Element::GetNumberOfElements();
106  if(numOfElm > 0) {
107  for(size_t i=0; i<numOfElm; ++i) {
108  G4int Z = G4int(((*theElmTable)[i])->GetZ());
109  if(Z < 1) { Z = 1; }
110  else if(Z > 100) { Z = 100; }
111  //G4cout << "Z= " << Z << G4endl;
112  // Initialisation
113  if(!dataSB[Z]) { ReadData(Z, path); }
114  }
115  }
116 
117  G4eBremsstrahlungRelModel95::Initialise(p, cuts);
118 }
119 
120 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
121 
122 void G4SeltzerBergerModel95::ReadData(size_t Z, const char* path)
123 {
124  // G4cout << "ReadData Z= " << Z << G4endl;
125  // G4cout << "Status for Z= " << dataSB[Z] << G4endl;
126  //if(path) { G4cout << path << G4endl; }
127  if(dataSB[Z]) { return; }
128  const char* datadir = path;
129 
130  if(!datadir) {
131  datadir = getenv("G4LEDATA");
132  if(!datadir) {
133  G4Exception("G4SeltzerBergerModel95::ReadData() Environment variable G4LEDATA not defined");
134  //G4Exception("G4SeltzerBergerModel95::ReadData()","em0006",FatalException,
135  // "Environment variable G4LEDATA not defined");
136  return;
137  }
138  }
139  std::ostringstream ost;
140  ost << datadir << "/brem_SB/br" << Z;
141  std::ifstream fin(ost.str().c_str());
142  if( !fin.is_open()) {
143  //G4ExceptionDescription ed;
144  G4cout << "Bremsstrahlung data file <" << ost.str().c_str()
145  << "> is not opened!" << G4endl;
146  G4Exception("G4SeltzerBergerModel95::ReadData() G4LEDATA version should be G4EMLOW6.23 or later.");
147  return;
148  }
149  //G4cout << "G4SeltzerBergerModel95 read from <" << ost.str().c_str()
150  // << ">" << G4endl;
151  G4Physics2DVector95* v = new G4Physics2DVector95();
152  if(v->Retrieve(fin)) { dataSB[Z] = v; }
153  else {
154  //G4ExceptionDescription ed;
155  G4cout << "Bremsstrahlung data file <" << ost.str().c_str()
156  << "> is not retrieved!" << G4endl;
157  G4Exception("G4SeltzerBergerModel95::ReadData() G4LEDATA version should be G4EMLOW6.23 or later.");
158  delete v;
159  }
160  // G4cout << dataSB[Z] << G4endl;
161 }
162 
163 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
164 
165 G4double G4SeltzerBergerModel95::ComputeDXSectionPerAtom(G4double gammaEnergy)
166 {
167 
168  if(gammaEnergy < 0.0 || kinEnergy <= 0.0) { return 0.0; }
169  G4double x = gammaEnergy/kinEnergy;
170  G4double y = log(kinEnergy/MeV);
171  G4int Z = G4int(currentZ);
172  //G4cout << "G4SeltzerBergerModel95::ComputeDXSectionPerAtom Z= " << Z
173  // << " x= " << x << " y= " << y << " " << dataSB[Z] << G4endl;
174  if(!dataSB[Z]) { ReadData(Z); }
175  G4double invb2 = totalEnergy*totalEnergy/(kinEnergy*(kinEnergy + 2*particleMass));
176  G4double cross = dataSB[Z]->Value(x,y)*invb2*millibarn/bremFactor;
177 
178  if(!isElectron) {
179  if(1 - x < 1.e-20) { cross = 0.0; }
180  else {
181  G4double invbeta1 = sqrt(invb2);
182  G4double e2 = kinEnergy - gammaEnergy;
183  G4double invbeta2 = (e2 + particleMass)/sqrt(e2*(e2 + 2*particleMass));
184  cross *= exp(twopi*fine_structure_const*currentZ*(invbeta1 - invbeta2));
185  }
186  }
187 
188  return cross;
189 }
190 
191 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
192 
193 void
194 G4SeltzerBergerModel95::SampleSecondaries(std::vector<G4DynamicParticle*>* vdp,
195  const G4MaterialCutsCouple* couple,
196  const G4DynamicParticle* dp,
197  G4double cutEnergy,
198  G4double maxEnergy)
199 {
200  G4double kineticEnergy = dp->GetKineticEnergy();
201  G4double cut = std::min(cutEnergy, kineticEnergy);
202  G4double emax = std::min(maxEnergy, kineticEnergy);
203  if(cut >= emax) { return; }
204 
205  SetupForMaterial(particle, couple->GetMaterial(), kineticEnergy);
206 
207  const G4Element* elm =
208  SelectRandomAtom(couple,particle,kineticEnergy,cut,emax);
209  SetCurrentElement(elm->GetZ());
210  G4int Z = G4int(currentZ);
211 
212  totalEnergy = kineticEnergy + particleMass;
213  densityCorr = densityFactor*totalEnergy*totalEnergy;
214  G4double totMomentum = sqrt(kineticEnergy*(totalEnergy + electron_mass_c2));
215  G4ThreeVector direction = dp->GetMomentumDirection();
216  /*
217  G4cout << "G4SeltzerBergerModel95::SampleSecondaries E(MeV)= "
218  << kineticEnergy/MeV
219  << " Z= " << Z << " cut(MeV)= " << cut/MeV
220  << " emax(MeV)= " << emax/MeV << " corr= " << densityCorr << G4endl;
221  */
222  G4double xmin = log(cut*cut + densityCorr);
223  G4double xmax = log(emax*emax + densityCorr);
224  G4double y = log(kineticEnergy/MeV);
225 
226  G4double gammaEnergy, v;
227 
228  // majoranta
229  G4double vmax = dataSB[Z]->Value(cut/kineticEnergy, y);
230  if(isElectron && Z > 12 && kineticEnergy < 100*keV) {
231  if((Z < 41 && kineticEnergy < 10*keV) ||
232  (Z >= 41 && Z < 61 && kineticEnergy < 50*keV) ||
233  (Z >= 61) )
234  {
235  v = 1.05*dataSB[Z]->Value(emax/kineticEnergy, y);
236  if(v > vmax) { vmax = v; }
237  }
238  }
239 
240  //G4cout<<"y= "<<y<<" xmin= "<<xmin<<" xmax= "<<xmax<<" vmax= "<<vmax<<G4endl;
241 
242  do {
243  G4double x = exp(xmin + G4UniformRand()*(xmax - xmin)) - densityCorr;
244  if(x < 0.0) { x = 0.0; }
245  gammaEnergy = sqrt(x);
246  G4double x1 = gammaEnergy/kineticEnergy;
247  v = dataSB[Z]->Value(x1, y);
248 
249  if(!isElectron) {
250  if(1 - x1 < 1.e-20) { v = 0.0; }
251  else {
252  G4double e1 = kineticEnergy - cut;
253  G4double invbeta1 = (e1 + particleMass)/sqrt(e1*(e1 + 2*particleMass));
254  G4double e2 = kineticEnergy - gammaEnergy;
255  G4double invbeta2 = (e2 + particleMass)/sqrt(e2*(e2 + 2*particleMass));
256  v *= exp(twopi*fine_structure_const*currentZ*(invbeta1 - invbeta2));
257  }
258  }
259 
260  if ( v > 1.5*vmax ) {
261  G4cout << "### G4SeltzerBergerModel95 Warning: Majoranta exceeded! "
262  << v << " > " << vmax
263  << " Egamma(MeV)= " << gammaEnergy
264  << " Ee(MeV)= " << kineticEnergy
265  << " Z= " << Z << " " << particle->GetParticleName()
266  << G4endl;
267  }
268 
269  } while (v < vmax*G4UniformRand());
270 
271  //
272  // angles of the emitted gamma. ( Z - axis along the parent particle)
273  // use general interface
274  //
275  G4double theta =
276  GetAngularDistribution()->PolarAngle(totalEnergy,totalEnergy-gammaEnergy,Z);
277 
278  G4double sint = sin(theta);
279  G4double phi = twopi * G4UniformRand();
280  G4ThreeVector gammaDirection(sint*cos(phi),sint*sin(phi), cos(theta));
281  gammaDirection.rotateUz(direction);
282 
283  // create G4DynamicParticle object for the Gamma
284  G4DynamicParticle* g =
285  new G4DynamicParticle(theGamma,gammaDirection,gammaEnergy);
286  vdp->push_back(g);
287 
288  G4ThreeVector dir = totMomentum*direction - gammaEnergy*gammaDirection;
289  direction = dir.unit();
290 
291  // energy of primary
292  G4double finalE = kineticEnergy - gammaEnergy;
293 
294  // stop tracking and create new secondary instead of primary
295  if(gammaEnergy > SecondaryThreshold()) {
296  fParticleChange->ProposeTrackStatus(fStopAndKill);
297  fParticleChange->SetProposedKineticEnergy(0.0);
298  G4DynamicParticle* el =
299  new G4DynamicParticle(const_cast<G4ParticleDefinition*>(particle),
300  direction, finalE);
301  vdp->push_back(el);
302 
303  // continue tracking
304  } else {
305  fParticleChange->SetProposedMomentumDirection(direction);
306  fParticleChange->SetProposedKineticEnergy(finalE);
307  }
308 }
309 
310 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
311 
312 
const double Z[kNumberCalorimeter]
int i
Definition: DBlmapReader.cc:9
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Geom::Theta< T > theta() const
#define min(a, b)
Definition: mlp_lapack.h:161
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e g
Definition: Activities.doc:4
const TString datadir("/src/Fireworks/Core/")
list path
Definition: scaleCards.py:51
bool isElectron(const Candidate &part)
Definition: pdgIdUtils.h:7
T sqrt(T t)
Definition: SSEVec.h:46
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
dbl *** dir
Definition: mlp_gen.cc:35
Definition: DDAxes.h:10
mathSSE::Vec4< T > v
Basic3DVector cross(const Basic3DVector &v) const
Vector product, or &quot;cross&quot; product, with a vector of same type.
Definition: DDAxes.h:10