source: trunk/source/processes/electromagnetic/standard/test/BremLPMTest.cc@ 1228

Last change on this file since 1228 was 1199, checked in by garnier, 16 years ago

nvx fichiers dans CVS

File size: 6.6 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: BremLPMTest.cc,v 1.2 2008/08/22 08:16:17 schaelic Exp $
27// GEANT4 tag $Name: geant4-09-03-cand-01 $
28//
29// ------------------------------------------------------------
30// Test Routine for eBremsstrahlung-Models
31// -- including LPM effect
32// ------------------------------------------------------------
33// run
34// make -f GNUmakefile.root G4TARGET=BremLPMTest
35// BremLPMTest
36// this creates a file eBremRel01.root which can be read using e.g.
37// python readBremLPM.py
38// ------------------------------------------------------------
39//
40// History
41// 21.08.08 basic test created (A.Schaelicke)
42//
43//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
44//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
45
46#include "TFile.h"
47#include "TGraph.h"
48#include "TH1F.h"
49#include "TTree.h"
50
51#include "G4NistManager.hh"
52#include "G4Material.hh"
53#include "G4UnitsTable.hh"
54#include "G4Electron.hh"
55
56// ------------------------------------------------------------
57// nasty trick to access private or protected functions
58// use only in tests and with care !!!
59#define protected public
60// ------------------------------------------------------------
61
62#include "G4eBremsstrahlungHEModel.hh"
63#include "G4eBremsstrahlungRelModel.hh"
64#include "G4eBremsstrahlungModel.hh"
65
66
67using namespace std;
68
69//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
70
71
72void CalcCrossSection()
73{
74 // initialize models
75 G4eBremsstrahlungModel * model1 = new G4eBremsstrahlungModel();
76 model1->SetLPMflag(false);
77
78 G4eBremsstrahlungRelModel * modelR = new G4eBremsstrahlungRelModel(); // rel. version
79 // modelR->SetLPMflag(false);
80
81 // define energy range and cut
82 const G4int nmax=80;
83 G4double emin=1.e2*MeV;
84 G4double emax=1.e6*MeV;
85 G4double kinEs[nmax+1];
86
87 G4double cut=10.*keV;
88 for (int j=0; j<=nmax; ++j)
89 kinEs[j]=exp(G4double(j)/nmax*log(emax/emin))*emin;
90
91
92 // writing information to root file
93 // Z , cut, emin, emax
94 G4int currentZ=0;
95
96 TTree tree("info","info");
97 tree.Branch("Z",&currentZ,"Z/I");
98 tree.Branch("cut",&cut,"cut/D");
99 tree.Branch("emin",&emin,"emin/D");
100 tree.Branch("emax",&emax,"emax/D");
101
102 // setup test materials
103 const G4int nElements=6;
104 G4int theZ[]={1,8,29,47,82,92};
105
106 G4Element* els[nElements];
107 G4Material* mats[nElements];
108 for (G4int i=0; i<nElements; i++) {
109 G4Element * el = G4NistManager::Instance()->FindOrBuildElement(theZ[i]);
110 std::vector<G4String> names = G4NistManager::Instance()->GetNistMaterialNames();
111 G4Material* mat = G4NistManager::Instance()->FindOrBuildMaterial(names[theZ[i]-1]);
112 els[i]=el;
113 mats[i]=mat;
114 }
115
116
117 // fill plots
118 G4double cross[nmax+1], cross1[nmax+1];
119 G4double xDEDX[nmax+1], xDEDX1[nmax+1];
120
121 for (G4int i=0; i<nElements; i++ ) {
122 currentZ=theZ[i];
123 tree.Fill();
124 G4cout<<" Z="<<theZ[i]<<" ("<<mats[i]->GetName()<<","<<els[i]->GetName()<<")"<<G4endl;
125 const G4double* theAtomicNumDensityVector = mats[i]->GetAtomicNumDensityVector();
126 G4double dndV = 0.0;
127 for (size_t j=0; j<mats[i]->GetNumberOfElements(); j++)
128 dndV += theAtomicNumDensityVector[j];
129
130 for (int j=0; j<=nmax; ++j) {
131
132 modelR->SetupForMaterial(0,mats[i]);
133 modelR->G4VEmModel::SetCurrentElement(els[i]);
134// cross[j] =
135// modelR->ComputeCrossSectionPerAtom( G4Electron::Electron(),
136// kinEs[j], theZ[i], dum, cut)/barn;
137// cross1[j] =
138// model1->ComputeCrossSectionPerAtom( G4Electron::Electron(),
139// kinEs[j], theZ[i], dum, cut)/barn;
140 xDEDX[j]= modelR->ComputeDEDXPerVolume(mats[i], G4Electron::Electron(), kinEs[j], cut)/dndV/barn;
141 xDEDX1[j] = model1->ComputeDEDXPerVolume(mats[i], G4Electron::Electron(), kinEs[j], cut)/dndV/barn;
142 cross[j] = modelR->CrossSectionPerVolume(mats[i], G4Electron::Electron(), kinEs[j], cut, kinEs[j])/dndV/barn;
143 cross1[j] = model1->CrossSectionPerVolume(mats[i], G4Electron::Electron(), kinEs[j], cut, kinEs[j])/dndV/barn;
144 }
145 TGraph gR(nmax,kinEs,cross);
146 gR.Draw("Alp");
147 gR.SetLineColor(2);
148 gR.GetHistogram()->SetXTitle("E_{kin}");
149 gR.GetHistogram()->SetYTitle("#sigma(E_{kin})");
150 gR.Write("modelR");
151 TGraph g1(nmax,kinEs,cross1);
152 g1.Draw("lp");
153 g1.SetLineColor(4);
154 gR.GetHistogram()->SetXTitle("E_{kin}");
155 gR.GetHistogram()->SetYTitle("#sigma(E_{kin})");
156 g1.Write("model1");
157
158 TGraph xR(nmax,kinEs,xDEDX);
159 xR.Draw("Alp");
160 xR.SetLineColor(2);
161 xR.GetHistogram()->SetXTitle("E_{kin}");
162 xR.GetHistogram()->SetYTitle("dE/dx");
163 xR.Write("xDEDXR");
164 TGraph x1(nmax,kinEs,xDEDX1);
165 x1.Draw("lp");
166 x1.SetLineColor(4);
167 xR.GetHistogram()->SetXTitle("E_{kin}");
168 xR.GetHistogram()->SetYTitle("dE/dx");
169 x1.Write("xDEDX1");
170 }
171 tree.Write();
172}
173
174int main(int argc,char** argv) {
175 G4UnitDefinition::BuildUnitsTable();
176
177 TFile tree("eBremRel01.root","RECREATE","eBremRel");
178
179 // make some tests
180 CalcCrossSection();
181
182 tree.Write();
183 tree.Close();
184
185 return 0;
186}
187
188//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
Note: See TracBrowser for help on using the repository browser.