Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
examples/extended/electromagnetic/TestEm7/src/RunAction.cc
이 파일의 문서화 페이지로 가기
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 101250 2016-11-10 08:54:02Z 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 
40 #include "G4Run.hh"
41 #include "G4RunManager.hh"
42 #include "G4UnitsTable.hh"
43 #include "G4SystemOfUnits.hh"
44 #include "G4ios.hh"
45 
46 #include "Randomize.hh"
47 
48 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
49 
52  : G4UserRunAction(),
53  fAnalysisManager(0), fDetector(det), fPhysics(phys), fKinematic(kin),
54  fTallyEdep(new G4double[kMaxTally]), fProjRange(0.), fProjRange2(0.),
55  fEdeptot(0.), fEniel(0.), fNbPrimarySteps(0), fRange(0)
56 {
57  // Book predefined histograms
58  BookHisto();
59 }
60 
61 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
62 
64 {
65  delete [] fTallyEdep;
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  if(!fAnalysisManager) { BookHisto(); }
75 
77 
78  //initialize projected range, tallies, Ebeam, and book histograms
79  //
80  fNbPrimarySteps = 0;
81  fRange = 0;
82  fProjRange = fProjRange2 = 0.;
83  fEdeptot = fEniel = 0.;
84  for (G4int j=0; j<kMaxTally; ++j) { fTallyEdep[j] = 0.; }
86 
87  if (fAnalysisManager->IsActive()) {
88  fAnalysisManager->OpenFile();
89 
90  // histogram "1" is defined by the length of the target
91  // zoomed histograms are defined by UI command
92  G4double length = fDetector->GetAbsorSizeX();
93  G4double stepMax = fPhysics->GetStepMaxProcess()->GetMaxStep();
94  G4int nbmin = 100;
95  G4int nbBins = (G4int)(0.5 + length/stepMax);
96  if (nbBins < nbmin) nbBins = nbmin;
97  fAnalysisManager->SetH1(1, nbBins, 0., length, "mm");
98  }
99 }
100 
101 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
102 
103 void RunAction::EndOfRunAction(const G4Run* aRun)
104 {
105  G4int nbofEvents = aRun->GetNumberOfEvent();
106  if (nbofEvents == 0) return;
107 
108  //run conditions
109  //
110  const G4Material* material = fDetector->GetAbsorMaterial();
111  G4double density = material->GetDensity();
112 
114  ->GetParticleName();
116  G4cout << "\n The run consists of " << nbofEvents << " "<< particle << " of "
117  << G4BestUnit(energy,"Energy") << " through "
118  << G4BestUnit(fDetector->GetAbsorSizeX(),"Length") << " of "
119  << material->GetName() << " (density: "
120  << G4BestUnit(density,"Volumic Mass") << ")" << G4endl;
121 
122  //compute projected range and straggling
123  //
124  if(fRange > 0) {
125  fProjRange /= fRange;
126  fProjRange2 /= fRange;
127  }
129  if (rms>0.) rms = std::sqrt(rms); else rms = 0.;
130 
131  G4double nstep = G4double(fNbPrimarySteps)/G4double(nbofEvents);
132 
133  G4cout.precision(6);
134  G4cout << "\n Projected Range= "<< G4BestUnit(fProjRange,"Length")
135  << " rms= " << G4BestUnit( rms,"Length")
136  << G4endl;
137  G4cout << " Mean number of primary steps = "<< nstep << G4endl;
138 
139  //compute energy deposition and niel
140  //
141  fEdeptot /= nbofEvents;
142  G4cout << " Total energy deposit= "<< G4BestUnit(fEdeptot,"Energy")
143  << G4endl;
144  fEniel /= nbofEvents;
145  G4cout << " niel energy deposit = "<< G4BestUnit(fEniel,"Energy")
146  << G4endl;
147 
148  //print dose in tallies
149  //
150  G4int tallyNumber = fDetector->GetTallyNumber();
151  if (tallyNumber > 0) {
152  G4double Ebeam = fKinematic->GetEbeamCumul();
153  G4cout << "\n---------------------------------------------------------\n";
154  G4cout << " Cumulated Doses : \tEdep \tEdep/Ebeam \tDose" << G4endl;
155  for (G4int j=0; j < tallyNumber; ++j) {
156  G4double Edep = fTallyEdep[j], ratio = 100*Edep/Ebeam;
157  G4double tallyMass = fDetector->GetTallyMass(j);
158  G4double Dose = Edep/tallyMass;
159  G4cout << " tally " << j << ": \t \t"
160  << G4BestUnit(Edep,"Energy") << "\t"
161  << ratio << " % \t"
162  << G4BestUnit(Dose,"Dose") << G4endl;
163  }
164  G4cout << "\n---------------------------------------------------------\n";
165  G4cout << G4endl;
166  }
167 
168  if (fAnalysisManager->IsActive() ) {
169  // normalize histograms
170  //
171  for (G4int j=1; j<3; ++j) {
172  G4double binWidth = fAnalysisManager->GetH1Width(j);
173  G4double fac = (mm/MeV)/(nbofEvents * binWidth);
174  fAnalysisManager->ScaleH1(j, fac);
175  }
176  fAnalysisManager->ScaleH1(3, 1./nbofEvents);
177 
178  // save histograms
179  fAnalysisManager->Write();
180  fAnalysisManager->CloseFile();
181  delete fAnalysisManager;
182  fAnalysisManager = 0;
183  }
184 
185  // show Rndm status
186  //
188 }
189 
190 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
191 
193 {
194  // Create or get analysis manager
195  // The choice of analysis technology is done via selection of a namespace
196  // in HistoManager.hh
197  fAnalysisManager = G4AnalysisManager::Instance();
198  fAnalysisManager->SetFileName("testem7");
199  fAnalysisManager->SetVerboseLevel(1);
200  fAnalysisManager->SetActivation(true); // enable inactivation of histograms
201 
202  // Define histograms start values
203  const G4int kMaxHisto = 4;
204  const G4String id[] = { "h0", "h1", "h2", "h3" };
205  const G4String title[] =
206  { "dummy", //0
207  "Edep (MeV/mm) along absorber ", //1
208  "Edep (MeV/mm) along absorber zoomed", //2
209  "projectile range" //3
210  };
211 
212  // Default values (to be reset via /analysis/h1/set command)
213  G4int nbins = 100;
214  G4double vmin = 0.;
215  G4double vmax = 100.;
216 
217  // Create all histograms as inactivated
218  // as we have not yet set nbins, vmin, vmax
219  for (G4int k=0; k<kMaxHisto; ++k) {
220  G4int ih = fAnalysisManager->CreateH1(id[k], title[k], nbins, vmin, vmax);
221  G4bool activ = false;
222  if (k == 1) activ = true;
223  fAnalysisManager->SetH1Activation(ih, activ);
224  }
225 }
226 
227 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
static constexpr double MeV
Definition: G4SIunits.hh:214
static constexpr double mm
Definition: G4SIunits.hh:115
#define G4endl
Definition: G4ios.hh:61
const G4String & GetParticleName() const
G4int GetNumberOfEvent() const
Definition: G4Run.hh:79
const G4String & GetName() const
Definition: G4Material.hh:179
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
G4int GetRunID() const
Definition: G4Run.hh:76
double energy
Definition: plottest35.C:25
The primary generator action class with particle gun.
Definition: G4Run.hh:46
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
int G4int
Definition: G4Types.hh:78
G4double GetParticleEnergy() const
static const G4double fac
G4GLOB_DLL std::ostream G4cout
G4ParticleDefinition * GetParticleDefinition() const
static void showEngineStatus()
Definition: Random.cc:303
Simple detector construction with a box volume placed in a world.
G4double GetDensity() const
Definition: G4Material.hh:181