Geant4  10.03.p03
 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: RunAction.cc 68036 2013-03-13 14:13:45Z gcosmo $
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 "RunActionMessenger.hh"
38 
39 #include "G4Run.hh"
40 #include "G4RunManager.hh"
41 #include "G4UnitsTable.hh"
42 #include "G4SystemOfUnits.hh"
43 #include "G4ios.hh"
44 
45 #include "Histo.hh"
46 #include "G4EmCalculator.hh"
47 
48 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
49 
51  :G4UserRunAction(),
52  fHisto(0),fDetector(det),fKinematic(kin),fRunActionMessenger(0)
53 {
54  fVerboseLevel = 1;
55  fProjRange = fProjRange2 = fBinLength = fOffsetX = 0.;
56  fHisto = new Histo();
57  fHisto->SetFileName("monopole");
58  // create commands for interactive definition of the detector
59  fRunActionMessenger = new RunActionMessenger(this);
60 }
61 
62 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
63 
65 {
66  delete fHisto;
67  delete fRunActionMessenger;
68 }
69 
70 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
71 
72 void RunAction::BeginOfRunAction(const G4Run* aRun)
73 {
74  if(GetVerbose() > 0) {
75  G4cout << "### Run " << aRun->GetRunID() << " start." << G4endl;
76  }
77 
78  //initialize projected range, tallies, Ebeam, and book histograms
79  fProjRange = fProjRange2 = 0.;
80 
81  G4double length = fDetector->GetAbsorSizeX();
82  if(0.0 == fBinLength) { fBinLength = 5 * mm; }
83  if(fBinLength > fDetector->GetMaxStepSize()) {
84  fBinLength = fDetector->GetMaxStepSize();
85  }
86  fOffsetX = -0.5 * length;
87 
88  G4int nbBins = G4lrint(length / fBinLength);
89 
90 
91  // Create histograms
92  fHisto->Add1D("1","Edep (MeV/mm) along absorber (mm)", nbBins, 0, length, mm);
93  fHisto->Add1D("2","DEDX (MeV/mm) of proton", 100, -3., 7.);
94  fHisto->Add1D("3","DEDX (MeV/mm) of monopole", 100, -3., 7.);
95  fHisto->Add1D("4","Range(mm) of proton", 100, -3., 7., mm);
96  fHisto->Add1D("5","Range(mm) of monopole", 100, -3., 7., mm);
97 
98  fHisto->Book();
99 }
100 
101 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
102 
103 void RunAction::EndOfRunAction(const G4Run* aRun)
104 {
105  G4int nEvents = aRun->GetNumberOfEvent();
106  if (nEvents == 0) { return; }
107 
108  //run conditions
109  //
110  const G4Material* material = fDetector->GetAbsorMaterial();
111  G4double density = material->GetDensity();
112  G4String matName = material->GetName();
113  const G4ParticleDefinition* part =
114  fKinematic->GetParticleGun()->GetParticleDefinition();
115  G4String particle = part->GetParticleName();
116  G4double energy = fKinematic->GetParticleGun()->GetParticleEnergy();
117 
118  if(GetVerbose() > 0){
119  G4cout << "\n The run consists of " << nEvents << " "<< particle << " of "
120  << G4BestUnit(energy,"Energy") << " through "
121  << G4BestUnit(fDetector->GetAbsorSizeX(),"Length") << " of "
122  << matName << " (density: "
123  << G4BestUnit(density,"Volumic Mass") << ")" << G4endl;
124  };
125 
126  //compute projected range and straggling
127 
128  fProjRange /= nEvents; fProjRange2 /= nEvents;
129  G4double rms = fProjRange2 - fProjRange*fProjRange;
130  if (rms>0.) { rms = std::sqrt(rms); }
131  else { rms = 0.; }
132 
133  if(GetVerbose() > 0){
134  G4cout.precision(5);
135  G4cout << "\n projected Range= " << G4BestUnit(fProjRange, "Length")
136  << " rms= " << G4BestUnit(rms, "Length")
137  << G4endl;
138  };
139 
140  G4double ekin[100], dedxproton[100], dedxmp[100];
141  G4EmCalculator calc;
142  calc.SetVerbose(0);
143  G4int i;
144  for(i = 0; i < 100; ++i) {
145  ekin[i] = std::pow(10., 0.1*G4double(i)) * keV;
146  dedxproton[i] =
147  calc.ComputeElectronicDEDX(ekin[i], "proton", matName);
148  dedxmp[i] =
149  calc.ComputeElectronicDEDX(ekin[i], "monopole", matName);
150  }
151 
152  if(GetVerbose() > 0){
153  G4cout << "### Stopping Powers" << G4endl;
154  for(i=0; i<100; i++) {
155  G4cout << " E(MeV)= " << ekin[i] << " dedxp= " << dedxproton[i]
156  << " dedxmp= " << dedxmp[i]
157  << G4endl;
158  }
159  }
160  G4cout << "### End of stopping power table" << G4endl;
161 
162  // normalize histogram
163  G4double fac = (mm/MeV) / (nEvents * fBinLength);
164  fHisto->ScaleH1(0,fac);
165 
166  if(GetVerbose() > 0){
167  G4cout << "Range table for " << matName << G4endl;
168  }
169 
170  for(i=0; i<100; ++i) {
171  G4double e = std::log10(ekin[i] / MeV) + 0.05;
172  fHisto->Fill(1, e, dedxproton[i]);
173  fHisto->Fill(2, e, dedxmp[i]);
174  fHisto->Fill(3, e, std::log10(calc.GetRange(ekin[i],"proton",matName)/mm));
175  fHisto->Fill(4, e, std::log10(calc.GetRange(ekin[i],"monopole",matName)/mm));
176  }
177 
178  // save and clean histo
179  fHisto->Save();
180 
182 }
183 
184 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
185 
187 {
188  if(GetVerbose() > 1) {
189  G4cout << "FillHisto " << ih << " x=" << x << " weight= " << weight
190  << G4endl;
191  }
192  fHisto->Fill(ih, x, weight);
193 }
194 
195 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
The primary generator action class with particle gun.
static constexpr double mm
Definition: G4SIunits.hh:115
void BeginOfRunAction(const G4Run *)
Definition: RunAction.cc:57
G4double GetRange(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *, const G4Region *r=nullptr)
Definition: Histo.hh:56
const G4String & GetName() const
Definition: G4Material.hh:178
G4double GetDensity() const
Definition: G4Material.hh:180
G4double ComputeElectronicDEDX(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *mat, G4double cut=DBL_MAX)
tuple x
Definition: test.py:50
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
void Book()
Definition: Histo.cc:78
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
void Fill(G4int, G4double, G4double)
Definition: Histo.cc:221
G4Material * GetAbsorMaterial(G4int i)
string material
Definition: eplot.py:19
G4GLOB_DLL std::ostream G4cout
G4int GetNumberOfEvent() const
Definition: G4Run.hh:79
void FillHisto(G4int id, G4double x, G4double weight=1.0)
Definition: RunAction.cc:186
void ScaleH1(G4int, G4double)
Definition: Histo.cc:238
void EndOfRunAction(const G4Run *)
Definition: RunAction.cc:260
G4int GetVerbose()
Definition: RunAction.hh:65
G4int GetRunID() const
Definition: G4Run.hh:76
Definition: G4Run.hh:46
void Add1D(const G4String &, const G4String &, G4int nb, G4double x1, G4double x2, G4double u=1.)
Definition: Histo.cc:161
static void showEngineStatus()
Definition: Random.cc:303
~RunAction()
Definition: RunAction.cc:52
int G4lrint(double ad)
Definition: templates.hh:163
G4ParticleGun * GetParticleGun()
G4double energy(const ThreeVector &p, const G4double m)
G4ParticleDefinition * GetParticleDefinition() const
static const G4double fac
Detector construction class to define materials and geometry.
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
void SetVerbose(G4int val)
double G4double
Definition: G4Types.hh:76
static constexpr double keV
Definition: G4SIunits.hh:216
void Save()
Definition: Histo.cc:137
G4double GetParticleEnergy() const
void SetFileName(const G4String &)
Definition: Histo.cc:350