Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/examples/extended/medical/fanoCavity2/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 medical/fanoCavity2/src/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 "G4Electron.hh"
 40 #include "G4EmCalculator.hh"
 41 #include "G4Run.hh"
 42 #include "G4RunManager.hh"
 43 #include "G4SystemOfUnits.hh"
 44 #include "G4UnitsTable.hh"
 45 #include "Randomize.hh"
 46 
 47 #include <iomanip>
 48 
 49 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 50 
 51 Run::Run(DetectorConstruction* det, PrimaryGeneratorAction* kin, bool isMaster)
 52   : G4Run(),
 53     fDetector(det),
 54     fKinematic(kin),
 55     fProcCounter(0),
 56     fEdepCavity(0.),
 57     fEdepCavity2(0.),
 58     fTrkSegmCavity(0.),
 59     fNbEventCavity(0),
 60     fStepWall(0.),
 61     fStepWall2(0.),
 62     fStepCavity(0.),
 63     fStepCavity2(0.),
 64     fNbStepWall(0),
 65     fNbStepCavity(0),
 66     fEnergyGun(0.),
 67     fMassWall(0.),
 68     fMassCavity(0.),
 69     fIsMaster(isMaster)
 70 {
 71   // run conditions
 72   //
 73   G4ParticleDefinition* particleGun = fKinematic->GetParticleGun()->GetParticleDefinition();
 74   G4String partName = particleGun->GetParticleName();
 75   fEnergyGun = fKinematic->GetParticleGun()->GetParticleEnergy();
 76 
 77   // geometry : effective wall volume
 78   //
 79   G4double cavityThickness = fDetector->GetCavityThickness();
 80   G4Material* mateCavity = fDetector->GetCavityMaterial();
 81   G4double densityCavity = mateCavity->GetDensity();
 82   fMassCavity = cavityThickness * densityCavity;
 83 
 84   G4double wallThickness = fDetector->GetWallThickness();
 85   G4Material* mateWall = fDetector->GetWallMaterial();
 86   G4double densityWall = mateWall->GetDensity();
 87 
 88   G4EmCalculator emCal;
 89   G4double RangeWall = emCal.GetCSDARange(fEnergyGun, particleGun, mateWall);
 90   G4double factor = 1.2;
 91   G4double effWallThick = factor * RangeWall;
 92   if ((effWallThick > wallThickness) || (effWallThick <= 0.)) effWallThick = wallThickness;
 93   fMassWall = 2 * effWallThick * densityWall;
 94 
 95   G4double massTotal = fMassWall + fMassCavity;
 96   G4double fMassWallRatio = fMassWall / massTotal;
 97   fKinematic->RunInitialisation(effWallThick, fMassWallRatio);
 98 
 99   G4double massRatio = fMassCavity / fMassWall;
100 
101   // check radius
102   //
103   G4double worldRadius = fDetector->GetWallRadius();
104   G4double RangeCavity = emCal.GetCSDARange(fEnergyGun, particleGun, mateCavity);
105 
106   // G4String partName    = particleGun->GetParticleName();
107 
108   std::ios::fmtflags mode = G4cout.flags();
109   G4cout.setf(std::ios::fixed, std::ios::floatfield);
110   G4int prec = G4cout.precision(3);
111 
112   G4cout << "\n ===================== run conditions =====================\n";
113 
114   G4cout << "\n The run will be " << numberOfEvent << " " << partName << " of "
115          << G4BestUnit(fEnergyGun, "Energy") << " through 2*" << G4BestUnit(effWallThick, "Length")
116          << " of " << mateWall->GetName()
117          << " (density: " << G4BestUnit(densityWall, "Volumic Mass")
118          << "); Mass/cm2 = " << G4BestUnit(fMassWall * cm2, "Mass")
119          << "\n csdaRange: " << G4BestUnit(RangeWall, "Length") << G4endl;
120 
121   G4cout << "\n the cavity is " << G4BestUnit(cavityThickness, "Length") << " of "
122          << mateCavity->GetName() << " (density: " << G4BestUnit(densityCavity, "Volumic Mass")
123          << "); Mass/cm2 = " << G4BestUnit(fMassCavity * cm2, "Mass")
124          << " --> massRatio = " << std::setprecision(6) << massRatio << G4endl;
125 
126   G4cout.precision(3);
127   G4cout << " Wall radius: " << G4BestUnit(worldRadius, "Length")
128          << "; range in cavity: " << G4BestUnit(RangeCavity, "Length") << G4endl;
129 
130   G4cout << "\n ==========================================================\n";
131 
132   // stopping power from EmCalculator
133   //
134   G4double dedxWall = emCal.GetDEDX(fEnergyGun, G4Electron::Electron(), mateWall);
135   dedxWall /= densityWall;
136   G4double dedxCavity = emCal.GetDEDX(fEnergyGun, G4Electron::Electron(), mateCavity);
137   dedxCavity /= densityCavity;
138 
139   G4cout << std::setprecision(4)
140          << "\n StoppingPower in wall   = " << G4BestUnit(dedxWall, "Energy*Surface/Mass")
141          << "\n               in cavity = " << G4BestUnit(dedxCavity, "Energy*Surface/Mass")
142          << G4endl;
143 
144   // process counter
145   //
146   fProcCounter = new ProcessesCount;
147 
148   // charged particles and energy flow in cavity
149   //
150   fPartFlowCavity[0] = fPartFlowCavity[1] = 0;
151   fEnerFlowCavity[0] = fEnerFlowCavity[1] = 0.;
152 
153   // total energy deposit and charged track segment in cavity
154   //
155   fEdepCavity = fEdepCavity2 = fTrkSegmCavity = 0.;
156   fNbEventCavity = 0;
157 
158   // stepLenth of charged particles
159   //
160   fStepWall = fStepWall2 = fStepCavity = fStepCavity2 = 0.;
161   fNbStepWall = fNbStepCavity = 0;
162 
163   // reset default formats
164   G4cout.setf(mode, std::ios::floatfield);
165   G4cout.precision(prec);
166 
167   // histograms
168   //
169   G4AnalysisManager* analysisManager = G4AnalysisManager::Instance();
170   if (analysisManager->IsActive()) {
171     analysisManager->OpenFile();
172   }
173 }
174 
175 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
176 
177 Run::~Run() {}
178 
179 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
180 
181 void Run::CountProcesses(G4String procName)
182 {
183   // does the process  already encounted ?
184   size_t nbProc = fProcCounter->size();
185   size_t i = 0;
186   while ((i < nbProc) && ((*fProcCounter)[i]->GetName() != procName))
187     i++;
188   if (i == nbProc) fProcCounter->push_back(new OneProcessCount(procName));
189 
190   (*fProcCounter)[i]->Count();
191 }
192 
193 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
194 
195 void Run::SurveyConvergence(G4int NbofEvents)
196 {
197   if (NbofEvents == 0) return;
198 
199   // beam fluence
200   //
201   G4int Nwall = fKinematic->GetWallCount();
202   G4int Ncavity = fKinematic->GetCavityCount();
203   G4double Iwall = Nwall / fMassWall;
204   G4double Icavity = Ncavity / fMassCavity;
205   G4double Iratio = Icavity / Iwall;
206   G4double Itot = NbofEvents / (fMassWall + fMassCavity);
207   G4double energyFluence = fEnergyGun * Itot;
208 
209   // total dose in cavity
210   //
211   G4double doseCavity = fEdepCavity / fMassCavity;
212   G4double ratio = doseCavity / energyFluence;
213   G4double err = 100 * (ratio - 1.);
214 
215   std::ios::fmtflags mode = G4cout.flags();
216   G4cout.setf(std::ios::fixed, std::ios::floatfield);
217   G4int prec = G4cout.precision(5);
218 
219   G4cout << "--->evntNb= " << NbofEvents << " Nwall= " << Nwall << " Ncav= " << Ncavity
220          << " Ic/Iw= " << Iratio << " Ne-_cav= " << fPartFlowCavity[0]
221          << " doseCavity/Ebeam= " << ratio << "  (100*(ratio-1) = " << err << " %) \n"
222          << G4endl;
223 
224   // reset default formats
225   G4cout.setf(mode, std::ios::floatfield);
226   G4cout.precision(prec);
227 }
228 
229 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
230 
231 void Run::EndOfRun()
232 {  // Only call by Master thread
233 
234   std::ios::fmtflags mode = G4cout.flags();
235   G4cout.setf(std::ios::fixed, std::ios::floatfield);
236   G4int prec = G4cout.precision(3);
237 
238   if (numberOfEvent == 0) return;
239 
240   // frequency of processes
241   //
242   G4cout << "\n Process calls frequency --->";
243   for (size_t i = 0; i < fProcCounter->size(); i++) {
244     G4String procName = (*fProcCounter)[i]->GetName();
245     G4int count = (*fProcCounter)[i]->GetCounter();
246     G4cout << "  " << procName << "= " << count;
247   }
248   G4cout << G4endl;
249 
250   // charged particle flow in cavity
251   //
252   G4cout << "\n Charged particle flow in cavity :"
253          << "\n      Enter --> nbParticles = " << fPartFlowCavity[0]
254          << "\t Energy = " << G4BestUnit(fEnerFlowCavity[0], "Energy")
255          << "\n      Exit  --> nbParticles = " << fPartFlowCavity[1]
256          << "\t Energy = " << G4BestUnit(fEnerFlowCavity[1], "Energy") << G4endl;
257 
258   if (fPartFlowCavity[0] == 0) return;
259 
260   G4int Nwall = fKinematic->GetWallCount();
261   G4int Ncavity = fKinematic->GetCavityCount();
262 
263   G4double Iwall = Nwall / fMassWall;
264   G4double Icavity = Ncavity / fMassCavity;
265   G4double Iratio = Icavity / Iwall;
266   G4double Itot = numberOfEvent / (fMassWall + fMassCavity);
267   G4double energyFluence = fEnergyGun * Itot;
268 
269   G4cout.precision(5);
270   G4cout << "\n beamFluence in wall = " << Nwall << "\t in cavity = " << Ncavity
271          << "\t Icav/Iwall = " << Iratio
272          << "\t energyFluence = " << energyFluence / (MeV * cm2 / mg) << " MeV*cm2/mg" << G4endl;
273 
274   // error on Edep in cavity
275   //
276   if (fNbEventCavity == 0) return;
277   G4double meanEdep = fEdepCavity / fNbEventCavity;
278   G4double meanEdep2 = fEdepCavity2 / fNbEventCavity;
279   G4double varianceEdep = meanEdep2 - meanEdep * meanEdep;
280   G4double dEoverE = 0.;
281   if (varianceEdep > 0.) dEoverE = std::sqrt(varianceEdep / fNbEventCavity) / meanEdep;
282 
283   // total dose in cavity
284   //
285   G4double doseCavity = fEdepCavity / fMassCavity;
286 
287   G4double ratio = doseCavity / energyFluence, error = ratio * dEoverE;
288 
289   G4cout << "\n Total edep in cavity = " << G4BestUnit(fEdepCavity, "Energy") << " +- "
290          << 100 * dEoverE << " %"
291          << "\n Total dose in cavity = " << doseCavity / (MeV * cm2 / mg) << " MeV*cm2/mg"
292          << " +- " << 100 * dEoverE << " %"
293          << "\n\n DoseCavity/EnergyFluence = " << ratio << " +- " << error << G4endl;
294 
295   // track length in cavity
296   G4double meantrack = fTrkSegmCavity / fPartFlowCavity[0];
297 
298   G4cout.precision(4);
299   G4cout << "\n Total charged trackLength in cavity = " << G4BestUnit(fTrkSegmCavity, "Length")
300          << "   (mean value = " << G4BestUnit(meantrack, "Length") << ")" << G4endl;
301 
302   // compute mean step size of charged particles
303   //
304   fStepWall /= fNbStepWall;
305   fStepWall2 /= fNbStepWall;
306   G4double rms = fStepWall2 - fStepWall * fStepWall;
307   if (rms > 0.)
308     rms = std::sqrt(rms);
309   else
310     rms = 0.;
311   G4double nbTrackWall = fKinematic->GetWallCount();
312 
313   G4cout << "\n StepSize of ch. tracks in wall   = " << G4BestUnit(fStepWall, "Length") << " +- "
314          << G4BestUnit(rms, "Length") << "\t (nbSteps/track = " << double(fNbStepWall) / nbTrackWall
315          << ")";
316 
317   fStepCavity /= fNbStepCavity;
318   fStepCavity2 /= fNbStepCavity;
319   rms = fStepCavity2 - fStepCavity * fStepCavity;
320   if (rms > 0.)
321     rms = std::sqrt(rms);
322   else
323     rms = 0.;
324 
325   G4cout << "\n StepSize of ch. tracks in cavity = " << G4BestUnit(fStepCavity, "Length") << " +- "
326          << G4BestUnit(rms, "Length")
327          << "\t (nbSteps/track = " << double(fNbStepCavity) / fPartFlowCavity[0] << ")";
328 
329   G4cout << G4endl;
330 
331   // reset default formats
332   G4cout.setf(mode, std::ios::floatfield);
333   G4cout.precision(prec);
334 
335   // delete and remove all contents in fProcCounter
336   while (fProcCounter->size() > 0) {
337     OneProcessCount* aProcCount = fProcCounter->back();
338     fProcCounter->pop_back();
339     delete aProcCount;
340   }
341   delete fProcCounter;
342 
343   // show Rndm status
344   CLHEP::HepRandom::showEngineStatus();
345 }
346 
347 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
348 
349 void Run::Merge(const G4Run* run)
350 {
351   const Run* localRun = static_cast<const Run*>(run);
352 
353   // Merge Run variables
354   fPartFlowCavity[0] += localRun->fPartFlowCavity[0];
355   fPartFlowCavity[1] += localRun->fPartFlowCavity[1];
356   fEnerFlowCavity[0] += localRun->fEnerFlowCavity[0];
357   fEnerFlowCavity[1] += localRun->fEnerFlowCavity[1];
358   fEdepCavity += localRun->fEdepCavity;
359   fEdepCavity2 += localRun->fEdepCavity2;
360   fTrkSegmCavity += localRun->fTrkSegmCavity;
361   fNbEventCavity += localRun->fNbEventCavity;
362   fStepWall += localRun->fStepWall;
363   fStepWall2 += localRun->fStepWall2;
364   fStepCavity += localRun->fStepCavity;
365   fStepCavity2 += localRun->fStepCavity2;
366   fNbStepWall += localRun->fNbStepWall;
367   fNbStepCavity += localRun->fNbStepCavity;
368 
369   // Merge PrimaryGenerator variables
370   fKinematic->AddWallCount(localRun->fKinematic->GetWallCount());
371   fKinematic->AddCavityCount(localRun->fKinematic->GetCavityCount());
372 
373   // Merge ProcessCount varaibles
374   std::vector<OneProcessCount*>::iterator it;
375   for (it = localRun->fProcCounter->begin(); it != localRun->fProcCounter->end(); it++) {
376     OneProcessCount* process = *it;
377     for (G4int i = 0; i < process->GetCounter(); i++)
378       this->CountProcesses(process->GetName());
379   }
380 
381   G4Run::Merge(run);
382 }
383 
384 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
385