Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4LivermoreBremsstrahlungModel.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 // $Id$
27 //
28 // Author: Luciano Pandola
29 // on base of G4LowEnergyBremsstrahlung developed by A.Forti and V.Ivanchenko
30 //
31 // History:
32 // --------
33 // 03 Mar 2009 L Pandola Migration from process to model
34 // 12 Apr 2009 V Ivanchenko Cleanup initialisation and generation of secondaries:
35 // - apply internal high-energy limit only in constructor
36 // - do not apply low-energy limit (default is 0)
37 // - added MinEnergyCut method
38 // - do not change track status
39 // - do not initialize element selectors
40 // - use cut value from the interface
41 // - fixed bug in sampling of angles between keV and MeV
42 // 19 May 2009 L Pandola Explicitely set to zero pointers deleted in
43 // Initialise(), since they might be checked later on
44 //
45 
47 #include "G4PhysicalConstants.hh"
48 #include "G4SystemOfUnits.hh"
49 #include "G4ParticleDefinition.hh"
50 #include "G4MaterialCutsCouple.hh"
51 
52 #include "G4DynamicParticle.hh"
53 #include "G4Element.hh"
54 #include "G4Gamma.hh"
55 #include "G4Electron.hh"
57 //
59 #include "G4ModifiedTsai.hh"
60 #include "G4Generator2BS.hh"
61 //#include "G4Generator2BN.hh"
62 //
64 //
65 #include "G4VEnergySpectrum.hh"
67 #include "G4VEMDataSet.hh"
68 
69 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
70 
71 
73  const G4String& nam)
74  :G4VEmModel(nam),fParticleChange(0),isInitialised(false),
75  crossSectionHandler(0),energySpectrum(0)
76 {
77  fIntrinsicLowEnergyLimit = 10.0*eV;
78  fIntrinsicHighEnergyLimit = 100.0*GeV;
79  fNBinEnergyLoss = 360;
80  // SetLowEnergyLimit(fIntrinsicLowEnergyLimit);
81  SetHighEnergyLimit(fIntrinsicHighEnergyLimit);
82  //
83  verboseLevel = 0;
85  //
86  //generatorName = "tsai";
87  //angularDistribution = new G4ModifiedTsai("TsaiGenerator"); //default generator
88  //
89  //TsaiAngularDistribution = new G4ModifiedTsai("TsaiGenerator");
90  //
91 }
92 
93 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
94 
96 {
97  if (crossSectionHandler) delete crossSectionHandler;
98  if (energySpectrum) delete energySpectrum;
99  energyBins.clear();
100  //delete angularDistribution;
101  //delete TsaiAngularDistribution;
102 }
103 
104 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
105 
107  const G4DataVector& cuts)
108 {
109  //Check that the Livermore Bremsstrahlung is NOT attached to e+
110  if (particle != G4Electron::Electron())
111  {
112  G4Exception("G4LivermoreBremsstrahlungModel::Initialise",
113  "em0002",FatalException,"Livermore Bremsstrahlung Model is applicable only to electrons");
114  }
115  //Prepare energy spectrum
116  if (energySpectrum)
117  {
118  delete energySpectrum;
119  energySpectrum = 0;
120  }
121 
122  energyBins.clear();
123  for(size_t i=0; i<15; i++)
124  {
125  G4double x = 0.1*((G4double)i);
126  if(i == 0) x = 0.01;
127  if(i == 10) x = 0.95;
128  if(i == 11) x = 0.97;
129  if(i == 12) x = 0.99;
130  if(i == 13) x = 0.995;
131  if(i == 14) x = 1.0;
132  energyBins.push_back(x);
133  }
134  const G4String dataName("/brem/br-sp.dat");
135  energySpectrum = new G4eBremsstrahlungSpectrum(energyBins,dataName);
136 
137  if (verboseLevel > 0)
138  G4cout << "G4eBremsstrahlungSpectrum is initialized" << G4endl;
139 
140  //Initialize cross section handler
141  if (crossSectionHandler)
142  {
143  delete crossSectionHandler;
144  crossSectionHandler = 0;
145  }
146  G4VDataSetAlgorithm* interpolation = 0;//new G4SemiLogInterpolation();
147  crossSectionHandler = new G4BremsstrahlungCrossSectionHandler(energySpectrum,interpolation);
148  crossSectionHandler->Initialise(0,LowEnergyLimit(),HighEnergyLimit(),
149  fNBinEnergyLoss);
150  crossSectionHandler->Clear();
151  crossSectionHandler->LoadShellData("brem/br-cs-");
152  //This is used to retrieve cross section values later on
153  G4VEMDataSet* p = crossSectionHandler->BuildMeanFreePathForMaterials(&cuts);
154  delete p;
155 
156  if (verboseLevel > 0)
157  {
158  G4cout << "Livermore Bremsstrahlung model is initialized " << G4endl
159  << "Energy range: "
160  << LowEnergyLimit() / keV << " keV - "
161  << HighEnergyLimit() / GeV << " GeV"
162  << G4endl;
163  }
164 
165  if (verboseLevel > 1)
166  {
167  G4cout << "Cross section data: " << G4endl;
168  crossSectionHandler->PrintData();
169  G4cout << "Parameters: " << G4endl;
170  energySpectrum->PrintData();
171  }
172 
173  if(isInitialised) return;
175  isInitialised = true;
176 }
177 
178 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
179 
181  const G4MaterialCutsCouple*)
182 {
183  return 250.*eV;
184 }
185 
186 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
187 
188 G4double
192  G4double cutEnergy,
193  G4double)
194 {
195  G4int iZ = (G4int) Z;
196  if (!crossSectionHandler)
197  {
198  G4Exception("G4LivermoreBremsstrahlungModel::ComputeCrossSectionPerAtom",
199  "em1007",FatalException,"The cross section handler is not correctly initialized");
200  return 0;
201  }
202 
203  //The cut is already included in the crossSectionHandler
204  G4double cs =
205  crossSectionHandler->GetCrossSectionAboveThresholdForElement(energy,cutEnergy,iZ);
206 
207  if (verboseLevel > 1)
208  {
209  G4cout << "G4LivermoreBremsstrahlungModel " << G4endl;
210  G4cout << "Cross section for gamma emission > " << cutEnergy/keV << " keV at " <<
211  energy/keV << " keV and Z = " << iZ << " --> " << cs/barn << " barn" << G4endl;
212  }
213  return cs;
214 }
215 
216 
217 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
218 
220  const G4ParticleDefinition* ,
221  G4double kineticEnergy,
222  G4double cutEnergy)
223 {
224  G4double sPower = 0.0;
225 
226  const G4ElementVector* theElementVector = material->GetElementVector();
227  size_t NumberOfElements = material->GetNumberOfElements() ;
228  const G4double* theAtomicNumDensityVector =
229  material->GetAtomicNumDensityVector();
230 
231  // loop for elements in the material
232  for (size_t iel=0; iel<NumberOfElements; iel++ )
233  {
234  G4int iZ = (G4int)((*theElementVector)[iel]->GetZ());
235  G4double e = energySpectrum->AverageEnergy(iZ, 0.0,cutEnergy,
236  kineticEnergy);
237  G4double cs= crossSectionHandler->FindValue(iZ,kineticEnergy);
238  sPower += e * cs * theAtomicNumDensityVector[iel];
239  }
240 
241  if (verboseLevel > 2)
242  {
243  G4cout << "G4LivermoreBremsstrahlungModel " << G4endl;
244  G4cout << "Stopping power < " << cutEnergy/keV << " keV at " <<
245  kineticEnergy/keV << " keV = " << sPower/(keV/mm) << " keV/mm" << G4endl;
246  }
247 
248  return sPower;
249 }
250 
251 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
252 
253 void G4LivermoreBremsstrahlungModel::SampleSecondaries(std::vector<G4DynamicParticle*>* fvect,
254  const G4MaterialCutsCouple* couple,
255  const G4DynamicParticle* aDynamicParticle,
256  G4double energyCut,
257  G4double)
258 {
259 
260  G4double kineticEnergy = aDynamicParticle->GetKineticEnergy();
261 
262  // this is neede for pathalogical cases of no ionisation
263  if (kineticEnergy <= fIntrinsicLowEnergyLimit)
264  {
267  return;
268  }
269 
270  //Sample material
271  G4int Z = crossSectionHandler->SelectRandomAtom(couple, kineticEnergy);
272 
273  //Sample gamma energy
274  G4double tGamma = energySpectrum->SampleEnergy(Z, energyCut, kineticEnergy, kineticEnergy);
275  //nothing happens
276  if (tGamma == 0.) { return; }
277 
278  G4double totalEnergy = kineticEnergy + electron_mass_c2;
279  G4double finalEnergy = kineticEnergy - tGamma; // electron final energy
280 
281  //Sample gamma direction
282  G4ThreeVector gammaDirection =
283  GetAngularDistribution()->SampleDirection(aDynamicParticle,
284  totalEnergy-tGamma,
285  Z,
286  couple->GetMaterial());
287 
288  G4ThreeVector electronDirection = aDynamicParticle->GetMomentumDirection();
289 
290  //Update the incident particle
291  if (finalEnergy < 0.)
292  {
293  // Kinematic problem
294  tGamma = kineticEnergy;
296  }
297  else
298  {
299  G4double momentum = std::sqrt((totalEnergy + electron_mass_c2)*kineticEnergy);
300  G4double finalX = momentum*electronDirection.x() - tGamma*gammaDirection.x();
301  G4double finalY = momentum*electronDirection.y() - tGamma*gammaDirection.y();
302  G4double finalZ = momentum*electronDirection.z() - tGamma*gammaDirection.z();
303  G4double norm = 1./std::sqrt(finalX*finalX + finalY*finalY + finalZ*finalZ);
304 
305  fParticleChange->ProposeMomentumDirection(finalX*norm, finalY*norm, finalZ*norm);
307  }
308 
309  //Generate the bremsstrahlung gamma
311  gammaDirection, tGamma);
312  fvect->push_back(aGamma);
313 
314  if (verboseLevel > 1)
315  {
316  G4cout << "-----------------------------------------------------------" << G4endl;
317  G4cout << "Energy balance from G4LivermoreBremsstrahlung" << G4endl;
318  G4cout << "Incoming primary energy: " << kineticEnergy/keV << " keV" << G4endl;
319  G4cout << "-----------------------------------------------------------" << G4endl;
320  G4cout << "Outgoing primary energy: " << finalEnergy/keV << " keV" << G4endl;
321  G4cout << "Gamma ray " << tGamma/keV << " keV" << G4endl;
322  G4cout << "Total final state: " << (finalEnergy+tGamma)/keV << " keV" << G4endl;
323  G4cout << "-----------------------------------------------------------" << G4endl;
324  }
325  if (verboseLevel > 0)
326  {
327  G4double energyDiff = std::fabs(finalEnergy+tGamma-kineticEnergy);
328  if (energyDiff > 0.05*keV)
329  G4cout << "G4LivermoreBremsstrahlung WARNING: problem with energy conservation: "
330  << (finalEnergy+tGamma)/keV << " keV (final) vs. "
331  << kineticEnergy/keV << " keV (initial)" << G4endl;
332  }
333  return;
334 }
335 
336 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
337 /*
338 void
339 G4LivermoreBremsstrahlungModel::SetAngularGenerator(G4VBremAngularDistribution* distribution)
340 {
341  if(angularDistribution == distribution) return;
342  if(angularDistribution) delete angularDistribution;
343  angularDistribution = distribution;
344  angularDistribution->PrintGeneratorInformation();
345 }
346 */
347 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
348  /*
349 void G4LivermoreBremsstrahlungModel::SetAngularGenerator(const G4String& theGenName)
350 {
351  if(theGenName == generatorName) return;
352  if (theGenName == "tsai")
353  {
354  delete angularDistribution;
355  angularDistribution = new G4ModifiedTsai("TsaiGenerator");
356  generatorName = theGenName;
357  }
358  else if (theGenName == "2bn")
359  {
360  delete angularDistribution;
361  angularDistribution = new G4Generator2BN("2BNGenerator");
362  generatorName = theGenName;
363  }
364  else if (theGenName == "2bs")
365  {
366  delete angularDistribution;
367  angularDistribution = new G4Generator2BS("2BSGenerator");
368  generatorName = theGenName;
369  }
370  else
371  {
372  G4cout << "### G4LivermoreBremsstrahlungModel::SetAngularGenerator WARNING:"
373  << " generator <" << theGenName << "> is not known" << G4endl;
374  return;
375 
376  }
377 
378  angularDistribution->PrintGeneratorInformation();
379 }
380  */