Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/examples/extended/hadronic/Hadr10/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 
 31 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 32 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 33 
 34 #include "Run.hh"
 35 
 36 #include "G4Run.hh"
 37 #include "G4RunManager.hh"
 38 #include "G4SystemOfUnits.hh"
 39 
 40 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 41 
 42 Run::Run()
 43   : G4Run(),
 44     fNumEvents(0),
 45     fPrimaryParticleId(0),
 46     fPrimaryParticleInitialKineticEnergy(0.0),
 47     fPrimaryParticleInitialTotalEnergy(0.0),
 48     fPrimaryParticleInitialMomentum(0.0),
 49     fPrimaryParticleInitialBeta(0.0),
 50     fPrimaryParticleInitialGamma(0.0),
 51     fPrimaryParticleInitial3Momentum(G4ThreeVector()),
 52     fPrimaryParticleInitialPosition(G4ThreeVector()),
 53     fToleranceEPviolations(0.0),
 54     fToleranceDeltaDecayRadius(0.0),
 55     fIsPreassignedDecayEnabled(true),
 56     fIsBoostToLabEnabled(true),
 57     fNumDecays(0),
 58     fNumBadDecays(0),
 59     fNumUnexpectedDecays(0),
 60     fNumEviolations(0),
 61     fNumPviolations(0),
 62     fNum_mc_truth_rPos_deltaMax_above(0),
 63     fNum_underestimated_mc_truth_rPos_delta_above(0),
 64     fNum_overestimated_mc_truth_rPos_delta_above(0),
 65     fNumLargeUnderestimates(0),
 66     fNumLargeOverestimates(0),
 67     fDecayT(0.0),
 68     fSumDecayT(0.0),
 69     fMinDecayT(999999.9),
 70     fMaxDecayT(-999999.9),
 71     fDecayR_mc_truth(0.0),
 72     fDecayR(0.0),
 73     fSumDecayR(0.0),
 74     fMinDecayR(999999.9),
 75     fMaxDecayR(-999999.9),
 76     fDecayX(0.0),
 77     fSumDecayX(0.0),
 78     fMinDecayX(999999.9),
 79     fMaxDecayX(-999999.9),
 80     fDecayY(0.0),
 81     fSumDecayY(0.0),
 82     fMinDecayY(999999.9),
 83     fMaxDecayY(-999999.9),
 84     fDecayZ(0.0),
 85     fSumDecayZ(0.0),
 86     fMinDecayZ(999999.9),
 87     fMaxDecayZ(-999999.9),
 88     fDeltaDecayR(0.0),
 89     fSumDeltaDecayR(0.0),
 90     fMinDeltaDecayR(999999.9),
 91     fMaxDeltaDecayR(-999999.9),
 92     fDeflectionAngle(0.0),
 93     fSumDeflectionAngle(0.0),
 94     fMinDeflectionAngle(999999.9),
 95     fMaxDeflectionAngle(-999999.9),
 96     fDeltaEkin(0.0),
 97     fSumDeltaEkin(0.0),
 98     fMinDeltaEkin(999999.9),
 99     fMaxDeltaEkin(-999999.9),
