Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/examples/extended/hadronic/Hadr06/src/Run.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 Run.cc
 27 /// \brief Implementation of the Run class
 28 //
 29 //
 30 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 31 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 32 
 33 #include "Run.hh"
 34 
 35 #include "DetectorConstruction.hh"
 36 #include "HistoManager.hh"
 37 #include "PrimaryGeneratorAction.hh"
 38 
 39 #include "G4SystemOfUnits.hh"
 40 #include "G4UnitsTable.hh"
 41 
 42 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 43 
 44 Run::Run(DetectorConstruction* det) : fDetector(det) {}
 45 
 46 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 47 
 48 void Run::SetPrimary(G4ParticleDefinition* particle, G4double energy)
 49 {
 50   fParticle = particle;
 51   fEkin = energy;
 52 }
 53 
 54 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 55 
 56 void Run::CountProcesses(const G4VProcess* process)
 57 {
 58   if (process == nullptr) return;
 59   G4String procName = process->GetProcessName();
 60   std::map<G4String, G4int>::iterator it = fProcCounter.find(procName);
 61   if (it == fProcCounter.end()) {
 62     fProcCounter[procName] = 1;
 63   }
 64   else {
 65     fProcCounter[procName]++;
 66   }
 67 }
 68 
 69 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 70 
 71 void Run::ParticleCount(G4String name, G4double Ekin, G4double meanLife)
 72 {
 73   std::map<G4String, ParticleData>::iterator it = fParticleDataMap1.find(name);
 74   if (it == fParticleDataMap1.end()) {
 75     fParticleDataMap1[name] = ParticleData(1, Ekin, Ekin, Ekin, meanLife);
 76   }
 77   else {
 78     ParticleData& data = it->second;
 79     data.fCount++;
 80     data.fEmean += Ekin;
 81     // update min max
 82     G4double emin = data.fEmin;
 83     if (Ekin < emin) data.fEmin = Ekin;
 84     G4double emax = data.fEmax;
 85     if (Ekin > emax) data.fEmax = Ekin;
 86     data.fTmean = meanLife;
 87   }
 88 }
 89 
 90 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 91 
 92 void Run::SumEnergies(G4double edep, G4double eflow, G4double etot)
 93 {
 94   fEnergyDeposit += edep;
 95   fEnergyDeposit2 += edep * edep;
 96 
 97   fEnergyFlow += eflow;
 98   fEnergyFlow2 += eflow * eflow;
 99 
100   fEnergyTotal += etot;
101   fEnergyTotal2 += etot * etot;
102 }
103 
104 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
105 
106 void Run::ParticleFlux(G4String name, G4double Ekin)
107 {
108   std::map<G4String, ParticleData>::iterator it = fParticleDataMap2.find(name);
109   if (it == fParticleDataMap2.end()) {
110     fParticleDataMap2[name] = ParticleData(1, Ekin, Ekin, Ekin, -1 * ns);
111   }
112   else {
113     ParticleData& data = it->second;
114     data.fCount++;
115     data.fEmean += Ekin;
116     // update min max
117     G4double emin = data.fEmin;
118     if (Ekin < emin) data.fEmin = Ekin;
119     G4double emax = data.fEmax;
120     if (Ekin > emax) data.fEmax = Ekin;
121     data.fTmean = -1 * ns;
122   }
123 }
124 
125 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
126 
127 void Run::Merge(const G4Run* run)
128 {
129   const Run* localRun = static_cast<const Run*>(run);
130 
131   // primary particle info
132   //
133   fParticle = localRun->fParticle;
134   fEkin = localRun->fEkin;
135 
136   // accumulate sums
137   //
138   fEnergyDeposit += localRun->fEnergyDeposit;
139   fEnergyDeposit2 += localRun->fEnergyDeposit2;
140   fEnergyFlow += localRun->fEnergyFlow;
141   fEnergyFlow2 += localRun->fEnergyFlow2;
142   fEnergyTotal += localRun->fEnergyTotal;
143   fEnergyTotal2 += localRun->fEnergyTotal2;
144 
145   // map: processes count
146   std::map<G4String, G4int>::const_iterator itp;
147   for (itp = localRun->fProcCounter.begin(); itp != localRun->fProcCounter.end(); ++itp) {
148     G4String procName = itp->first;
149     G4int localCount = itp->second;
150     if (fProcCounter.find(procName) == fProcCounter.end()) {
151       fProcCounter[procName] = localCount;
152     }
153     else {
154       fProcCounter[procName] += localCount;
155     }
156   }
157 
158   // map: created particles count
159   std::map<G4String, ParticleData>::const_iterator itc;
160   for (itc = localRun->fParticleDataMap1.begin(); itc != localRun->fParticleDataMap1.end(); ++itc) {
161     G4String name = itc->first;
162     const ParticleData& localData = itc->second;
163     if (fParticleDataMap1.find(name) == fParticleDataMap1.end()) {
164       fParticleDataMap1[name] = ParticleData(localData.fCount, localData.fEmean, localData.fEmin,
165                                              localData.fEmax, localData.fTmean);
166     }
167     else {
168       ParticleData& data = fParticleDataMap1[name];
169       data.fCount += localData.fCount;
170       data.fEmean += localData.fEmean;
171       G4double emin = localData.fEmin;
172       if (emin < data.fEmin) data.fEmin = emin;
173       G4double emax = localData.fEmax;
174       if (emax > data.fEmax) data.fEmax = emax;
175       data.fTmean = localData.fTmean;
176     }
177   }
178 
179   // map: particles flux count
180   std::map<G4String, ParticleData>::const_iterator itn;
181   for (itn = localRun->fParticleDataMap2.begin(); itn != localRun->fParticleDataMap2.end(); ++itn) {
182     G4String name = itn->first;
183     const ParticleData& localData = itn->second;
184     if (fParticleDataMap2.find(name) == fParticleDataMap2.end()) {
185       fParticleDataMap2[name] = ParticleData(localData.fCount, localData.fEmean, localData.fEmin,
186                                              localData.fEmax, localData.fTmean);
187     }
188     else {
189       ParticleData& data = fParticleDataMap2[name];
190       data.fCount += localData.fCount;
191       data.fEmean += localData.fEmean;
192       G4double emin = localData.fEmin;
193       if (emin < data.fEmin) data.fEmin = emin;
194       G4double emax = localData.fEmax;
195       if (emax > data.fEmax) data.fEmax = emax;
196       data.fTmean = localData.fTmean;
197     }
198   }
199 
200   G4Run::Merge(run);
201 }
202 
203 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
204 
205 void Run::EndOfRun()
206 {
207   G4int prec = 5, wid = prec + 2;
208   G4int dfprec = G4cout.precision(prec);
209 
210   // run condition
211   //
212   G4Material* material = fDetector->GetMaterial();
213   G4double density = material->GetDensity();
214 
215   G4String Particle = fParticle->GetParticleName();
216   G4cout << "\n The run is " << numberOfEvent << " " << Particle << " of "
217          << G4BestUnit(fEkin, "Energy") << " through "
218          << G4BestUnit(fDetector->GetRadius(), "Length") << " of " << material->GetName()
219          << " (density: " << G4BestUnit(density, "Volumic Mass") << ")" << G4endl;
220 
221   if (numberOfEvent == 0) {
222     G4cout.precision(dfprec);
223     return;
224   }
225 
226   // frequency of processes
227   //
228   G4cout << "\n Process calls frequency :" << G4endl;
229   G4int index = 0;
230   std::map<G4String, G4int>::iterator it;
231   for (it = fProcCounter.begin(); it != fProcCounter.end(); it++) {
232     G4String procName = it->first;
233     G4int count = it->second;
234     G4String space = " ";
235     if (++index % 3 == 0) space = "\n";
236     G4cout << " " << std::setw(20) << procName << "=" << std::setw(7) << count << space;
237   }
238   G4cout << G4endl;
239 
240   // compute mean Energy deposited and rms
241   //
242   G4int TotNbofEvents = numberOfEvent;
243   fEnergyDeposit /= TotNbofEvents;
244   fEnergyDeposit2 /= TotNbofEvents;
245   G4double rmsEdep = fEnergyDeposit2 - fEnergyDeposit * fEnergyDeposit;
246   if (rmsEdep > 0.)
247     rmsEdep = std::sqrt(rmsEdep);
248   else
249     rmsEdep = 0.;
250 
251   G4cout << "\n Mean energy deposit per event = " << G4BestUnit(fEnergyDeposit, "Energy")
252          << ";  rms = " << G4BestUnit(rmsEdep, "Energy") << G4endl;
253 
254   // compute mean Energy leakage and rms
255   //
256   fEnergyFlow /= TotNbofEvents;
257   fEnergyFlow2 /= TotNbofEvents;
258   G4double rmsEflow = fEnergyFlow2 - fEnergyFlow * fEnergyFlow;
259   if (rmsEflow > 0.)
260     rmsEflow = std::sqrt(rmsEflow);
261   else
262     rmsEflow = 0.;
263 
264   G4cout << " Mean energy leakage per event = " << G4BestUnit(fEnergyFlow, "Energy")
265          << ";  rms = " << G4BestUnit(rmsEflow, "Energy") << G4endl;
266 
267   // energy balance
268   //
269   fEnergyTotal /= TotNbofEvents;
270   fEnergyTotal2 /= TotNbofEvents;
271   G4double rmsEtotal = fEnergyTotal2 - fEnergyTotal * fEnergyTotal;
272   if (rmsEtotal > 0.)
273     rmsEtotal = std::sqrt(rmsEtotal);
274   else
275     rmsEflow = 0.;
276 
277   G4cout << "\n Mean energy total   per event = " << G4BestUnit(fEnergyTotal, "Energy")
278          << ";  rms = " << G4BestUnit(rmsEtotal, "Energy") << G4endl;
279 
280   // particles at creation
281   //
282   if (fParticleDataMap1.size() > 0) {
283     G4cout << "\n List of particles at creation :" << G4endl;
284     std::map<G4String, ParticleData>::iterator itc;
285     for (itc = fParticleDataMap1.begin(); itc != fParticleDataMap1.end(); itc++) {
286       G4String name = itc->first;
287       ParticleData data = itc->second;
288       G4int count = data.fCount;
289       G4double eMean = data.fEmean / count;
290       G4double eMin = data.fEmin;
291       G4double eMax = data.fEmax;
292       G4double meanLife = data.fTmean;
293 
294       G4cout << "  " << std::setw(13) << name << ": " << std::setw(7) << count
295              << "  Emean = " << std::setw(wid) << G4BestUnit(eMean, "Energy") << "\t( "
296              << G4BestUnit(eMin, "Energy") << " --> " << G4BestUnit(eMax, "Energy") << ")";
297       if (meanLife >= 0.)
298         G4cout << "\tmean life = " << G4BestUnit(meanLife, "Time") << G4endl;
299       else
300         G4cout << "\tstable" << G4endl;
301     }
302   }
303 
304   // emerging particles
305   //
306   G4cout << "\n List of particles emerging from the absorber :" << G4endl;
307 
308   std::map<G4String, ParticleData>::iterator itn;
309   for (itn = fParticleDataMap2.begin(); itn != fParticleDataMap2.end(); itn++) {
310     G4String name = itn->first;
311     ParticleData data = itn->second;
312     G4int count = data.fCount;
313     G4double eMean = data.fEmean / count;
314     G4double eMin = data.fEmin;
315     G4double eMax = data.fEmax;
316     G4double Eflow = data.fEmean / TotNbofEvents;
317 
318     G4cout << "  " << std::setw(13) << name << ": " << std::setw(7) << count
319            << "  Emean = " << std::setw(wid) << G4BestUnit(eMean, "Energy") << "\t( "
320            << G4BestUnit(eMin, "Energy") << " --> " << G4BestUnit(eMax, "Energy")
321            << ") \tEleak/event = " << G4BestUnit(Eflow, "Energy") << G4endl;
322   }
323 
324   // normalize histograms
325   G4AnalysisManager* analysisManager = G4AnalysisManager::Instance();
326 
327   G4int ih = 2;
328   G4double binWidth = analysisManager->GetH1Width(ih);
329   G4double fac = (1. / (numberOfEvent * binWidth)) * (mm / MeV);
330   analysisManager->ScaleH1(ih, fac);
331 
332   // remove all contents in fProcCounter, fCount
333   fProcCounter.clear();
334   fParticleDataMap2.clear();
335 
336   // restore default format
337   G4cout.precision(dfprec);
338 }
339 
340 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
341