Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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$
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 
40 #include "G4Run.hh"
41 #include "G4RunManager.hh"
42 #include "G4UnitsTable.hh"
43 #include "G4EmCalculator.hh"
44 
45 #include "Randomize.hh"
46 #include "G4SystemOfUnits.hh"
47 #include <iomanip>
48 
49 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
50 
53 :fDetector(det),fPhysics(phys),fKinematic(kin)
54 {
55  // Book predefined histograms
56  fHistoManager = new HistoManager();
57 }
58 
59 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
60 
62 {
63  delete fHistoManager;
64 }
65 
66 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
67 
68 void RunAction::BeginOfRunAction(const G4Run* aRun)
69 {
70  G4cout << "### Run " << aRun->GetRunID() << " start." << G4endl;
71 
72  // save Rndm status
75 
76  //initialize total energy deposit
77  //
78  fEdeposit = fEdeposit2 = 0.;
79 
80  //initialize track legth of primary
81  //
82  fTrackLen = fTrackLen2 = 0.;
83 
84  //initialize projected range
85  //
86  fProjRange = fProjRange2 = 0.;
87 
88  //initialize mean step size
89  //
90  fNbOfSteps = fNbOfSteps2 = 0; fStepSize = fStepSize2 = 0.;
91 
92  //initialize track status
93  //
94  fStatus[0] = fStatus[1] = fStatus[2] = 0;
95 
96  //normalized binwidth
97  //
98  for (G4int i=0; i< MaxAbsor; i++) {
99  fCsdaRange[i] = fXfrontNorm[i] = 0.;
100  }
101 
102  //histograms
103  //
104  G4AnalysisManager* analysisManager = G4AnalysisManager::Instance();
105  if ( analysisManager->IsActive() ) {
106  analysisManager->OpenFile();
107  }
108 
109  //set StepMax from histos 1 and 8
110  //e
111  G4double stepMax = DBL_MAX;
112  G4int ih = 1;
113  if (analysisManager->GetH1Activation(ih))
114  stepMax = analysisManager->GetH1Width(ih);
115  //
116  ih = 8;
117  G4ParticleDefinition* particle = fKinematic->GetParticleGun()
119  if (particle->GetPDGCharge() != 0.) {
120  G4double width = analysisManager->GetH1Width(ih);
121  if (width == 0.) width = 1.;
122  G4EmCalculator emCalculator;
123  G4double energy = fKinematic->GetParticleGun()->GetParticleEnergy();
124  G4int NbOfAbsor = fDetector->GetNbOfAbsor();
125  for (G4int i=1; i<= NbOfAbsor; i++) {
126  G4Material* material = fDetector->GetAbsorMaterial(i);
127  fCsdaRange[i] = emCalculator.GetCSDARange(energy,particle,material);
128  if (analysisManager->GetH1Activation(ih))
129  stepMax = std::min(stepMax, width*fCsdaRange[i]);
130  if (i>1) {
131  G4double thickness = fDetector->GetAbsorThickness(i-1);
132  fXfrontNorm[i] = fXfrontNorm[i-1] + thickness/fCsdaRange[i-1];
133  }
134  }
135  }
136  fPhysics->GetStepMaxProcess()->SetMaxStep2(stepMax);
137 }
138 
139 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
140 
141 void RunAction::EndOfRunAction(const G4Run* aRun)
142 {
143  std::ios::fmtflags mode = G4cout.flags();
144  G4cout.setf(std::ios::fixed,std::ios::floatfield);
145 
146  G4int NbofEvents = aRun->GetNumberOfEvent();
147  if (NbofEvents == 0) return;
148  G4double fNbofEvents = double(NbofEvents);
149 
150  //run conditions
151  //
152  G4ParticleDefinition* particle = fKinematic->GetParticleGun()
154  G4String partName = particle->GetParticleName();
155  G4double energy = fKinematic->GetParticleGun()->GetParticleEnergy();
156 
157  G4int NbOfAbsor = fDetector->GetNbOfAbsor();
158 
159  G4cout << "\n ======================== run summary ======================\n";
160 
161  G4int prec = G4cout.precision(2);
162 
163  G4cout
164  << "\n The run consists of " << NbofEvents << " "<< partName << " of "
165  << G4BestUnit(energy,"Energy")
166  << " through " << NbOfAbsor << " absorbers: \n";
167  for (G4int i=1; i<= NbOfAbsor; i++) {
168  G4Material* material = fDetector->GetAbsorMaterial(i);
169  G4double thickness = fDetector->GetAbsorThickness(i);
170  G4double density = material->GetDensity();
171  G4cout << std::setw(20) << G4BestUnit(thickness,"Length") << " of "
172  << material->GetName() << " (density: "
173  << G4BestUnit(density,"Volumic Mass") << ")" << G4endl;
174  }
175  G4cout << "\n ============================================================\n";
176 
177  //compute total energy deposit
178  //
179  fEdeposit /= NbofEvents; fEdeposit2 /= NbofEvents;
180  G4double rms = fEdeposit2 - fEdeposit*fEdeposit;
181  if (rms>0.) rms = std::sqrt(rms); else rms = 0.;
182 
183  G4cout.precision(3);
184  G4cout
185  << "\n Total Energy deposited = " << G4BestUnit(fEdeposit,"Energy")
186  << " +- " << G4BestUnit( rms,"Energy")
187  << G4endl;
188 
189  //compute track length of primary track
190  //
191  fTrackLen /= NbofEvents; fTrackLen2 /= NbofEvents;
192  rms = fTrackLen2 - fTrackLen*fTrackLen;
193  if (rms>0.) rms = std::sqrt(rms); else rms = 0.;
194 
195  G4cout.precision(3);
196  G4cout
197  << "\n Track length of primary track = " << G4BestUnit(fTrackLen,"Length")
198  << " +- " << G4BestUnit( rms,"Length");
199 
200  //compare with csda range
201  //
202  if (NbOfAbsor == 1) {
203  G4cout
204  << "\n Range from EmCalculator = " << G4BestUnit(fCsdaRange[1],"Length")
205  << " (from full dE/dx)" << G4endl;
206  }
207 
208  //compute projected range of primary track
209  //
210  fProjRange /= NbofEvents; fProjRange2 /= NbofEvents;
211  rms = fProjRange2 - fProjRange*fProjRange;
212  if (rms>0.) rms = std::sqrt(rms); else rms = 0.;
213 
214  G4cout
215  << "\n Projected range = " << G4BestUnit(fProjRange,"Length")
216  << " +- " << G4BestUnit( rms,"Length")
217  << G4endl;
218 
219  //nb of steps and step size of primary track
220  //
221  G4double fNbSteps = fNbOfSteps/fNbofEvents, fNbSteps2 = fNbOfSteps2/fNbofEvents;
222  rms = fNbSteps2 - fNbSteps*fNbSteps;
223  if (rms>0.) rms = std::sqrt(rms); else rms = 0.;
224 
225  G4cout.precision(2);
226  G4cout << "\n Nb of steps of primary track = " << fNbSteps << " +- " << rms;
227 
228  fStepSize /= NbofEvents; fStepSize2 /= NbofEvents;
229  rms = fStepSize2 - fStepSize*fStepSize;
230  if (rms>0.) rms = std::sqrt(rms); else rms = 0.;
231 
232  G4cout.precision(3);
233  G4cout
234  << "\t Step size= " << G4BestUnit(fStepSize,"Length")
235  << " +- " << G4BestUnit( rms,"Length")
236  << G4endl;
237 
238  //transmission coefficients
239  //
240  G4double absorbed = 100.*fStatus[0]/fNbofEvents;
241  G4double transmit = 100.*fStatus[1]/fNbofEvents;
242  G4double reflected = 100.*fStatus[2]/fNbofEvents;
243 
244  G4cout.precision(2);
245  G4cout
246  << "\n absorbed = " << absorbed << " %"
247  << " transmit = " << transmit << " %"
248  << " reflected = " << reflected << " %" << G4endl;
249 
250  // normalize histograms of longitudinal energy profile
251  //
252  G4AnalysisManager* analysisManager = G4AnalysisManager::Instance();
253  G4int ih = 1;
254  G4double binWidth = analysisManager->GetH1Width(ih);
255  G4double fac = (1./(NbofEvents*binWidth))*(mm/MeV);
256  analysisManager->ScaleH1(ih,fac);
257 
258  ih = 8;
259  binWidth = analysisManager->GetH1Width(ih);
260  fac = (1./(NbofEvents*binWidth))*(g/(MeV*cm2));
261  analysisManager->ScaleH1(ih,fac);
262 
263  // reset default formats
264  G4cout.setf(mode,std::ios::floatfield);
265  G4cout.precision(prec);
266 
267  // save histograms
268  if ( analysisManager->IsActive() ) {
269  analysisManager->Write();
270  analysisManager->CloseFile();
271  }
272 
273  // show Rndm status
275 }
276 
277 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......