Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
examples/extended/electromagnetic/TestEm18/src/RunAction.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 /// \file electromagnetic/TestEm18/src/RunAction.cc
27 /// \brief Implementation of the RunAction class
28 //
29 // $Id: RunAction.cc 67268 2013-02-13 11:38:40Z ihrivnac $
30 //
31 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
32 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
33 
34 #include "RunAction.hh"
35 #include "DetectorConstruction.hh"
36 #include "PrimaryGeneratorAction.hh"
37 #include "HistoManager.hh"
38 
39 #include "G4Run.hh"
40 #include "G4RunManager.hh"
41 #include "G4UnitsTable.hh"
42 #include "G4EmCalculator.hh"
43 
44 #include "Randomize.hh"
45 #include <iomanip>
46 
47 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
48 
50 :G4UserRunAction(),fDetector(det), fPrimary(kin), fHistoManager(0)
51 {
52  fHistoManager = new HistoManager();
53 }
54 
55 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
56 
58 {
59  delete fHistoManager;
60 }
61 
62 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
63 
64 void RunAction::BeginOfRunAction(const G4Run* run)
65 {
66  G4cout << "### Run " << run->GetRunID() << " start." << G4endl;
67 
68  //initialisation
69  //
70  fEnergyDeposit = 0.;
71 
72  fNbCharged = fNbNeutral = 0;
73  fEnergyCharged = fEnergyNeutral = 0.;
74  fEmin[0] = fEmin[1] = DBL_MAX;
75  fEmax[0] = fEmax[1] = 0.;
76 
77  fNbSteps = 0;
78  fTrackLength = 0.;
79 
80  //histograms
81  //
82  G4AnalysisManager* analysisManager = G4AnalysisManager::Instance();
83  if ( analysisManager->IsActive() ) {
84  analysisManager->OpenFile();
85  }
86 
87  // do not save Rndm status
90 }
91 
92 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
93 
94 void RunAction::EndOfRunAction(const G4Run* aRun)
95 {
96  G4int nbEvents = aRun->GetNumberOfEvent();
97  if (nbEvents == 0) return;
98 
99  G4Material* material = fDetector->GetMaterial();
100  G4double length = fDetector->GetSize();
101  G4double density = material->GetDensity();
102 
103  G4ParticleDefinition* particle = fPrimary->GetParticleGun()
105  G4String partName = particle->GetParticleName();
106  G4double ePrimary = fPrimary->GetParticleGun()->GetParticleEnergy();
107 
108  G4int prec = G4cout.precision(3);
109  G4cout << "\n ======================== run summary ======================\n";
110  G4cout << "\n The run was " << nbEvents << " " << partName << " of "
111  << G4BestUnit(ePrimary,"Energy") << " through "
112  << G4BestUnit(length,"Length") << " of "
113  << material->GetName() << " (density: "
114  << G4BestUnit(density,"Volumic Mass") << ")";
115  G4cout << "\n ===========================================================\n";
116  G4cout << G4endl;
117 
118  //save histograms
119  G4AnalysisManager* analysisManager = G4AnalysisManager::Instance();
120  if ( analysisManager->IsActive() ) {
121  analysisManager->Write();
122  analysisManager->CloseFile();
123  }
124 
125  if (particle->GetPDGCharge() == 0.) return;
126 
127  G4cout.precision(5);
128 
129  //track length
130  //
131  G4double trackLPerEvent = fTrackLength/nbEvents;
132  G4double nbStepPerEvent = double(fNbSteps)/nbEvents;
133  G4double stepSize = fTrackLength/fNbSteps;
134 
135  G4cout
136  << "\n TrackLength= "
137  << G4BestUnit(trackLPerEvent, "Length")
138  << "\t nb of steps= " << nbStepPerEvent
139  << " stepSize= " << G4BestUnit(stepSize, "Length")
140  << G4endl;
141 
142  //charged secondaries (ionization, direct pair production)
143  //
144  G4double energyPerEvent = fEnergyCharged/nbEvents;
145  G4double nbPerEvent = double(fNbCharged)/nbEvents;
146  G4double meanEkin = 0.;
147  if (fNbCharged) meanEkin = fEnergyCharged/fNbCharged;
148 
149  G4cout
150  << "\n d-rays : eLoss/primary= "
151  << G4BestUnit(energyPerEvent, "Energy")
152  << "\t nb of d-rays= " << nbPerEvent
153  << " <Tkin>= " << G4BestUnit(meanEkin, "Energy")
154  << " Tmin= " << G4BestUnit(fEmin[0], "Energy")
155  << " Tmax= " << G4BestUnit(fEmax[0], "Energy")
156  << G4endl;
157 
158  //neutral secondaries (bremsstrahlung, pixe)
159  //
160  energyPerEvent = fEnergyNeutral/nbEvents;
161  nbPerEvent = double(fNbNeutral)/nbEvents;
162  meanEkin = 0.;
163  if (fNbNeutral) meanEkin = fEnergyNeutral/fNbNeutral;
164 
165  G4cout
166  << "\n gamma : eLoss/primary= "
167  << G4BestUnit(energyPerEvent, "Energy")
168  << "\t nb of gammas= " << nbPerEvent
169  << " <Tkin>= " << G4BestUnit(meanEkin, "Energy")
170  << " Tmin= " << G4BestUnit(fEmin[1], "Energy")
171  << " Tmax= " << G4BestUnit(fEmax[1], "Energy")
172  << G4endl;
173 
174 
175  G4EmCalculator emCal;
176 
177  //local energy deposit
178  //
179  energyPerEvent = fEnergyDeposit/nbEvents;
180  //
181  G4double r0 = emCal.GetRangeFromRestricteDEDX(ePrimary,particle,material);
182  G4double r1 = r0 - trackLPerEvent;
183  G4double etry = ePrimary - energyPerEvent;
184  G4double efinal = 0.;
185  if (r1 > 0.) efinal = GetEnergyFromRestrictedRange(r1,particle,material,etry);
186  G4double dEtable = ePrimary - efinal;
187  G4double ratio = 0.;
188  if (dEtable > 0.) ratio = energyPerEvent/dEtable;
189 
190  G4cout
191  << "\n deposit : eLoss/primary= "
192  << G4BestUnit(energyPerEvent, "Energy")
193  << "\t <dEcut > table= "
194  << G4BestUnit(dEtable, "Energy")
195  << " ---> simul/reference= " << ratio
196  << G4endl;
197 
198  //total energy transferred
199  //
200  G4double energyTotal = fEnergyDeposit + fEnergyCharged + fEnergyNeutral;
201  energyPerEvent = energyTotal/nbEvents;
202  //
203  r0 = emCal.GetCSDARange(ePrimary,particle,material);
204  r1 = r0 - trackLPerEvent;
205  etry = ePrimary - energyPerEvent;
206  efinal = 0.;
207  if (r1 > 0.) efinal = GetEnergyFromCSDARange(r1,particle,material,etry);
208  dEtable = ePrimary - efinal;
209  ratio = 0.;
210  if (dEtable > 0.) ratio = energyPerEvent/dEtable;
211 
212  G4cout
213  << "\n total : eLoss/primary= "
214  << G4BestUnit(energyPerEvent, "Energy")
215  << "\t <dEfull> table= "
216  << G4BestUnit(dEtable, "Energy")
217  << " ---> simul/reference= " << ratio
218  << G4endl;
219 
220  G4cout.precision(prec);
221 
222  // show Rndm status
224 }
225 
226 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
227 
229  G4ParticleDefinition* particle, G4Material* material, G4double Etry)
230 {
231  G4EmCalculator emCal;
232 
233  G4double Energy = Etry, dE = 0., dEdx;
234  G4double r, dr;
235  G4double err = 1., errmax = 0.00001;
236  G4int iter = 0 , itermax = 10;
237  while (err > errmax && iter < itermax) {
238  iter++;
239  Energy -= dE;
240  r = emCal.GetRangeFromRestricteDEDX(Energy,particle,material);
241  dr = r - range;
242  dEdx = emCal.GetDEDX(Energy,particle,material);
243  dE = dEdx*dr;
244  err = std::abs(dE)/Energy;
245  }
246  if (iter == itermax) {
247  G4cout
248  << "\n ---> warning: RunAction::GetEnergyFromRestRange() did not converge"
249  << " Etry = " << G4BestUnit(Etry,"Energy")
250  << " Energy = " << G4BestUnit(Energy,"Energy")
251  << " err = " << err
252  << " iter = " << iter << G4endl;
253  }
254 
255  return Energy;
256 }
257 
258 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
259 
261  G4ParticleDefinition* particle, G4Material* material, G4double Etry)
262 {
263  G4EmCalculator emCal;
264 
265  G4double Energy = Etry, dE = 0., dEdx;
266  G4double r, dr;
267  G4double err = 1., errmax = 0.00001;
268  G4int iter = 0 , itermax = 10;
269  while (err > errmax && iter < itermax) {
270  iter++;
271  Energy -= dE;
272  r = emCal.GetCSDARange(Energy,particle,material);
273  dr = r - range;
274  dEdx = emCal.ComputeTotalDEDX(Energy,particle,material);
275  dE = dEdx*dr;
276  err = std::abs(dE)/Energy;
277  }
278  if (iter == itermax) {
279  G4cout
280  << "\n ---> warning: RunAction::GetEnergyFromCSDARange() did not converge"
281  << " Etry = " << G4BestUnit(Etry,"Energy")
282  << " Energy = " << G4BestUnit(Energy,"Energy")
283  << " err = " << err
284  << " iter = " << iter << G4endl;
285  }
286 
287  return Energy;
288 }
289 
290 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
G4double GetEnergyFromCSDARange(G4double, G4ParticleDefinition *, G4Material *, G4double)
const G4String & GetName() const
Definition: G4Material.hh:176
G4double GetDensity() const
Definition: G4Material.hh:178
G4double ComputeTotalDEDX(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *, G4double cut=DBL_MAX)
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
void SetRandomNumberStore(G4bool flag)
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
G4double GetDEDX(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *, const G4Region *r=0)
string material
Definition: eplot.py:19
G4double GetCSDARange(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *, const G4Region *r=0)
G4double density
Definition: TRTMaterials.hh:39
G4GLOB_DLL std::ostream G4cout
G4int GetNumberOfEvent() const
Definition: G4Run.hh:79
G4double GetRangeFromRestricteDEDX(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *, const G4Region *r=0)
G4int GetRunID() const
Definition: G4Run.hh:76
Definition: G4Run.hh:46
ExG4HbookAnalysisManager G4AnalysisManager
Definition: g4hbook_defs.hh:46
static void showEngineStatus()
Definition: Random.cc:203
static G4RunManager * GetRunManager()
Definition: G4RunManager.cc:74
G4ParticleDefinition * GetParticleDefinition() const
G4double GetEnergyFromRestrictedRange(G4double, G4ParticleDefinition *, G4Material *, G4double)
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
G4double GetPDGCharge() const
#define DBL_MAX
Definition: templates.hh:83
G4double GetParticleEnergy() const