Geant4  10.00.p02
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 76258 2013-11-08 11:36:51Z gcosmo $
30 //
31 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
32 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
33 
34 #include "RunAction.hh"
35 #include "PhysicsList.hh"
36 #include "StepMax.hh"
37 #include "PrimaryGeneratorAction.hh"
38 #include "HistoManager.hh"
39 #include "Run.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 
54 :G4UserRunAction(),fDetector(det),fPhysics(phys),fKinematic(kin),fRun(0),
55  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 
71 {
72  fRun = new Run(fDetector);
73  return fRun;
74 }
75 
76 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
77 
79 {
80  // save Rndm status
82  G4Random::showEngineStatus();
83 
84  //histograms
85  //
86  G4AnalysisManager* analysisManager = G4AnalysisManager::Instance();
87  if ( analysisManager->IsActive() ) {
88  analysisManager->OpenFile();
89  }
90 
91  //set StepMax from histos 1 and 8
92  //e
93  G4double stepMax = DBL_MAX;
94  G4int ih = 1;
95  if (analysisManager->GetH1Activation(ih)) {
96  stepMax = analysisManager->GetH1Width(ih);
97  }
98 
99  if ( ! fKinematic ) return;
100 
101  //
102  ih = 8;
105  if (particle->GetPDGCharge() != 0.) {
106  G4double width = analysisManager->GetH1Width(ih);
107  if (width == 0.) width = 1.;
108  G4EmCalculator emCalculator;
110  G4int nbOfAbsor = fDetector->GetNbOfAbsor();
111  for (G4int i=1; i<= nbOfAbsor; i++) {
112  G4Material* material = fDetector->GetAbsorMaterial(i);
113  G4double newCsdaRange
114  = emCalculator.GetCSDARange(energy,particle,material);
115  fRun->SetCsdaRange(i, newCsdaRange);
116  if (analysisManager->GetH1Activation(ih))
117  stepMax = std::min(stepMax, width*newCsdaRange);
118  if (i>1) {
119  G4double thickness = fDetector->GetAbsorThickness(i-1);
120  G4double xfrontNorm = fRun->GetXfrontNorm(i-1);
121  G4double csdaRange = fRun->GetCsdaRange(i-1);
122  G4double newXfrontNorm = xfrontNorm + thickness/csdaRange;
123  fRun->SetXfrontNorm(i, newXfrontNorm);
124  }
125  }
126  }
127  fPhysics->GetStepMaxProcess()->SetMaxStep2(stepMax);
128 }
129 
130 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
131 
132 void RunAction::EndOfRunAction(const G4Run* run)
133 {
134  std::ios::fmtflags mode = G4cout.flags();
135  G4cout.setf(std::ios::fixed,std::ios::floatfield);
136  G4int prec = G4cout.precision(2);
137 
138  G4int nbofEvents = run->GetNumberOfEvent();
139 
140  if ( fKinematic && nbofEvents ) {
141 
142  //run conditions
143  //
146  G4String partName = particle->GetParticleName();
148 
149  G4int nbOfAbsor = fDetector->GetNbOfAbsor();
150 
151  G4cout << "\n ======================== run summary =====================\n";
152 
153  G4cout
154  << "\n The run consists of " << nbofEvents << " "<< partName << " of "
155  << G4BestUnit(energy,"Energy")
156  << " through " << nbOfAbsor << " absorbers: \n";
157  for (G4int i=1; i<= nbOfAbsor; i++) {
158  G4Material* material = fDetector->GetAbsorMaterial(i);
159  G4double thickness = fDetector->GetAbsorThickness(i);
160  G4double density = material->GetDensity();
161  G4cout << std::setw(20) << G4BestUnit(thickness,"Length") << " of "
162  << material->GetName() << " (density: "
163  << G4BestUnit(density,"Volumic Mass") << ")" << G4endl;
164  }
165  G4cout << "\n ==========================================================\n";
166  }
167 
168  //compute and print total energy deposit
169  //
170  if ( isMaster ) {
172 
173  // normalize histograms of longitudinal energy profile
174  //
175  G4AnalysisManager* analysisManager = G4AnalysisManager::Instance();
176  G4int ih = 1;
177  G4double binWidth = analysisManager->GetH1Width(ih);
178  G4double fac = (1./(nbofEvents*binWidth))*(mm/MeV);
179  analysisManager->ScaleH1(ih,fac);
180 
181  ih = 8;
182  binWidth = analysisManager->GetH1Width(ih);
183  fac = (1./(nbofEvents*binWidth))*(g/(MeV*cm2));
184  analysisManager->ScaleH1(ih,fac);
185  }
186 
187  // reset default formats
188  G4cout.setf(mode,std::ios::floatfield);
189  G4cout.precision(prec);
190 
191  // save histograms
192  G4AnalysisManager* analysisManager = G4AnalysisManager::Instance();
193  if ( analysisManager->IsActive() ) {
194  analysisManager->Write();
195  analysisManager->CloseFile();
196  }
197 
198  // show Rndm status
199  G4Random::showEngineStatus();
200 }
201 
202 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
G4double GetCsdaRange(G4int i)
Definition: Run.hh:66
G4StepLimiter * GetStepMaxProcess()
Definition: PhysicsList.hh:64
static const double MeV
Definition: G4SIunits.hh:193
static const double cm2
Definition: G4SIunits.hh:107
void BeginOfRunAction(const G4Run *)
Definition: RunAction.cc:57
static const G4double fac
const G4String & GetName() const
Definition: G4Material.hh:176
void SetCsdaRange(G4int i, G4double value)
Definition: Run.hh:63
G4double GetDensity() const
Definition: G4Material.hh:178
#define width
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
virtual G4Run * GenerateRun()
Definition: RunAction.cc:65
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)
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
void SetXfrontNorm(G4int i, G4double value)
Definition: Run.hh:64
void EndOfRunAction(const G4Run *)
Definition: RunAction.cc:260
Definition: G4Run.hh:46
PhysicsList * fPhysics
Definition: RunAction.hh:61
ExG4HbookAnalysisManager G4AnalysisManager
Definition: g4hbook_defs.hh:46
~RunAction()
Definition: RunAction.cc:52
G4double GetAbsorThickness(G4int i)
void ComputeStatistics()
Definition: Run.cc:97
G4double GetXfrontNorm(G4int i)
Definition: Run.hh:67
G4ParticleGun * GetParticleGun()
G4double energy(const ThreeVector &p, const G4double m)
static const double g
Definition: G4SIunits.hh:162
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
double G4double
Definition: G4Types.hh:76
Definition: Run.hh:46
G4double GetPDGCharge() const
#define DBL_MAX
Definition: templates.hh:83
static const double mm
Definition: G4SIunits.hh:102
G4double GetParticleEnergy() const
Run * fRun
Definition: RunAction.hh:63