100     fDecayEkin(0.0),
101     fSumDecayEkin(0.0),
102     fMinDecayEkin(999999.9),
103     fMaxDecayEkin(-999999.9),
104     fDecayPx(0.0),
105     fSumDecayPx(0.0),
106     fMinDecayPx(999999.9),
107     fMaxDecayPx(-999999.9),
108     fDecayPy(0.0),
109     fSumDecayPy(0.0),
110     fMinDecayPy(999999.9),
111     fMaxDecayPy(-999999.9),
112     fDecayPz(0.0),
113     fSumDecayPz(0.0),
114     fMinDecayPz(999999.9),
115     fMaxDecayPz(-999999.9),
116     fDecayEtotViolation(0.0),
117     fSumDecayEtotViolation(0.0),
118     fMinDecayEtotViolation(999999.9),
119     fMaxDecayEtotViolation(-999999.9),
120     fDecayPxViolation(0.0),
121     fSumDecayPxViolation(0.0),
122     fMinDecayPxViolation(999999.9),
123     fMaxDecayPxViolation(-999999.9),
124     fDecayPyViolation(0.0),
125     fSumDecayPyViolation(0.0),
126     fMinDecayPyViolation(999999.9),
127     fMaxDecayPyViolation(-999999.9),
128     fDecayPzViolation(0.0),
129     fSumDecayPzViolation(0.0),
130     fMinDecayPzViolation(999999.9),
131     fMaxDecayPzViolation(-999999.9),
132     fMaxEkin_deltaMax(0.0),
133     fMaxEtot_deltaMax(0.0),
134     fMaxP_deltaMax(0.0),
135     fMaxPdir_deltaMax(0.0),
136     fMaxMass_deltaMax1(0.0),
137     fMaxMass_deltaMax2(0.0),
138     fSumMass_deltaMax3(0.0),
139     fMaxMass_deltaMax3(0.0),
140     fMaxBeta_deltaMax1(0.0),
141     fMaxBeta_deltaMax2(0.0),
142     fMaxGamma_deltaMax1(0.0),
143     fMaxGamma_deltaMax2(0.0),
144     fMaxGamma_deltaMax3(0.0),
145     fMaxT_proper_deltaMax(0.0),
146     fMaxT_lab_deltaMax(0.0),
147     fSumMc_truth_rPos_deltaMax(0.0),
148     fMaxMc_truth_rPos_deltaMax(0.0),
149     fSumUnderestimated_mc_truth_rPos_delta(0.0),
150     fMinUnderestimated_mc_truth_rPos_delta(999999.9),
151     fSumOverestimated_mc_truth_rPos_delta(0.0),
152     fMaxOverestimated_mc_truth_rPos_delta(-999999.9),
153     fSumUnderestimated_rDeltaPos(0.0),
154     fMinUnderestimated_rDeltaPos(999999.9),
155     fSumOverestimated_rDeltaPos(0.0),
156     fMaxOverestimated_rDeltaPos(-999999.9),
157     fMaxFloat_rDeltaPos_deltaMax(-999999.9)
158 {}
159 
160 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
161 
162 Run::~Run() {}
163 
164 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
165 
166 void Run::RecordEvent(const G4Event* anEvent)
167 {
168   // This method is called automatically by the Geant4 kernel (not by the user!)
169   // at the end of each event : in MT-mode, it is called only for the Working thread
170   // that handled that event.
171   G4int nEvt = anEvent->GetEventID();
172   if (nEvt % 1000 == 0) {
173     G4cout << std::setprecision(6) << " Event#=" << nEvt << " ; t[ns]=" << fDecayT
174            << " ; r[mm]=" << fDecayR
175            << " ; deltaR[mum]=" << (fDecayR_mc_truth - fDecayR) / CLHEP::micrometer
176            << " ; deltaEkin[MeV]=" << fDeltaEkin << " ; deltaAngle(deg)=" << fDeflectionAngle
177            << G4endl;
178   }
179   G4Run::RecordEvent(anEvent);
180 }
181 
182 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
183 
184 void Run::Merge(const G4Run* aRun)
185 {
186   // This method is called automatically by the Geant4 kernel (not by the user!)
187   // only in the case of multithreaded mode and only for Working threads.
188   const Run* localRun = static_cast<const Run*>(aRun);
189   fNumEvents += localRun->GetNumberOfEvent();
190 
191   fPrimaryParticleId = localRun->GetPrimaryParticleId();
192   fPrimaryParticleInitialKineticEnergy = localRun->GetPrimaryParticleInitialKineticEnergy();
193   fPrimaryParticleInitialTotalEnergy = localRun->GetPrimaryParticleInitialTotalEnergy();
194   fPrimaryParticleInitialMomentum = localRun->GetPrimaryParticleInitialMomentum();
195   fPrimaryParticleInitialBeta = localRun->GetPrimaryParticleInitialBeta();
196   fPrimaryParticleInitialGamma = localRun->GetPrimaryParticleInitialGamma();
197   fPrimaryParticleInitial3Momentum = localRun->GetPrimaryParticleInitial3Momentum();
198   fPrimaryParticleInitialPosition = localRun->GetPrimaryParticleInitialPosition();
199   fToleranceEPviolations = localRun->GetToleranceEPviolations();
200   fToleranceDeltaDecayRadius = localRun->GetToleranceDeltaDecayRadius();
201   fIsPreassignedDecayEnabled = localRun->GetIsPreassignedDecayEnabled();
202   fIsBoostToLabEnabled = localRun->GetIsBoostToLabEnabled();
203 
204   fNumDecays += localRun->GetNumberDecays();
205   fNumBadDecays += localRun->GetNumberBadDecays();
206   fNumUnexpectedDecays += localRun->GetNumberUnexpectedDecays();
207   fNumEviolations += localRun->GetNumberEviolations();
208   fNumPviolations += localRun->GetNumberPviolations();
209   fNum_mc_truth_rPos_deltaMax_above += localRun->GetNumber_mc_truth_rPos_deltaMax_above();
210   fNum_underestimated_mc_truth_rPos_delta_above +=
211     localRun->GetNumberUnderestimated_mc_truth_rPos_delta_above();
212   fNum_overestimated_mc_truth_rPos_delta_above +=
213     localRun->GetNumberOverestimated_mc_truth_rPos_delta_above();
214   fNumLargeUnderestimates += localRun->GetNumberLargeUnderestimates();
215   fNumLargeOverestimates += localRun->GetNumberLargeOverestimates();
216 
217   fSumDecayT += localRun->GetSumDecayT();
218   fMinDecayT = std::min(fMinDecayT, localRun->GetMinDecayT());
219   fMaxDecayT = std::max(fMaxDecayT, localRun->GetMaxDecayT());
220   fSumDecayR += localRun->GetSumDecayR();
221   fMinDecayR = std::min(fMinDecayR, localRun->GetMinDecayR());
222   fMaxDecayR = std::max(fMaxDecayR, localRun->GetMaxDecayR());
223   fSumDecayX += localRun->GetSumDecayX();
224   fMinDecayX = std::min(fMinDecayX, localRun->GetMinDecayX());
225   fMaxDecayX = std::max(fMaxDecayX, localRun->GetMaxDecayX());
226   fSumDecayY += localRun->GetSumDecayY();
227   fMinDecayY = std::min(fMinDecayY, localRun->GetMinDecayY());
228   fMaxDecayY = std::max(fMaxDecayY, localRun->GetMaxDecayY());
229   fSumDecayZ += localRun->GetSumDecayZ();
230   fMinDecayZ = std::min(fMinDecayZ, localRun->GetMinDecayZ());
231   fMaxDecayZ = std::max(fMaxDecayZ, localRun->GetMaxDecayZ());
232   fSumDeltaDecayR += localRun->GetSumDeltaDecayR();
233   fMinDeltaDecayR = std::min(fMinDeltaDecayR, localRun->GetMinDeltaDecayR());
234   fMaxDeltaDecayR = std::max(fMaxDeltaDecayR, localRun->GetMaxDeltaDecayR());
235   fSumDeflectionAngle += localRun->GetSumDeflectionAngle();
236   fMinDeflectionAngle = std::min(fMinDeflectionAngle, localRun->GetMinDeflectionAngle());
237   fMaxDeflectionAngle = std::max(fMaxDeflectionAngle, localRun->GetMaxDeflectionAngle());
238   fSumDeltaEkin += localRun->GetSumDeltaEkin();
239   fMinDeltaEkin = std::min(fMinDeltaEkin, localRun->GetMinDeltaEkin());
240   fMaxDeltaEkin = std::max(fMaxDeltaEkin, localRun->GetMaxDeltaEkin());
241   fSumDecayEkin += localRun->GetSumDecayEkin();
242   fMinDecayEkin = std::min(fMinDecayEkin, localRun->GetMinDecayEkin());
243   fMaxDecayEkin = std::max(fMaxDecayEkin, localRun->GetMaxDecayEkin());
244   fSumDecayPx += localRun->GetSumDecayPx();
245   fMinDecayPx = std::min(fMinDecayPx, localRun->GetMinDecayPx());
246   fMaxDecayPx = std::max(fMaxDecayPx, localRun->GetMaxDecayPx());
247   fSumDecayPy += localRun->GetSumDecayPy();
248   fMinDecayPy = std::min(fMinDecayPy, localRun->GetMinDecayPy());
249   fMaxDecayPy = std::max(fMaxDecayPy, localRun->GetMaxDecayPy());
250   fSumDecayPz += localRun->GetSumDecayPz();
251   fMinDecayPz = std::min(fMinDecayPz, localRun->GetMinDecayPz());
252   fMaxDecayPz = std::max(fMaxDecayPz, localRun->GetMaxDecayPz());
253   fSumDecayEtotViolation += localRun->GetSumDecayEtotViolation();
254   fMinDecayEtotViolation = std::min(fMinDecayEtotViolation, localRun->GetMinDecayEtotViolation());
255   fMaxDecayEtotViolation = std::max(fMaxDecayEtotViolation, localRun->GetMaxDecayEtotViolation());
256   fSumDecayPxViolation += localRun->GetSumDecayPxViolation();
257   fMinDecayPxViolation = std::min(fMinDecayPxViolation, localRun->GetMinDecayPxViolation());
258   fMaxDecayPxViolation = std::max(fMaxDecayPxViolation, localRun->GetMaxDecayPxViolation());
259   fSumDecayPyViolation += localRun->GetSumDecayPyViolation();
260   fMinDecayPyViolation = std::min(fMinDecayPyViolation, localRun->GetMinDecayPyViolation());
261   fMaxDecayPyViolation = std::max(fMaxDecayPyViolation, localRun->GetMaxDecayPyViolation());
262   fSumDecayPzViolation += localRun->GetSumDecayPzViolation();
263   fMinDecayPzViolation = std::min(fMinDecayPzViolation, localRun->GetMinDecayPzViolation());
264   fMaxDecayPzViolation = std::max(fMaxDecayPzViolation, localRun->GetMaxDecayPzViolation());
265 
266   fMaxEkin_deltaMax = std::max(fMaxEkin_deltaMax, localRun->GetMaxEkin_deltaMax());
267   fMaxEtot_deltaMax = std::max(fMaxEtot_deltaMax, localRun->GetMaxEtot_deltaMax());
268   fMaxP_deltaMax = std::max(fMaxP_deltaMax, localRun->GetMaxP_deltaMax());
269   fMaxPdir_deltaMax = std::max(fMaxPdir_deltaMax, localRun->GetMaxPdir_deltaMax());
270   fMaxMass_deltaMax1 = std::max(fMaxMass_deltaMax1, localRun->GetMaxMass_deltaMax1());
271   fMaxMass_deltaMax2 = std::max(fMaxMass_deltaMax2, localRun->GetMaxMass_deltaMax2());
272   fMaxMass_deltaMax3 = std::max(fMaxMass_deltaMax3, localRun->GetMaxMass_deltaMax3());
273   fSumMass_deltaMax3 += localRun->GetSumMass_deltaMax3();
274   fMaxBeta_deltaMax1 = std::max(fMaxBeta_deltaMax1, localRun->GetMaxBeta_deltaMax1());
275   fMaxBeta_deltaMax2 = std::max(fMaxBeta_deltaMax2, localRun->GetMaxBeta_deltaMax2());
276   fMaxGamma_deltaMax1 = std::max(fMaxGamma_deltaMax1, localRun->GetMaxGamma_deltaMax1());
277   fMaxGamma_deltaMax2 = std::max(fMaxGamma_deltaMax2, localRun->GetMaxGamma_deltaMax2());
278   fMaxGamma_deltaMax3 = std::max(fMaxGamma_deltaMax3, localRun->GetMaxGamma_deltaMax3());
279   fMaxT_proper_deltaMax = std::max(fMaxT_proper_deltaMax, localRun->GetMaxT_proper_deltaMax());
280   fMaxT_lab_deltaMax = std::max(fMaxT_lab_deltaMax, localRun->GetMaxT_lab_deltaMax());
281   fMaxMc_truth_rPos_deltaMax =
282     std::max(fMaxMc_truth_rPos_deltaMax, localRun->GetMaxMc_truth_rPos_deltaMax());
283   fSumMc_truth_rPos_deltaMax += localRun->GetSumMc_truth_rPos_deltaMax();
284 
285   fMinUnderestimated_mc_truth_rPos_delta = std::min(
286     fMinUnderestimated_mc_truth_rPos_delta, localRun->GetMinUnderestimated_mc_truth_rPos_delta());
287   fSumUnderestimated_mc_truth_rPos_delta += localRun->GetSumUnderestimated_mc_truth_rPos_delta();
288   fMaxOverestimated_mc_truth_rPos_delta = std::max(
289     fMaxOverestimated_mc_truth_rPos_delta, localRun->GetMaxOverestimated_mc_truth_rPos_delta());
290   fSumOverestimated_mc_truth_rPos_delta += localRun->GetSumOverestimated_mc_truth_rPos_delta();
291   fMinUnderestimated_rDeltaPos =
292     std::min(fMinUnderestimated_rDeltaPos, localRun->GetMinUnderestimated_rDeltaPos());
293   fSumUnderestimated_rDeltaPos += localRun->GetSumUnderestimated_rDeltaPos();
294   fMaxOverestimated_rDeltaPos =
295     std::max(fMaxOverestimated_rDeltaPos, localRun->GetMaxOverestimated_rDeltaPos());
296   fSumOverestimated_rDeltaPos += localRun->GetSumOverestimated_rDeltaPos();
297 
298   fMaxFloat_rDeltaPos_deltaMax =
299     std::max(fMaxFloat_rDeltaPos_deltaMax, localRun->GetMaxFloat_rDeltaPos_deltaMax());
300 
301   G4Run::Merge(aRun);
302 }
303 
304 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
305 
306 void Run::PrintInfo() const
307 {
308   // This method is called by RunAction::EndOfRunAction.
309   // In MT-mode, only the master thread calls it.
310   const G4double femtosecond = 0.001 * CLHEP::picosecond;
311   const G4double NN = (fNumDecays > 0 ? G4double(fNumDecays) : 1.0);
312   G4cout
313     << std::setprecision(6) << G4endl << G4endl
314     << " ============ Run::printInfo() =============== \t RunID = " << GetRunID() << G4endl
315     << " primary PDG code                     = " << fPrimaryParticleId << G4endl
316     << " primary initial kinetic energy [GeV] = "
317     << fPrimaryParticleInitialKineticEnergy / CLHEP::GeV << G4endl
318     << " primary initial total energy [GeV]   = " << fPrimaryParticleInitialTotalEnergy / CLHEP::GeV
319     << G4endl
320     << " primary initial momentum [GeV]       = " << fPrimaryParticleInitialMomentum / CLHEP::GeV
321     << G4endl << " primary initial Lorentz beta         = " << fPrimaryParticleInitialBeta << G4endl
322     << " primary initial Lorentz gamma        = " << fPrimaryParticleInitialGamma << G4endl
323     << " primary initial 3 momentum [GeV]     = " << fPrimaryParticleInitial3Momentum / CLHEP::GeV
324     << G4endl << " primary initial position [mm]        = " << fPrimaryParticleInitialPosition
325     << G4endl << " toleranceEPviolations [eV]           = " << fToleranceEPviolations / CLHEP::eV
326     << G4endl
327     << " toleranceDeltaDecayRadius [mum]      = " << fToleranceDeltaDecayRadius / CLHEP::micrometer
328     << G4endl << " isPreassignedDecayEnabled            = " << fIsPreassignedDecayEnabled << G4endl
329     << " isBoostToLabEnabled                  = " << fIsBoostToLabEnabled << G4endl
330     << " # Events               = " << (fNumEvents > 0 ? fNumEvents : GetNumberOfEvent()) << G4endl
331     << " # Decays               = " << fNumDecays << G4endl
332     << "   (# Bad decays        = " << fNumBadDecays << " )" << G4endl
333     << "   (# Unexpected decays = " << fNumUnexpectedDecays << " )" << G4endl
334     << " # E violations         = " << fNumEviolations << G4endl
335     << " # P violations         = " << fNumPviolations << G4endl
336     << " decay T [ns] : min=" << fMinDecayT << "\t mean=" << fSumDecayT / NN
337     << "\t max=" << fMaxDecayT << G4endl << " decay R [mm] : min=" << fMinDecayR
338     << "\t mean=" << fSumDecayR / NN << "\t max=" << fMaxDecayR << G4endl
339     << " decay X [mm] : min=" << fMinDecayX << "\t mean=" << fSumDecayX / NN
340     << "\t max=" << fMaxDecayX << G4endl << " decay Y [mm] : min=" << fMinDecayY
341     << "\t mean=" << fSumDecayY / NN << "\t max=" << fMaxDecayY << G4endl
342     << " decay Z [mm] : min=" << fMinDecayZ << "\t mean=" << fSumDecayZ / NN
343     << "\t max=" << fMaxDecayZ << G4endl << " Delta decay R [mm] : min=" << fMinDeltaDecayR
344     << "\t mean=" << fSumDeltaDecayR / NN << "\t max=" << fMaxDeltaDecayR << G4endl
345     << " deflection angle [deg] : min=" << fMinDeflectionAngle
346     << "\t mean=" << fSumDeflectionAngle / NN << "\t max=" << fMaxDeflectionAngle << G4endl
347     << " Delta Ekin [MeV] : min=" << fMinDeltaEkin << "\t mean=" << fSumDeltaEkin / NN
348     << "\t max=" << fMaxDeltaEkin << G4endl
349     << " decay Ekin [GeV] : min=" << fMinDecayEkin / CLHEP::GeV
350     << "\t mean=" << fSumDecayEkin / NN / CLHEP::GeV << "\t max=" << fMaxDecayEkin / CLHEP::GeV
351     << G4endl << " decay Px [GeV]   : min=" << fMinDecayPx / CLHEP::GeV
352     << "\t mean=" << fSumDecayPx / NN / CLHEP::GeV << "\t max=" << fMaxDecayPx / CLHEP::GeV
353     << G4endl << " decay Py [GeV]   : min=" << fMinDecayPy / CLHEP::GeV
354     << "\t mean=" << fSumDecayPy / NN / CLHEP::GeV << "\t max=" << fMaxDecayPy / CLHEP::GeV
355     << G4endl << " decay Pz [GeV]   : min=" << fMinDecayPz / CLHEP::GeV
356     << "\t mean=" << fSumDecayPz / NN / CLHEP::GeV << "\t max=" << fMaxDecayPz / CLHEP::GeV
357     << G4endl << " decay Etot violation [MeV] : min=" << fMinDecayEtotViolation
358     << "\t mean=" << fSumDecayEtotViolation / NN << "\t max=" << fMaxDecayEtotViolation << G4endl
359     << " decay Px violation [MeV]   : min=" << fMinDecayPxViolation
360     << "\t mean=" << fSumDecayPxViolation / NN << "\t max=" << fMaxDecayPxViolation << G4endl
361     << " decay Py violation [MeV]   : min=" << fMinDecayPyViolation
362     << "\t mean=" << fSumDecayPyViolation / NN << "\t max=" << fMaxDecayPyViolation << G4endl
363     << " decay Pz violation [MeV]   : min=" << fMinDecayPzViolation
364     << "\t mean=" << fSumDecayPzViolation / NN << "\t max=" << fMaxDecayPzViolation << G4endl
365     << " --- Consistency checks --- " << G4endl
366     << " maxEkin_deltaMax [eV]           = " << fMaxEkin_deltaMax / CLHEP::eV << G4endl
367     << " maxEtot_deltaMax [eV]           = " << fMaxEtot_deltaMax / CLHEP::eV << G4endl
368     << " maxP_deltaMax [eV]              = " << fMaxP_deltaMax / CLHEP::eV << G4endl
369     << " maxPdir_deltaMax                = " << fMaxPdir_deltaMax << G4endl
370     << " maxMass_deltaMax{1,2,3} [eV]    = " << fMaxMass_deltaMax1 / CLHEP::eV << " , "
371     << fMaxMass_deltaMax2 / CLHEP::eV << " , " << fMaxMass_deltaMax3 / CLHEP::eV
372     << " (mean=" << fSumMass_deltaMax3 / NN / CLHEP::eV << ")" << G4endl
373     << " maxBeta_deltaMax{1,2}           = " << fMaxBeta_deltaMax1 << " , " << fMaxBeta_deltaMax2
374     << G4endl << " maxGamma_deltaMax{1,2,3}        = " << fMaxGamma_deltaMax1 << " , "
375     << fMaxGamma_deltaMax2 << " , " << fMaxGamma_deltaMax3 << G4endl
376     << " maxT_proper_deltaMax [fs]       = " << fMaxT_proper_deltaMax / femtosecond << G4endl
377     << " maxT_lab_deltaMax [fs]          = " << fMaxT_lab_deltaMax / femtosecond << G4endl
378     << " maxMc_truth_rPos_deltaMax [mum] = " << fMaxMc_truth_rPos_deltaMax / CLHEP::micrometer
379     << " (mean=" << fSumMc_truth_rPos_deltaMax / NN / CLHEP::micrometer
380     << ")\t (# above threshold = " << fNum_mc_truth_rPos_deltaMax_above << ")" << G4endl
381     << " --- Extra checks --- " << G4endl << " minUnderestimated_mc_truth_rPos_delta [mum] = "
382     << fMinUnderestimated_mc_truth_rPos_delta / CLHEP::micrometer
383     << " (mean=" << fSumUnderestimated_mc_truth_rPos_delta / NN / CLHEP::micrometer
384     << ")\t (#above threshold = " << fNum_underestimated_mc_truth_rPos_delta_above << ")" << G4endl
385     << " maxOverestimated_mc_truth_rPos_delta [mum]  = "
386     << fMaxOverestimated_mc_truth_rPos_delta / CLHEP::micrometer
387     << " (mean=" << fSumOverestimated_mc_truth_rPos_delta / NN / CLHEP::micrometer
388     << ")\t (#above threshold = " << fNum_overestimated_mc_truth_rPos_delta_above << ")" << G4endl
389     << " minUnderestimated_rDeltaPos [mum]           = "
390     << fMinUnderestimated_rDeltaPos / CLHEP::micrometer
391     << " (mean=" << fSumUnderestimated_rDeltaPos / NN / CLHEP::micrometer
392     << ")\t (#above threshold = " << fNumLargeUnderestimates << ")" << G4endl
393     << " maxOverestimated_rDeltaPos [mum]            = "
394     << fMaxOverestimated_rDeltaPos / CLHEP::micrometer
395     << " (mean=" << fSumOverestimated_rDeltaPos / NN / CLHEP::micrometer
396     << ")\t (#above threshold = " << fNumLargeOverestimates << ")" << G4endl
397     << " --- float  instead of  double --- " << G4endl
398     << " fMaxFloat_rDeltaPos_deltaMax [mum] = " << fMaxFloat_rDeltaPos_deltaMax / CLHEP::micrometer
399     << G4endl << " ============================================= " << G4endl << G4endl;
400 }
401 
402 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
403 
404 void Run::SetDecayT(const G4double inputValue)
405 {
406   fDecayT = inputValue;
407   fSumDecayT += inputValue;
408   fMinDecayT = std::min(fMinDecayT, inputValue);
409   fMaxDecayT = std::max(fMaxDecayT, inputValue);
410 }
411 
412 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
413 
414 void Run::SetDecayR_mc_truth(const G4double inputValue)
415 {
416   fDecayR_mc_truth = inputValue;
417 }
418 
419 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
420 
421 void Run::SetDecayR(const G4double inputValue)
422 {
423   fDecayR = inputValue;
424   fSumDecayR += inputValue;
425   fMinDecayR = std::min(fMinDecayR, inputValue);
426   fMaxDecayR = std::max(fMaxDecayR, inputValue);
427 }
428 
429 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
430 
431 void Run::SetDecayX(const G4double inputValue)
432 {
433   fDecayX = inputValue;
434   fSumDecayX += inputValue;
435   fMinDecayX = std::min(fMinDecayX, inputValue);
436   fMaxDecayX = std::max(fMaxDecayX, inputValue);
437 }
438 
439 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
440 
441 void Run::SetDecayY(const G4double inputValue)
442 {
443   fDecayY = inputValue;
444   fSumDecayY += inputValue;
445   fMinDecayY = std::min(fMinDecayY, inputValue);
446   fMaxDecayY = std::max(fMaxDecayY, inputValue);
447 }
448 
449 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
450 
451 void Run::SetDecayZ(const G4double inputValue)
452 {
453   fDecayZ = inputValue;
454   fSumDecayZ += inputValue;
455   fMinDecayZ = std::min(fMinDecayZ, inputValue);
456   fMaxDecayZ = std::max(fMaxDecayZ, inputValue);
457 }
458 
459 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
460 
461 void Run::SetDeltaDecayR(const G4double inputValue)
462 {
463   fDeltaDecayR = inputValue;
464   fSumDeltaDecayR += inputValue;
465   fMinDeltaDecayR = std::min(fMinDeltaDecayR, inputValue);
466   fMaxDeltaDecayR = std::max(fMaxDeltaDecayR, inputValue);
467 }
468 
469 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
470 
471 void Run::SetDeflectionAngle(const G4double inputValue)
472 {
473   fDeflectionAngle = inputValue;
474   fSumDeflectionAngle += inputValue;
475   fMinDeflectionAngle = std::min(fMinDeflectionAngle, inputValue);
476   fMaxDeflectionAngle = std::max(fMaxDeflectionAngle, inputValue);
477 }
478 
479 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
480 
481 void Run::SetDeltaEkin(const G4double inputValue)
482 {
483   fDeltaEkin = inputValue;
484   fSumDeltaEkin += inputValue;
485   fMinDeltaEkin = std::min(fMinDeltaEkin, inputValue);
486   fMaxDeltaEkin = std::max(fMaxDeltaEkin, inputValue);
487 }
488 
489 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
490 
491 void Run::SetDecayEkin(const G4double inputValue)
492 {
493   fDecayEkin = inputValue;
494   fSumDecayEkin += inputValue;
495   fMinDecayEkin = std::min(fMinDecayEkin, inputValue);
496   fMaxDecayEkin = std::max(fMaxDecayEkin, inputValue);
497 }
498 
499 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
500 
501 void Run::SetDecayPx(const G4double inputValue)
502 {
503   fDecayPx = inputValue;
504   fSumDecayPx += inputValue;
505   fMinDecayPx = std::min(fMinDecayPx, inputValue);
506   fMaxDecayPx = std::max(fMaxDecayPx, inputValue);
507 }
508 
509 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
510 
511 void Run::SetDecayPy(const G4double inputValue)
512 {
513   fDecayPy = inputValue;
514   fSumDecayPy += inputValue;
515   fMinDecayPy = std::min(fMinDecayPy, inputValue);
516   fMaxDecayPy = std::max(fMaxDecayPy, inputValue);
517 }
518 
519 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
520 
521 void Run::SetDecayPz(const G4double inputValue)
522 {
523   fDecayPz = inputValue;
524   fSumDecayPz += inputValue;
525   fMinDecayPz = std::min(fMinDecayPz, inputValue);
526   fMaxDecayPz = std::max(fMaxDecayPz, inputValue);
527 }
528 
529 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
530 
531 void Run::SetDecayEtotViolation(const G4double inputValue)
532 {
533   fDecayEtotViolation = inputValue;
534   fSumDecayEtotViolation += inputValue;
535   fMinDecayEtotViolation = std::min(fMinDecayEtotViolation, inputValue);
536   fMaxDecayEtotViolation = std::max(fMaxDecayEtotViolation, inputValue);
537 }
538 
539 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
540 
541 void Run::SetDecayPxViolation(const G4double inputValue)
542 {
543   fDecayPxViolation = inputValue;
544   fSumDecayPxViolation += inputValue;
545   fMinDecayPxViolation = std::min(fMinDecayPxViolation, inputValue);
546   fMaxDecayPxViolation = std::max(fMaxDecayPxViolation, inputValue);
547 }
548 
549 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
550 
551 void Run::SetDecayPyViolation(const G4double inputValue)
552 {
553   fDecayPyViolation = inputValue;
554   fSumDecayPyViolation += inputValue;
555   fMinDecayPyViolation = std::min(fMinDecayPyViolation, inputValue);
556   fMaxDecayPyViolation = std::max(fMaxDecayPyViolation, inputValue);
557 }
558 
559 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
560 
561 void Run::SetDecayPzViolation(const G4double inputValue)
562 {
563   fDecayPzViolation = inputValue;
564   fSumDecayPzViolation += inputValue;
565   fMinDecayPzViolation = std::min(fMinDecayPzViolation, inputValue);
566   fMaxDecayPzViolation = std::max(fMaxDecayPzViolation, inputValue);
567 }
568 
569 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
570 
571 void Run::SetMaxEkin_deltaMax(const G4double inputValue)
572 {
573   fMaxEkin_deltaMax = std::max(fMaxEkin_deltaMax, inputValue);
574 }
575 
576 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
577 
578 void Run::SetMaxEtot_deltaMax(const G4double inputValue)
579 {
580   fMaxEtot_deltaMax = std::max(fMaxEtot_deltaMax, inputValue);
581 }
582 
583 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
584 
585 void Run::SetMaxP_deltaMax(const G4double inputValue)
586 {
587   fMaxP_deltaMax = std::max(fMaxP_deltaMax, inputValue);
588 }
589 
590 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
591 
592 void Run::SetMaxPdir_deltaMax(const G4double inputValue)
593 {
594   fMaxPdir_deltaMax = std::max(fMaxPdir_deltaMax, inputValue);
595 }
596 
597 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
598 
599 void Run::SetMaxMass_deltaMax1(const G4double inputValue)
600 {
601   fMaxMass_deltaMax1 = std::max(fMaxMass_deltaMax1, inputValue);
602 }
603 
604 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
605 
606 void Run::SetMaxMass_deltaMax2(const G4double inputValue)
607 {
608   fMaxMass_deltaMax2 = std::max(fMaxMass_deltaMax2, inputValue);
609 }
610 
611 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
612 
613 void Run::SetMaxMass_deltaMax3(const G4double inputValue)
614 {
615   fSumMass_deltaMax3 += inputValue;
616   fMaxMass_deltaMax3 = std::max(fMaxMass_deltaMax3, inputValue);
617 }
618 
619 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
620 
621 void Run::SetMaxBeta_deltaMax1(const G4double inputValue)
622 {
623   fMaxBeta_deltaMax1 = std::max(fMaxBeta_deltaMax1, inputValue);
624 }
625 
626 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
627 
628 void Run::SetMaxBeta_deltaMax2(const G4double inputValue)
629 {
630   fMaxBeta_deltaMax2 = std::max(fMaxBeta_deltaMax2, inputValue);
631 }
632 
633 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
634 
635 void Run::SetMaxGamma_deltaMax1(const G4double inputValue)
636 {
637   fMaxGamma_deltaMax1 = std::max(fMaxGamma_deltaMax1, inputValue);
638 }
639 
640 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
641 
642 void Run::SetMaxGamma_deltaMax2(const G4double inputValue)
643 {
644   fMaxGamma_deltaMax2 = std::max(fMaxGamma_deltaMax2, inputValue);
645 }
646 
647 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
648 
649 void Run::SetMaxGamma_deltaMax3(const G4double inputValue)
650 {
651   fMaxGamma_deltaMax3 = std::max(fMaxGamma_deltaMax3, inputValue);
652 }
653 
654 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
655 
656 void Run::SetMaxT_proper_deltaMax(const G4double inputValue)
657 {
658   fMaxT_proper_deltaMax = std::max(fMaxT_proper_deltaMax, inputValue);
659 }
660 
661 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
662 
663 void Run::SetMaxT_lab_deltaMax(const G4double inputValue)
664 {
665   fMaxT_lab_deltaMax = std::max(fMaxT_lab_deltaMax, inputValue);
666 }
667 
668 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
669 
670 void Run::SetMaxMc_truth_rPos_deltaMax(const G4double inputValue)
671 {
672   fSumMc_truth_rPos_deltaMax += inputValue;
673   fMaxMc_truth_rPos_deltaMax = std::max(fMaxMc_truth_rPos_deltaMax, inputValue);
674 }
675 
676 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
677 
678 void Run::SetMinUnderestimated_mc_truth_rPos_delta(const G4double inputValue)
679 {
680   fSumUnderestimated_mc_truth_rPos_delta += inputValue;
681   fMinUnderestimated_mc_truth_rPos_delta =
682     std::min(fMinUnderestimated_mc_truth_rPos_delta, inputValue);
683 }
684 
685 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
686 
687 void Run::SetMaxOverestimated_mc_truth_rPos_delta(const G4double inputValue)
688 {
689   fSumOverestimated_mc_truth_rPos_delta += inputValue;
690   fMaxOverestimated_mc_truth_rPos_delta =
691     std::max(fMaxOverestimated_mc_truth_rPos_delta, inputValue);
692 }
693 
694 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
695 
696 void Run::SetMinUnderestimated_rDeltaPos(const G4double inputValue)
697 {
698   fSumUnderestimated_rDeltaPos += inputValue;
699   fMinUnderestimated_rDeltaPos = std::min(fMinUnderestimated_rDeltaPos, inputValue);
700 }
701 
702 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
703 
704 void Run::SetMaxOverestimated_rDeltaPos(const G4double inputValue)
705 {
706   fSumOverestimated_rDeltaPos += inputValue;
707   fMaxOverestimated_rDeltaPos = std::max(fMaxOverestimated_rDeltaPos, inputValue);
708 }
709 
710 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
711 
712 void Run::SetMaxFloat_rDeltaPos_deltaMax(const G4double inputValue)
713 {
714   fMaxFloat_rDeltaPos_deltaMax = std::max(fMaxFloat_rDeltaPos_deltaMax, inputValue);
715 }
716 
717 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
718