Geant4  10.00.p01
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 //
28 //
29 // $Id: RunAction.cc 69108 2013-04-18 13:10:10Z gcosmo $
30 //
31 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
32 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
33 
34 #include "RunAction.hh"
35 #include "DetectorConstruction.hh"
36 #include "PhysicsList.hh"
37 #include "StepMax.hh"
38 #include "PrimaryGeneratorAction.hh"
39 #include "HistoManager.hh"
40 
41 #include "G4Run.hh"
42 #include "G4RunManager.hh"
43 #include "G4UnitsTable.hh"
44 #include "G4EmCalculator.hh"
45 
46 #include "Randomize.hh"
47 #include "G4SystemOfUnits.hh"
48 #include <iomanip>
49 
50 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
51 
55  fDetector(det),fPhysics(phys),fKinematic(kin),fHistoManager(0)
56 {
57  // Book predefined histograms
59  }
60 
61 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
62 
64 {
65  delete fHistoManager;
66 }
67 
68 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
69 
70 void RunAction::BeginOfRunAction(const G4Run* aRun)
71 {
72  G4cout << "### Run " << aRun->GetRunID() << " start." << G4endl;
73 
74  // save Rndm status
76  CLHEP::HepRandom::showEngineStatus();
77 
78  //initialize total energy deposit
79  //
80  fEdeposit = fEdeposit2 = 0.;
81 
82  //initialize track legth of primary
83  //
84  fTrackLen = fTrackLen2 = 0.;
85 
86  //initialize projected range
87  //
88  fProjRange = fProjRange2 = 0.;
89 
90  //initialize mean step size
91  //
93 
94  //get fCsdaRange from EmCalculator
95  //
96  G4EmCalculator emCalculator;
97  G4Material* material = fDetector->GetAbsorMaterial();
102  if (particle->GetPDGCharge() != 0.)
103  fCsdaRange = emCalculator.GetCSDARange(energy,particle,material);
104 
105  //histograms
106  //
107  G4AnalysisManager* analysisManager = G4AnalysisManager::Instance();
108  if ( analysisManager->IsActive() ) {
109  analysisManager->OpenFile();
110  }
111 
112  //set StepMax from histos 1 and 8
113  //
114  G4double stepMax = fCsdaRange;
115  G4int ih = 1;
116  if (analysisManager->GetH1Activation(ih))
117  stepMax = analysisManager->GetH1Width(ih);
118  ih = 8;
119  if (analysisManager->GetH1Activation(ih)) {
120  G4double width = analysisManager->GetH1Width(ih);
121  stepMax = std::min(stepMax, width*fCsdaRange);
122  }
123  fPhysics->GetStepMaxProcess()->SetMaxStep2(stepMax);
124  G4cout << "\n---> stepMax from histos 1 and 8 = "
125  << G4BestUnit(stepMax,"Length") << G4endl;
126 }
127 
128 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
129 
130 void RunAction::EndOfRunAction(const G4Run* aRun)
131 {
132  std::ios::fmtflags mode = G4cout.flags();
133  G4cout.setf(std::ios::fixed,std::ios::floatfield);
134 
135  G4int NbofEvents = aRun->GetNumberOfEvent();
136  if (NbofEvents == 0) return;
137  G4double fNbofEvents = double(NbofEvents);
138 
139  //run conditions
140  //
141  G4Material* material = fDetector->GetAbsorMaterial();
142  G4double density = material->GetDensity();
143 
146  G4String partName = particle->GetParticleName();
148 
149  G4cout << "\n ======================== run summary ======================\n";
150 
151  G4int prec = G4cout.precision(2);
152 
153  G4cout << "\n The run consists of " << NbofEvents << " "<< partName << " of "
154  << G4BestUnit(energy,"Energy") << " through "
155  << G4BestUnit(fDetector->GetAbsorRadius(),"Length") << " of "
156  << material->GetName() << " (density: "
157  << G4BestUnit(density,"Volumic Mass") << ")" << G4endl;
158 
159  G4cout << "\n ============================================================\n";
160 
161  //compute total energy deposit
162  //
163  fEdeposit /= NbofEvents; fEdeposit2 /= NbofEvents;
165  if (rms>0.) rms = std::sqrt(rms); else rms = 0.;
166 
167  G4cout.precision(3);
168  G4cout
169  << "\n Total Energy deposited = " << G4BestUnit(fEdeposit,"Energy")
170  << " +- " << G4BestUnit( rms,"Energy")
171  << G4endl;
172 
173  //compute track length of primary track
174  //
175  fTrackLen /= NbofEvents; fTrackLen2 /= NbofEvents;
176  rms = fTrackLen2 - fTrackLen*fTrackLen;
177  if (rms>0.) rms = std::sqrt(rms); else rms = 0.;
178 
179  G4cout.precision(3);
180  G4cout
181  << "\n Track length of primary track = " << G4BestUnit(fTrackLen,"Length")
182  << " +- " << G4BestUnit( rms,"Length");
183 
184  //compare with csda range
185  //
186  //G4EmCalculator emCalculator;
187  //G4double fCsdaRange = 0.;
188  //if (particle->GetPDGCharge() != 0.)
189  // fCsdaRange = emCalculator.GetCSDARange(energy,particle,material);
190  G4cout
191  << "\n Range from EmCalculator = " << G4BestUnit(fCsdaRange,"Length")
192  << " (from full dE/dx)" << G4endl;
193 
194  //compute projected range of primary track
195  //
196  fProjRange /= NbofEvents; fProjRange2 /= NbofEvents;
198  if (rms>0.) rms = std::sqrt(rms); else rms = 0.;
199 
200  G4cout
201  << "\n Projected range = " << G4BestUnit(fProjRange,"Length")
202  << " +- " << G4BestUnit( rms,"Length")
203  << G4endl;
204 
205  //nb of steps and step size of primary track
206  //
207  G4double fNbSteps = fNbOfSteps/fNbofEvents, fNbSteps2 = fNbOfSteps2/fNbofEvents;
208  rms = fNbSteps2 - fNbSteps*fNbSteps;
209  if (rms>0.) rms = std::sqrt(rms); else rms = 0.;
210 
211  G4cout.precision(2);
212  G4cout << "\n Nb of steps of primary track = " << fNbSteps << " +- " << rms;
213 
214  fStepSize /= NbofEvents; fStepSize2 /= NbofEvents;
215  rms = fStepSize2 - fStepSize*fStepSize;
216  if (rms>0.) rms = std::sqrt(rms); else rms = 0.;
217 
218  G4cout.precision(3);
219  G4cout
220  << "\t Step size= " << G4BestUnit(fStepSize,"Length")
221  << " +- " << G4BestUnit( rms,"Length")
222  << G4endl;
223 
224  // normalize histogram of longitudinal energy profile
225  //
226  G4AnalysisManager* analysisManager = G4AnalysisManager::Instance();
227  G4int ih = 1;
228  G4double binWidth = analysisManager->GetH1Width(ih);
229  G4double fac = (1./(NbofEvents*binWidth))*(mm/MeV);
230  analysisManager->ScaleH1(ih,fac);
231 
232  // normalize histogram d(E/E0)/d(r/r0)
233  //
234  ih = 8;
235  binWidth = analysisManager->GetH1Width(ih);
236  fac = 1./(NbofEvents*binWidth*energy);
237  analysisManager->ScaleH1(ih,fac);
238 
239  // reset default formats
240  G4cout.setf(mode,std::ios::floatfield);
241  G4cout.precision(prec);
242 
243  // save histograms
244  if ( analysisManager->IsActive() ) {
245  analysisManager->Write();
246  analysisManager->CloseFile();
247  }
248 
249  // show Rndm status
250  CLHEP::HepRandom::showEngineStatus();
251 }
252 
253 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
G4double fCsdaRange
Definition: RunAction.hh:80
G4StepLimiter * GetStepMaxProcess()
Definition: PhysicsList.hh:64
static const double MeV
Definition: G4SIunits.hh:193
void BeginOfRunAction(const G4Run *)
Definition: RunAction.cc:57
static const G4double fac
G4double fEdeposit2
Definition: RunAction.hh:74
const G4String & GetName() const
Definition: G4Material.hh:176
G4double fTrackLen2
Definition: RunAction.hh:75
G4double GetDensity() const
Definition: G4Material.hh:178
#define width
G4double fProjRange
Definition: RunAction.hh:76
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
int G4int
Definition: G4Types.hh:78
DetectorConstruction * fDetector
Definition: RunAction.hh:87
const G4String & GetParticleName() const
HistoManager * fHistoManager
Definition: RunAction.hh:60
G4Material * GetAbsorMaterial(G4int i)
G4long fNbSteps
Definition: RunAction.hh:94
G4double GetCSDARange(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *, const G4Region *r=0)
G4double density
Definition: TRTMaterials.hh:39
PrimaryGeneratorAction * fKinematic
Definition: RunAction.hh:62
static const double prec
Definition: RanecuEngine.cc:51
G4GLOB_DLL std::ostream G4cout
G4int GetNumberOfEvent() const
Definition: G4Run.hh:79
G4double fStepSize
Definition: RunAction.hh:78
void EndOfRunAction(const G4Run *)
Definition: RunAction.cc:260
G4int GetRunID() const
Definition: G4Run.hh:76
Definition: G4Run.hh:46
G4double fTrackLen
Definition: RunAction.hh:75
PhysicsList * fPhysics
Definition: RunAction.hh:61
ExG4HbookAnalysisManager G4AnalysisManager
Definition: g4hbook_defs.hh:46
~RunAction()
Definition: RunAction.cc:52
G4int fNbOfSteps2
Definition: RunAction.hh:77
G4ParticleGun * GetParticleGun()
G4double energy(const ThreeVector &p, const G4double m)
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
G4ParticleDefinition * GetParticleDefinition() const
Detector construction class to demonstrate various ways of placement.
#define G4endl
Definition: G4ios.hh:61
G4int fNbOfSteps
Definition: RunAction.hh:77
double G4double
Definition: G4Types.hh:76
G4double fStepSize2
Definition: RunAction.hh:78
G4double GetPDGCharge() const
#define DBL_MAX
Definition: templates.hh:83
static const double mm
Definition: G4SIunits.hh:102
G4double fEdeposit
Definition: RunAction.hh:74
G4double GetParticleEnergy() const
G4double fProjRange2
Definition: RunAction.hh:76