Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/examples/extended/electromagnetic/TestEm7/src/RunAction.cc

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  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/TestEm7/src/RunAction.cc
 27 /// \brief Implementation of the RunAction class
 28 //
 29 //
 30 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 31 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 32 
 33 #include "RunAction.hh"
 34 
 35 #include "DetectorConstruction.hh"
 36 #include "PhysicsList.hh"
 37 #include "PrimaryGeneratorAction.hh"
 38 #include "StepMax.hh"
 39 
 40 #include "G4Run.hh"
 41 #include "G4RunManager.hh"
 42 #include "G4SystemOfUnits.hh"
 43 #include "G4UnitsTable.hh"
 44 #include "G4ios.hh"
 45 #include "Randomize.hh"
 46 
 47 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 48 
 49 RunAction::RunAction(DetectorConstruction* det, PhysicsList* phys, PrimaryGeneratorAction* kin)
 50   : G4UserRunAction(),
 51     fAnalysisManager(0),
 52     fDetector(det),
 53     fPhysics(phys),
 54     fKinematic(kin),
 55     fTallyEdep(new G4double[kMaxTally]),
 56     fProjRange(0.),
 57     fProjRange2(0.),
 58     fEdeptot(0.),
 59     fEniel(0.),
 60     fNbPrimarySteps(0),
 61     fRange(0)
 62 {
 63   // Book predefined histograms
 64   BookHisto();
 65 }
 66 
 67 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 68 
 69 RunAction::~RunAction()
 70 {
 71   delete[] fTallyEdep;
 72 }
 73 
 74 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 75 
 76 void RunAction::BeginOfRunAction(const G4Run* aRun)
 77 {
 78   G4cout << "### Run " << aRun->GetRunID() << " start." << G4endl;
 79 
 80   if (!fAnalysisManager) {
 81     BookHisto();
 82   }
 83 
 84   CLHEP::HepRandom::showEngineStatus();
 85 
 86   // initialize projected range, tallies, Ebeam, and book histograms
 87   //
 88   fNbPrimarySteps = 0;
 89   fRange = 0;
 90   fProjRange = fProjRange2 = 0.;
 91   fEdeptot = fEniel = 0.;
 92   for (G4int j = 0; j < kMaxTally; ++j) {
 93     fTallyEdep[j] = 0.;
 94   }
 95   fKinematic->ResetEbeamCumul();
 96 
 97   if (fAnalysisManager->IsActive()) {
 98     fAnalysisManager->OpenFile();
 99 
100     // histogram "1" is defined by the length of the target
101     // zoomed histograms are defined by UI command
102     G4double length = fDetector->GetAbsorSizeX();
103     G4double stepMax = fPhysics->GetStepMaxProcess()->GetMaxStep();
104     G4int nbBins = 100;
105     if (stepMax < DBL_MAX) {
106       G4int nb = (G4int)(0.5 + length / stepMax);
107       nbBins = std::min(std::max(nbBins, nb), 10000);
108     }
109     fAnalysisManager->SetH1(1, nbBins, 0., length, "mm");
110   }
111 }
112 
113 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
114 
115 void RunAction::EndOfRunAction(const G4Run* aRun)
116 {
117   G4int nbofEvents = aRun->GetNumberOfEvent();
118   if (nbofEvents == 0) return;
119 
120   // run conditions
121   //
122   const G4Material* material = fDetector->GetAbsorMaterial();
123   G4double density = material->GetDensity();
124 
125   G4String particle = fKinematic->GetParticleGun()->GetParticleDefinition()->GetParticleName();
126   G4double energy = fKinematic->GetParticleGun()->GetParticleEnergy();
127   G4cout << "\n The run consists of " << nbofEvents << " " << particle << " of "
128          << G4BestUnit(energy, "Energy") << " through "
129          << G4BestUnit(fDetector->GetAbsorSizeX(), "Length") << " of " << material->GetName()
130          << " (density: " << G4BestUnit(density, "Volumic Mass") << ")" << G4endl;
131 
132   // compute projected range and straggling
133   //
134   if (fRange > 0) {
135     fProjRange /= fRange;
136     fProjRange2 /= fRange;
137   }
138   G4double rms = fProjRange2 - fProjRange * fProjRange;
139   if (rms > 0.)
140     rms = std::sqrt(rms);
141   else
142     rms = 0.;
143 
144   G4double nstep = G4double(fNbPrimarySteps) / G4double(nbofEvents);
145 
146   G4cout.precision(6);
147   G4cout << "\n Projected Range= " << G4BestUnit(fProjRange, "Length")
148          << "   rms= " << G4BestUnit(rms, "Length") << G4endl;
149   G4cout << " Mean number of primary steps = " << nstep << G4endl;
150 
151   // compute energy deposition and niel
152   //
153   fEdeptot /= nbofEvents;
154   G4cout << " Total energy deposit= " << G4BestUnit(fEdeptot, "Energy") << G4endl;
155   fEniel /= nbofEvents;
156   G4cout << " niel energy deposit = " << G4BestUnit(fEniel, "Energy") << G4endl;
157 
158   // print dose in tallies
159   //
160   G4int tallyNumber = fDetector->GetTallyNumber();
161   if (tallyNumber > 0) {
162     G4double Ebeam = fKinematic->GetEbeamCumul();
163     G4cout << "\n---------------------------------------------------------\n";
164     G4cout << " Cumulated Doses : \tEdep      \tEdep/Ebeam \tDose" << G4endl;
165     for (G4int j = 0; j < tallyNumber; ++j) {
166       G4double Edep = fTallyEdep[j], ratio = 100 * Edep / Ebeam;
167       G4double tallyMass = fDetector->GetTallyMass(j);
168       G4double Dose = Edep / tallyMass;
169       G4cout << " tally " << j << ": \t \t" << G4BestUnit(Edep, "Energy") << "\t" << ratio
170              << " % \t" << G4BestUnit(Dose, "Dose") << G4endl;
171     }
172     G4cout << "\n---------------------------------------------------------\n";
173     G4cout << G4endl;
174   }
175 
176   if (fAnalysisManager->IsActive()) {
177     // normalize histograms
178     //
179     for (G4int j = 1; j < 3; ++j) {
180       G4double binWidth = fAnalysisManager->GetH1Width(j);
181       G4double fac = (mm / MeV) / (nbofEvents * binWidth);
182       fAnalysisManager->ScaleH1(j, fac);
183     }
184     fAnalysisManager->ScaleH1(3, 1. / nbofEvents);
185 
186     // save histograms
187     fAnalysisManager->Write();
188     fAnalysisManager->CloseFile();
189   }
190 
191   // show Rndm status
192   //
193   CLHEP::HepRandom::showEngineStatus();
194 }
195 
196 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
197 
198 void RunAction::BookHisto()
199 {
200   // Create or get analysis manager
201   // The choice of analysis technology is done via selection of a namespace
202   // in HistoManager.hh
203   fAnalysisManager = G4AnalysisManager::Instance();
204   fAnalysisManager->SetDefaultFileType("root");
205   fAnalysisManager->SetFileName("testem7");
206   fAnalysisManager->SetVerboseLevel(1);
207   fAnalysisManager->SetActivation(true);  // enable inactivation of histograms
208 
209   // Define histograms start values
210   const G4int kMaxHisto = 4;
211   const G4String id[] = {"h0", "h1", "h2", "h3"};
212   const G4String title[] = {
213     "dummy",  // 0
214     "Edep (MeV/mm) along absorber ",  // 1
215     "Edep (MeV/mm) along absorber zoomed",  // 2
216     "projectile range"  // 3
217   };
218 
219   // Default values (to be reset via /analysis/h1/set command)
220   G4int nbins = 100;
221   G4double vmin = 0.;
222   G4double vmax = 100.;
223 
224   // Create all histograms as inactivated
225   // as we have not yet set nbins, vmin, vmax
226   for (G4int k = 0; k < kMaxHisto; ++k) {
227     G4int ih = fAnalysisManager->CreateH1(id[k], title[k], nbins, vmin, vmax);
228     G4bool activ = false;
229     if (k == 1) activ = true;
230     fAnalysisManager->SetH1Activation(ih, activ);
231   }
232 }
233 
234 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
235