File indexing completed on 2025-04-04 08:05:20
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027
0028
0029
0030
0031
0032
0033
0034 #include "SteppingAction.hh"
0035
0036 #include "Run.hh"
0037
0038 #include "G4DecayProducts.hh"
0039 #include "G4DecayTable.hh"
0040 #include "G4LossTableManager.hh"
0041 #include "G4ParticleDefinition.hh"
0042 #include "G4ParticleTypes.hh"
0043 #include "G4Step.hh"
0044 #include "G4StepPoint.hh"
0045 #include "G4SystemOfUnits.hh"
0046 #include "G4TouchableHistory.hh"
0047 #include "G4Track.hh"
0048 #include "G4VDecayChannel.hh"
0049 #include "G4VPhysicalVolume.hh"
0050 #include "G4VTouchable.hh"
0051
0052
0053
0054 SteppingAction::SteppingAction() : G4UserSteppingAction()
0055 {
0056 Initialize();
0057 }
0058
0059
0060
0061 SteppingAction::~SteppingAction() {}
0062
0063
0064
0065 void SteppingAction::Initialize()
0066 {
0067
0068 fRunPtr = nullptr;
0069 fToleranceEPviolations = 1.0 * CLHEP::eV;
0070 fPrimaryParticleId = 0;
0071 fPrimaryParticleInitialKineticEnergy = 0.0;
0072 fPrimaryParticleInitialTotalEnergy = 0.0;
0073 fPrimaryParticleInitialMomentum = 0.0;
0074 fPrimaryParticleInitialBeta = 1.0;
0075 fPrimaryParticleInitialGamma = 1.0;
0076 fPrimaryParticleInitial3Momentum = G4ThreeVector(0.0, 0.0, 0.0);
0077 fPrimaryParticleInitialPosition = G4ThreeVector(0.0, 0.0, 0.0);
0078 fMaxEkin_deltaMax = 0.0;
0079 fMaxEtot_deltaMax = 0.0;
0080 fMaxP_deltaMax = 0.0;
0081 fMaxPdir_deltaMax = 0.0;
0082 fMaxMass_deltaMax1 = 0.0;
0083 fMaxMass_deltaMax2 = 0.0;
0084 fMaxMass_deltaMax3 = 0.0;
0085 fMeanMass_deltaMax3 = 0.0;
0086 fMaxBeta_deltaMax1 = 0.0;
0087 fMaxBeta_deltaMax2 = 0.0;
0088 fMaxGamma_deltaMax1 = 0.0;
0089 fMaxGamma_deltaMax2 = 0.0;
0090 fMaxGamma_deltaMax3 = 0.0;
0091 fMaxT_proper_deltaMax = 0.0;
0092 fMaxT_lab_deltaMax = 0.0;
0093 fMaxMc_truth_rPos_deltaMax = 0.0;
0094 fMeanMc_truth_rPos_deltaMax = 0.0;
0095 fMeanDeltaR_primaryDecay = 0.0;
0096 fMinDeltaR_primaryDecay = 9999999.9;
0097 fMaxDeltaR_primaryDecay = -9999999.9;
0098 fMeanR_primaryDecay = 0.0;
0099 fMinR_primaryDecay = 9999999.9;
0100 fMaxR_primaryDecay = -9999999.9;
0101 fMeanX_primaryDecay = 0.0;
0102 fMinX_primaryDecay = 9999999.9;
0103 fMaxX_primaryDecay = -9999999.9;
0104 fMeanY_primaryDecay = 0.0;
0105 fMinY_primaryDecay = 9999999.9;
0106 fMaxY_primaryDecay = -9999999.9;
0107 fMeanZ_primaryDecay = 0.0;
0108 fMinZ_primaryDecay = 9999999.9;
0109 fMaxZ_primaryDecay = -9999999.9;
0110 fMeanDeltaAngle_primaryDecay = 0.0;
0111 fMinDeltaAngle_primaryDecay = 9999999.9;
0112 fMaxDeltaAngle_primaryDecay = -9999999.9;
0113 fMeanDeltaEkin_primaryDecay = 0.0;
0114 fMinDeltaEkin_primaryDecay = 9999999.9;
0115 fMaxDeltaEkin_primaryDecay = -9999999.9;
0116 fMeanEkin_primaryDecay = 0.0;
0117 fMinEkin_primaryDecay = 9999999.9;
0118 fMaxEkin_primaryDecay = -9999999.9;
0119 fMeanPx_primaryDecay = 0.0;
0120 fMinPx_primaryDecay = 9999999.9;
0121 fMaxPx_primaryDecay = -9999999.9;
0122 fMeanPy_primaryDecay = 0.0;
0123 fMinPy_primaryDecay = 9999999.9;
0124 fMaxPy_primaryDecay = -9999999.9;
0125 fMeanPz_primaryDecay = 0.0;
0126 fMinPz_primaryDecay = 9999999.9;
0127 fMaxPz_primaryDecay = -9999999.9;
0128 fMinUnderestimated_mc_truth_rPos_delta = 9999999.9;
0129 fMaxOverestimated_mc_truth_rPos_delta = -9999999.9;
0130 fMeanUnderestimated_mc_truth_rPos_delta = 0.0;
0131 fMeanOverestimated_mc_truth_rPos_delta = 0.0;
0132 fMinUnderestimated_rDeltaPos = 9999999.9;
0133 fMaxOverestimated_rDeltaPos = -9999999.9;
0134 fMeanUnderestimated_rDeltaPos = 0.0;
0135 fMeanOverestimated_rDeltaPos = 0.0;
0136 fMaxFloat_rDeltaPos_deltaMax = -9999999.9;
0137 fMeanViolationE_primaryDecay = 0.0;
0138 fMinViolationE_primaryDecay = 9999999.9;
0139 fMaxViolationE_primaryDecay = -9999999.9;
0140 fMeanViolationPx_primaryDecay = 0.0;
0141 fMinViolationPx_primaryDecay = 9999999.9;
0142 fMaxViolationPx_primaryDecay = -9999999.9;
0143 fMeanViolationPy_primaryDecay = 0.0;
0144 fMinViolationPy_primaryDecay = 9999999.9;
0145 fMaxViolationPy_primaryDecay = -9999999.9;
0146 fMeanViolationPz_primaryDecay = 0.0;
0147 fMinViolationPz_primaryDecay = 9999999.9;
0148 fMaxViolationPz_primaryDecay = -9999999.9;
0149 }
0150
0151
0152
0153 void SteppingAction::UserSteppingAction(const G4Step* theStep)
0154 {
0155
0156
0157
0158 if (theStep->GetTrack()->GetParentID() == 0 && theStep->GetTrack()->GetCurrentStepNumber() == 1) {
0159 fPrimaryParticleId = theStep->GetTrack()->GetDefinition()->GetPDGEncoding();
0160 fPrimaryParticleInitialKineticEnergy = theStep->GetPreStepPoint()->GetKineticEnergy();
0161 fPrimaryParticleInitialTotalEnergy = theStep->GetPreStepPoint()->GetTotalEnergy();
0162 fPrimaryParticleInitial3Momentum = theStep->GetPreStepPoint()->GetMomentum();
0163 fPrimaryParticleInitialMomentum = fPrimaryParticleInitial3Momentum.mag();
0164 fPrimaryParticleInitialPosition = theStep->GetPreStepPoint()->GetPosition();
0165 fPrimaryParticleInitialBeta = theStep->GetPreStepPoint()->GetBeta();
0166 fPrimaryParticleInitialGamma = theStep->GetPreStepPoint()->GetGamma();
0167
0168
0169 if (fToleranceEPviolations < fPrimaryParticleInitialKineticEnergy * 1.0e-9) {
0170 fToleranceEPviolations = fPrimaryParticleInitialKineticEnergy * 1.0e-9;
0171 }
0172
0173 if (fRunPtr) {
0174 fRunPtr->SetPrimaryParticleId(fPrimaryParticleId);
0175 fRunPtr->SetPrimaryParticleInitialKineticEnergy(fPrimaryParticleInitialKineticEnergy);
0176 fRunPtr->SetPrimaryParticleInitialTotalEnergy(fPrimaryParticleInitialTotalEnergy);
0177 fRunPtr->SetPrimaryParticleInitialMomentum(fPrimaryParticleInitialMomentum);
0178 fRunPtr->SetPrimaryParticleInitialBeta(fPrimaryParticleInitialBeta);
0179 fRunPtr->SetPrimaryParticleInitialGamma(fPrimaryParticleInitialGamma);
0180 fRunPtr->SetPrimaryParticleInitial3Momentum(fPrimaryParticleInitial3Momentum);
0181 fRunPtr->SetPrimaryParticleInitialPosition(fPrimaryParticleInitialPosition);
0182 fRunPtr->SetToleranceEPviolations(ToleranceEPviolations());
0183 fRunPtr->SetToleranceDeltaDecayRadius(ToleranceDeltaDecayRadius());
0184 fRunPtr->SetIsPreassignedDecayEnabled(IsPreassignedDecayEnabled());
0185 fRunPtr->SetIsBoostToLabEnabled(IsBoostToLabEnabled());
0186 }
0187
0188 if (IsPreassignedDecayEnabled() && (!theStep->GetTrack()->GetDefinition()->GetPDGStable())) {
0189 G4DynamicParticle* dynamicParent =
0190 const_cast<G4DynamicParticle*>(theStep->GetTrack()->GetDynamicParticle());
0191 if (dynamicParent != nullptr) {
0192 G4DecayProducts* decayProducts =
0193 (G4DecayProducts*)(dynamicParent->GetPreAssignedDecayProducts());
0194 if (decayProducts == nullptr) {
0195 G4ParticleDefinition* parentDef = theStep->GetTrack()->GetDefinition();
0196 G4DecayTable* decayTable = (parentDef == nullptr ? nullptr : parentDef->GetDecayTable());
0197 if (decayTable != nullptr) {
0198 G4double parentMass = dynamicParent->GetMass();
0199 G4VDecayChannel* decayChannel = decayTable->SelectADecayChannel(parentMass);
0200 if (decayChannel != nullptr) {
0201 decayProducts = decayChannel->DecayIt(parentMass);
0202 if (!decayProducts->IsChecked()) decayProducts->DumpInfo();
0203 if (IsBoostToLabEnabled()) {
0204
0205 decayProducts->Boost(dynamicParent->GetTotalEnergy(),
0206 dynamicParent->GetMomentumDirection());
0207 }
0208 }
0209 else {
0210 decayProducts = new G4DecayProducts(*dynamicParent);
0211 }
0212 dynamicParent->SetPreAssignedDecayProducts(decayProducts);
0213 }
0214 }
0215 else {
0216 G4cout << "WARNING : already present preassign decay !" << G4endl;
0217 }
0218 }
0219 }
0220 }
0221
0222
0223
0224
0225 if (theStep->GetTrack()->GetParentID() == 0
0226 && theStep->GetPostStepPoint()->GetProcessDefinedStep() != nullptr
0227 && theStep->GetPostStepPoint()->GetProcessDefinedStep()->GetProcessName().find("Decay")
0228 != std::string::npos)
0229 {
0230
0231
0232
0233
0234 const G4double ekin_dynamicParticle =
0235 theStep->GetTrack()->GetDynamicParticle()->GetKineticEnergy();
0236 const G4double ekin_track = theStep->GetTrack()->GetKineticEnergy();
0237 const G4double ekin_postStepPoint = theStep->GetPostStepPoint()->GetKineticEnergy();
0238 const G4double ekin_deltaMax = std::max(std::abs(ekin_dynamicParticle - ekin_track),
0239 std::abs(ekin_dynamicParticle - ekin_postStepPoint));
0240
0241 const G4double ekin_val = ekin_dynamicParticle;
0242
0243 const G4double etot_dynamicParticle =
0244 theStep->GetTrack()->GetDynamicParticle()->GetTotalEnergy();
0245 const G4double etot_track = theStep->GetTrack()->GetTotalEnergy();
0246 const G4double etot_postStepPoint = theStep->GetPostStepPoint()->GetTotalEnergy();
0247 const G4double etot_deltaMax = std::max(std::abs(etot_dynamicParticle - etot_track),
0248 std::abs(etot_dynamicParticle - etot_postStepPoint));
0249
0250 const G4double etot_val = etot_dynamicParticle;
0251
0252 const G4double p_dynamicParticle =
0253 theStep->GetTrack()->GetDynamicParticle()->GetMomentum().mag();
0254 const G4double p_track = theStep->GetTrack()->GetMomentum().mag();
0255 const G4double p_postStepPoint = theStep->GetPostStepPoint()->GetMomentum().mag();
0256 const G4double p_deltaMax = std::max(std::abs(p_dynamicParticle - p_track),
0257 std::abs(p_dynamicParticle - p_postStepPoint));
0258
0259 const G4double p_val = p_dynamicParticle;
0260
0261 const G4ThreeVector pdir_dynamicParticle =
0262 theStep->GetTrack()->GetDynamicParticle()->GetMomentumDirection();
0263 const G4ThreeVector pdir_track = theStep->GetTrack()->GetMomentumDirection();
0264 const G4ThreeVector pdir_postStepPoint = theStep->GetPostStepPoint()->GetMomentumDirection();
0265 const G4double pdir_x_deltaMax =
0266 std::max(std::abs(pdir_dynamicParticle.x() - pdir_track.x()),
0267 std::abs(pdir_dynamicParticle.x() - pdir_postStepPoint.x()));
0268 const G4double pdir_y_deltaMax =
0269 std::max(std::abs(pdir_dynamicParticle.y() - pdir_track.y()),
0270 std::abs(pdir_dynamicParticle.y() - pdir_postStepPoint.y()));
0271 const G4double pdir_z_deltaMax =
0272 std::max(std::abs(pdir_dynamicParticle.z() - pdir_track.z()),
0273 std::abs(pdir_dynamicParticle.z() - pdir_postStepPoint.z()));
0274 const G4double pdir_deltaMax =
0275 std::max(std::max(pdir_x_deltaMax, pdir_y_deltaMax), pdir_z_deltaMax);
0276
0277
0278 const G4double mass_dynamicParticle = theStep->GetTrack()->GetDynamicParticle()->GetMass();
0279 const G4double mass_preStepPoint = theStep->GetPreStepPoint()->GetMass();
0280 const G4double mass_postStepPoint = theStep->GetPostStepPoint()->GetMass();
0281 const G4double mass_from_etot_ekin = etot_val - ekin_val;
0282 const G4double mass_from4mom = std::sqrt(etot_val * etot_val - p_val * p_val);
0283 G4double mass_deltaMax1 = std::max(std::abs(mass_dynamicParticle - mass_preStepPoint),
0284 std::abs(mass_dynamicParticle - mass_postStepPoint));
0285 G4double mass_deltaMax2 = std::abs(mass_dynamicParticle - mass_from_etot_ekin);
0286 G4double mass_deltaMax3 = std::abs(mass_dynamicParticle - mass_from4mom);
0287 fMeanMass_deltaMax3 += mass_deltaMax3;
0288
0289
0290 const G4double mass_val = mass_dynamicParticle;
0291
0292
0293 const G4double beta_dynamicParticle = theStep->GetTrack()->GetDynamicParticle()->GetBeta();
0294 const G4double beta_postStepPoint = theStep->GetPostStepPoint()->GetBeta();
0295
0296 const G4double beta_velocity_track = theStep->GetTrack()->GetVelocity() / CLHEP::c_light;
0297 const G4double beta_velocity_postStepPoint =
0298 theStep->GetPostStepPoint()->GetVelocity() / CLHEP::c_light;
0299 const G4double beta_p_over_etot = p_val / etot_val;
0300 G4double beta_deltaMax1 = std::max(std::abs(beta_dynamicParticle - beta_postStepPoint),
0301 std::abs(beta_dynamicParticle - beta_velocity_track));
0302 beta_deltaMax1 =
0303 std::max(beta_deltaMax1, std::abs(beta_dynamicParticle - beta_velocity_postStepPoint));
0304 const G4double beta_deltaMax2 = std::abs(beta_dynamicParticle - beta_p_over_etot);
0305
0306 const G4double beta_val = beta_dynamicParticle;
0307
0308 const G4double gamma_postStepPoint = theStep->GetPostStepPoint()->GetGamma();
0309 const G4double gamma_from_e_over_m = etot_val / mass_val;
0310 const G4double gamma_deltaMax1 = std::abs(gamma_postStepPoint - gamma_from_e_over_m);
0311 G4double gamma_from_beta = 0.0;
0312 G4double gamma_deltaMax2 = 0.0;
0313 G4double gamma_deltaMax3 = 0.0;
0314 if (beta_val < 1.0) {
0315 gamma_from_beta = 1.0 / std::sqrt(1.0 - beta_val * beta_val);
0316 gamma_deltaMax2 = std::abs(gamma_postStepPoint - gamma_from_beta);
0317 gamma_deltaMax3 = std::abs(gamma_from_e_over_m - gamma_from_beta);
0318 }
0319 const G4double gamma_val = gamma_postStepPoint;
0320
0321
0322
0323
0324
0325 const G4double t_proper_track = theStep->GetTrack()->GetProperTime();
0326 const G4double t_proper_postStepPoint = theStep->GetPostStepPoint()->GetProperTime();
0327 const G4double t_proper_deltaMax = std::abs(t_proper_track - t_proper_postStepPoint);
0328
0329 const G4double t_proper_val = t_proper_track;
0330
0331
0332
0333
0334
0335
0336 const G4double t_local_track = theStep->GetTrack()->GetLocalTime();
0337 const G4double t_local_postStepPoint = theStep->GetPostStepPoint()->GetLocalTime();
0338 const G4double t_global_track = theStep->GetTrack()->GetGlobalTime();
0339 const G4double t_global_postStepPoint = theStep->GetPostStepPoint()->GetGlobalTime();
0340 G4double t_lab_deltaMax = std::max(std::abs(t_local_track - t_local_postStepPoint),
0341 std::abs(t_local_track - t_global_track));
0342 t_lab_deltaMax = std::max(t_lab_deltaMax, std::abs(t_local_track - t_global_postStepPoint));
0343
0344 const G4double t_lab_val = t_local_track;
0345
0346
0347
0348 const G4double primaryBeta =
0349 fPrimaryParticleInitialMomentum / fPrimaryParticleInitialTotalEnergy;
0350 const G4double mc_truth_rPos1 = t_lab_val * fPrimaryParticleInitialBeta * CLHEP::c_light;
0351 const G4double mc_truth_rPos2 = t_lab_val * primaryBeta * CLHEP::c_light;
0352 const G4double mc_truth_rPos_deltaMax = std::abs(mc_truth_rPos1 - mc_truth_rPos2);
0353 fMeanMc_truth_rPos_deltaMax += mc_truth_rPos_deltaMax;
0354
0355
0356 if (mc_truth_rPos_deltaMax > ToleranceDeltaDecayRadius()) {
0357
0358
0359
0360
0361 if (fRunPtr) fRunPtr->IncrementNumber_mc_truth_rPos_deltaMax_above();
0362 }
0363 const G4double mc_truth_rPos_val = mc_truth_rPos1;
0364
0365 fMaxEkin_deltaMax = std::max(fMaxEkin_deltaMax, ekin_deltaMax);
0366 fMaxEtot_deltaMax = std::max(fMaxEtot_deltaMax, etot_deltaMax);
0367 fMaxP_deltaMax = std::max(fMaxP_deltaMax, p_deltaMax);
0368 fMaxPdir_deltaMax = std::max(fMaxPdir_deltaMax, pdir_deltaMax);
0369 fMaxMass_deltaMax1 = std::max(fMaxMass_deltaMax1, mass_deltaMax1);
0370 fMaxMass_deltaMax2 = std::max(fMaxMass_deltaMax2, mass_deltaMax2);
0371 fMaxMass_deltaMax3 = std::max(fMaxMass_deltaMax3, mass_deltaMax3);
0372 fMaxBeta_deltaMax1 = std::max(fMaxBeta_deltaMax1, beta_deltaMax1);
0373 fMaxBeta_deltaMax2 = std::max(fMaxBeta_deltaMax2, beta_deltaMax2);
0374 fMaxGamma_deltaMax1 = std::max(fMaxGamma_deltaMax1, gamma_deltaMax1);
0375 fMaxGamma_deltaMax2 = std::max(fMaxGamma_deltaMax2, gamma_deltaMax2);
0376 fMaxGamma_deltaMax3 = std::max(fMaxGamma_deltaMax3, gamma_deltaMax3);
0377 fMaxT_lab_deltaMax = std::max(fMaxT_lab_deltaMax, t_lab_deltaMax);
0378 fMaxT_proper_deltaMax = std::max(fMaxT_proper_deltaMax, t_proper_deltaMax);
0379 fMaxMc_truth_rPos_deltaMax = std::max(fMaxMc_truth_rPos_deltaMax, mc_truth_rPos_deltaMax);
0380
0381
0382
0383 const G4double xPos = theStep->GetPostStepPoint()->GetPosition().x();
0384 const G4double yPos = theStep->GetPostStepPoint()->GetPosition().y();
0385 const G4double zPos = theStep->GetPostStepPoint()->GetPosition().z();
0386 const G4double rPos = std::sqrt(xPos * xPos + yPos * yPos + zPos * zPos);
0387
0388
0389
0390 G4double tPos = theStep->GetPostStepPoint()->GetLocalTime();
0391
0392
0393 const G4double mc_truth_rPos = tPos * fPrimaryParticleInitialBeta * CLHEP::c_light;
0394 const G4double rDeltaPos = mc_truth_rPos - rPos;
0395 const G4double eKin = theStep->GetPostStepPoint()->GetKineticEnergy();
0396 const G4double xMom = theStep->GetPostStepPoint()->GetMomentum().x();
0397 const G4double yMom = theStep->GetPostStepPoint()->GetMomentum().y();
0398 const G4double zMom = theStep->GetPostStepPoint()->GetMomentum().z();
0399
0400
0401 G4double xDirection = std::min(theStep->GetPostStepPoint()->GetMomentumDirection().x(), 1.0);
0402 if (xDirection < -1.0) xDirection = -1.0;
0403 const G4double deflection_angle_in_degrees = 57.29 * std::acos(xDirection);
0404 const G4double delta_ekin = fPrimaryParticleInitialKineticEnergy - eKin;
0405
0406
0407
0408
0409
0410
0411
0412
0413
0414
0415 if (std::abs(rDeltaPos) > ToleranceDeltaDecayRadius()) {
0416
0417
0418 if (rDeltaPos < 0.0) {
0419
0420 if (fRunPtr) fRunPtr->IncrementNumberUnexpectedDecays();
0421 }
0422
0423 }
0424 fMeanDeltaR_primaryDecay += rDeltaPos;
0425 fMinDeltaR_primaryDecay = std::min(fMinDeltaR_primaryDecay, rDeltaPos);
0426 fMaxDeltaR_primaryDecay = std::max(fMaxDeltaR_primaryDecay, rDeltaPos);
0427 fMeanR_primaryDecay += rPos;
0428 fMinR_primaryDecay = std::min(fMinR_primaryDecay, rPos);
0429 fMaxR_primaryDecay = std::max(fMaxR_primaryDecay, rPos);
0430 fMeanX_primaryDecay += xPos;
0431 fMinX_primaryDecay = std::min(fMinX_primaryDecay, xPos);
0432 fMaxX_primaryDecay = std::max(fMaxX_primaryDecay, xPos);
0433 fMeanY_primaryDecay += yPos;
0434 fMinY_primaryDecay = std::min(fMinY_primaryDecay, yPos);
0435 fMaxY_primaryDecay = std::max(fMaxY_primaryDecay, yPos);
0436 fMeanZ_primaryDecay += zPos;
0437 fMinZ_primaryDecay = std::min(fMinZ_primaryDecay, zPos);
0438 fMaxZ_primaryDecay = std::max(fMaxZ_primaryDecay, zPos);
0439 fMeanDeltaAngle_primaryDecay += deflection_angle_in_degrees;
0440 fMinDeltaAngle_primaryDecay =
0441 std::min(fMinDeltaAngle_primaryDecay, deflection_angle_in_degrees);
0442 fMaxDeltaAngle_primaryDecay =
0443 std::max(fMaxDeltaAngle_primaryDecay, deflection_angle_in_degrees);
0444 fMeanDeltaEkin_primaryDecay += delta_ekin;
0445 fMinDeltaEkin_primaryDecay = std::min(fMinDeltaEkin_primaryDecay, delta_ekin);
0446 fMaxDeltaEkin_primaryDecay = std::max(fMaxDeltaEkin_primaryDecay, delta_ekin);
0447 fMeanEkin_primaryDecay += eKin;
0448 fMinEkin_primaryDecay = std::min(fMinEkin_primaryDecay, eKin);
0449 fMaxEkin_primaryDecay = std::max(fMaxEkin_primaryDecay, eKin);
0450 fMeanPx_primaryDecay += xMom;
0451 fMinPx_primaryDecay = std::min(fMinPx_primaryDecay, xMom);
0452 fMaxPx_primaryDecay = std::max(fMaxPx_primaryDecay, xMom);
0453 fMeanPy_primaryDecay += yMom;
0454 fMinPy_primaryDecay = std::min(fMinPy_primaryDecay, yMom);
0455 fMaxPy_primaryDecay = std::max(fMaxPy_primaryDecay, yMom);
0456 fMeanPz_primaryDecay += zMom;
0457 fMinPz_primaryDecay = std::min(fMinPz_primaryDecay, zMom);
0458 fMaxPz_primaryDecay = std::max(fMaxPz_primaryDecay, zMom);
0459
0460
0461
0462
0463
0464
0465
0466
0467
0468
0469 const G4double underestimated_mc_truth_rPos =
0470 t_proper_val * gamma_val * beta_val * CLHEP::c_light;
0471 const G4double overestimated_mc_truth_rPos =
0472 t_proper_val * fPrimaryParticleInitialGamma * fPrimaryParticleInitialBeta * CLHEP::c_light;
0473 const G4double underestimated_mc_truth_rPos_delta =
0474 underestimated_mc_truth_rPos - mc_truth_rPos_val;
0475 const G4double overestimated_mc_truth_rPos_delta =
0476 overestimated_mc_truth_rPos - mc_truth_rPos_val;
0477 fMeanUnderestimated_mc_truth_rPos_delta += underestimated_mc_truth_rPos_delta;
0478 fMeanOverestimated_mc_truth_rPos_delta += overestimated_mc_truth_rPos_delta;
0479
0480
0481
0482
0483 if (-underestimated_mc_truth_rPos_delta > ToleranceDeltaDecayRadius()) {
0484
0485
0486
0487
0488
0489 if (fRunPtr) fRunPtr->IncrementNumber_underestimated_mc_truth_rPos_delta_above();
0490 }
0491 if (overestimated_mc_truth_rPos_delta > ToleranceDeltaDecayRadius()) {
0492
0493
0494
0495
0496
0497 if (fRunPtr) fRunPtr->IncrementNumber_overestimated_mc_truth_rPos_delta_above();
0498 }
0499 const G4double underestimated_rDeltaPos = underestimated_mc_truth_rPos - rPos;
0500 const G4double overestimated_rDeltaPos = overestimated_mc_truth_rPos - rPos;
0501 fMeanUnderestimated_rDeltaPos += underestimated_rDeltaPos;
0502 fMeanOverestimated_rDeltaPos += overestimated_rDeltaPos;
0503
0504
0505
0506
0507 if (-underestimated_rDeltaPos > ToleranceDeltaDecayRadius()) {
0508 if (fRunPtr) fRunPtr->IncrementNumberLargeUnderestimates();
0509 }
0510 if (overestimated_rDeltaPos > ToleranceDeltaDecayRadius()) {
0511 if (fRunPtr) fRunPtr->IncrementNumberLargeOverestimates();
0512 }
0513
0514 fMinUnderestimated_mc_truth_rPos_delta =
0515 std::min(fMinUnderestimated_mc_truth_rPos_delta, underestimated_mc_truth_rPos_delta);
0516 fMaxOverestimated_mc_truth_rPos_delta =
0517 std::max(fMaxOverestimated_mc_truth_rPos_delta, overestimated_mc_truth_rPos_delta);
0518 fMinUnderestimated_rDeltaPos = std::min(fMinUnderestimated_rDeltaPos, underestimated_rDeltaPos);
0519 fMaxOverestimated_rDeltaPos = std::max(fMaxOverestimated_rDeltaPos, overestimated_rDeltaPos);
0520
0521
0522
0523
0524 const G4float float_xPos = static_cast<G4float>(theStep->GetPostStepPoint()->GetPosition().x());
0525 const G4float float_yPos = static_cast<G4float>(theStep->GetPostStepPoint()->GetPosition().y());
0526 const G4float float_zPos = static_cast<G4float>(theStep->GetPostStepPoint()->GetPosition().z());
0527 const G4float float_rPos =
0528 std::sqrt(float_xPos * float_xPos + float_yPos * float_yPos + float_zPos * float_zPos);
0529 const G4float float_tPos = static_cast<G4float>(theStep->GetPostStepPoint()->GetLocalTime());
0530 const G4float float_initialBeta1 = static_cast<G4float>(fPrimaryParticleInitialBeta);
0531 const G4float float_initialBeta2 = static_cast<G4float>(fPrimaryParticleInitialMomentum)
0532 / static_cast<G4float>(fPrimaryParticleInitialTotalEnergy);
0533 const G4float float_initialGamma = static_cast<G4float>(fPrimaryParticleInitialGamma);
0534 const G4float float_initialBeta3 =
0535 std::sqrt(float_initialGamma * float_initialGamma - 1.0) / float_initialGamma;
0536 const G4float float_c_light = static_cast<G4float>(CLHEP::c_light);
0537 const G4float float_mc_truth_rPos1 = float_tPos * float_initialBeta1 * float_c_light;
0538 const G4float float_mc_truth_rPos2 = float_tPos * float_initialBeta2 * float_c_light;
0539 const G4float float_mc_truth_rPos3 = float_tPos * float_initialBeta3 * float_c_light;
0540 const G4float float_rDeltaPos_0 = static_cast<G4float>(rDeltaPos);
0541 const G4float float_rDeltaPos_1 = float_mc_truth_rPos1 - float_rPos;
0542 const G4float float_rDeltaPos_2 = float_mc_truth_rPos2 - float_rPos;
0543 const G4float float_rDeltaPos_3 = float_mc_truth_rPos3 - float_rPos;
0544 const G4float float_rDeltaPos_4 = static_cast<G4float>(mc_truth_rPos) - float_rPos;
0545 const G4float float_rDeltaPos_5 = float_mc_truth_rPos1 - static_cast<G4float>(rPos);
0546 const G4float float_rDeltaPos_6 = float_mc_truth_rPos2 - static_cast<G4float>(rPos);
0547 const G4float float_rDeltaPos_7 = float_mc_truth_rPos3 - static_cast<G4float>(rPos);
0548 G4double rDeltaPos_deltaMax =
0549 std::max(std::abs(float_rDeltaPos_0 - rDeltaPos), std::abs(float_rDeltaPos_1 - rDeltaPos));
0550 rDeltaPos_deltaMax = std::max(rDeltaPos_deltaMax, std::abs(float_rDeltaPos_2 - rDeltaPos));
0551 rDeltaPos_deltaMax = std::max(rDeltaPos_deltaMax, std::abs(float_rDeltaPos_3 - rDeltaPos));
0552 rDeltaPos_deltaMax = std::max(rDeltaPos_deltaMax, std::abs(float_rDeltaPos_4 - rDeltaPos));
0553 rDeltaPos_deltaMax = std::max(rDeltaPos_deltaMax, std::abs(float_rDeltaPos_5 - rDeltaPos));
0554 rDeltaPos_deltaMax = std::max(rDeltaPos_deltaMax, std::abs(float_rDeltaPos_6 - rDeltaPos));
0555 rDeltaPos_deltaMax = std::max(rDeltaPos_deltaMax, std::abs(float_rDeltaPos_7 - rDeltaPos));
0556
0557
0558 fMaxFloat_rDeltaPos_deltaMax = std::max(fMaxFloat_rDeltaPos_deltaMax, rDeltaPos_deltaMax);
0559
0560
0561
0562 std::size_t nSec = theStep->GetNumberOfSecondariesInCurrentStep();
0563 const std::vector<const G4Track*>* ptrVecSecondaries = theStep->GetSecondaryInCurrentStep();
0564 G4double deltaE = 0.0, deltaPx = 0.0, deltaPy = 0.0, deltaPz = 0.0;
0565 if (nSec > 0 && ptrVecSecondaries != nullptr) {
0566 G4double sumEsecondaries = 0.0;
0567 G4ThreeVector sumPsecondaries(0.0, 0.0, 0.0);
0568 for (std::size_t i = 0; i < nSec; ++i) {
0569 if ((*ptrVecSecondaries)[i]) {
0570 sumEsecondaries += (*ptrVecSecondaries)[i]->GetTotalEnergy();
0571 sumPsecondaries += (*ptrVecSecondaries)[i]->GetMomentum();
0572 }
0573 }
0574 deltaE = sumEsecondaries - theStep->GetPostStepPoint()->GetTotalEnergy();
0575 fMeanViolationE_primaryDecay += deltaE;
0576 fMinViolationE_primaryDecay = std::min(fMinViolationE_primaryDecay, deltaE);
0577 fMaxViolationE_primaryDecay = std::max(fMaxViolationE_primaryDecay, deltaE);
0578 if (std::abs(deltaE) > ToleranceEPviolations()) {
0579 if (fRunPtr) fRunPtr->IncrementNumberEviolations();
0580 }
0581 deltaPx = sumPsecondaries.x() - xMom;
0582 fMeanViolationPx_primaryDecay += deltaPx;
0583 fMinViolationPx_primaryDecay = std::min(fMinViolationPx_primaryDecay, deltaPx);
0584 fMaxViolationPx_primaryDecay = std::max(fMaxViolationPx_primaryDecay, deltaPx);
0585 deltaPy = sumPsecondaries.y() - yMom;
0586 fMeanViolationPy_primaryDecay += deltaPy;
0587 fMinViolationPy_primaryDecay = std::min(fMinViolationPy_primaryDecay, deltaPy);
0588 fMaxViolationPy_primaryDecay = std::max(fMaxViolationPy_primaryDecay, deltaPy);
0589 deltaPz = sumPsecondaries.z() - zMom;
0590 fMeanViolationPz_primaryDecay += deltaPz;
0591 fMinViolationPz_primaryDecay = std::min(fMinViolationPz_primaryDecay, deltaPz);
0592 fMaxViolationPz_primaryDecay = std::max(fMaxViolationPz_primaryDecay, deltaPz);
0593 if (std::abs(deltaPx) > ToleranceEPviolations() || std::abs(deltaPy) > ToleranceEPviolations()
0594 || std::abs(deltaPz) > ToleranceEPviolations())
0595 {
0596 if (fRunPtr) fRunPtr->IncrementNumberPviolations();
0597 }
0598 }
0599 else {
0600 if (fRunPtr) fRunPtr->IncrementNumberBadPrimaryDecays();
0601 }
0602
0603 if (fRunPtr) {
0604 fRunPtr->IncrementNumberDecays();
0605 fRunPtr->SetDecayT(tPos);
0606 fRunPtr->SetDecayR_mc_truth(mc_truth_rPos);
0607 fRunPtr->SetDecayR(rPos);
0608 fRunPtr->SetDecayX(xPos);
0609 fRunPtr->SetDecayY(yPos);
0610 fRunPtr->SetDecayZ(zPos);
0611 fRunPtr->SetDeltaDecayR(rDeltaPos);
0612 fRunPtr->SetDeflectionAngle(deflection_angle_in_degrees);
0613 fRunPtr->SetDeltaEkin(delta_ekin);
0614 fRunPtr->SetDecayEkin(eKin);
0615 fRunPtr->SetDecayPx(xMom);
0616 fRunPtr->SetDecayPy(yMom);
0617 fRunPtr->SetDecayPz(zMom);
0618 fRunPtr->SetDecayEtotViolation(deltaE);
0619 fRunPtr->SetDecayPxViolation(deltaPx);
0620 fRunPtr->SetDecayPyViolation(deltaPy);
0621 fRunPtr->SetDecayPzViolation(deltaPz);
0622 fRunPtr->SetMaxEkin_deltaMax(ekin_deltaMax);
0623 fRunPtr->SetMaxEtot_deltaMax(etot_deltaMax);
0624 fRunPtr->SetMaxP_deltaMax(p_deltaMax);
0625 fRunPtr->SetMaxPdir_deltaMax(pdir_deltaMax);
0626 fRunPtr->SetMaxMass_deltaMax1(mass_deltaMax1);
0627 fRunPtr->SetMaxMass_deltaMax2(mass_deltaMax2);
0628 fRunPtr->SetMaxMass_deltaMax3(mass_deltaMax3);
0629 fRunPtr->SetMaxBeta_deltaMax1(beta_deltaMax1);
0630 fRunPtr->SetMaxBeta_deltaMax2(beta_deltaMax2);
0631 fRunPtr->SetMaxGamma_deltaMax1(gamma_deltaMax1);
0632 fRunPtr->SetMaxGamma_deltaMax2(gamma_deltaMax2);
0633 fRunPtr->SetMaxGamma_deltaMax3(gamma_deltaMax3);
0634 fRunPtr->SetMaxT_proper_deltaMax(t_proper_deltaMax);
0635 fRunPtr->SetMaxT_lab_deltaMax(t_lab_deltaMax);
0636 fRunPtr->SetMaxMc_truth_rPos_deltaMax(mc_truth_rPos_deltaMax);
0637 fRunPtr->SetMinUnderestimated_mc_truth_rPos_delta(underestimated_mc_truth_rPos_delta);
0638 fRunPtr->SetMaxOverestimated_mc_truth_rPos_delta(overestimated_mc_truth_rPos_delta);
0639 fRunPtr->SetMinUnderestimated_rDeltaPos(underestimated_rDeltaPos);
0640 fRunPtr->SetMaxOverestimated_rDeltaPos(overestimated_rDeltaPos);
0641 fRunPtr->SetMaxFloat_rDeltaPos_deltaMax(fMaxFloat_rDeltaPos_deltaMax);
0642 }
0643 }
0644 }
0645
0646