File indexing completed on 2023-03-27 01:58:44
0001 #include "SimG4Core/Application/interface/SteppingAction.h"
0002
0003 #include "SimG4Core/Notification/interface/TrackInformation.h"
0004 #include "SimG4Core/Notification/interface/CMSSteppingVerbose.h"
0005
0006 #include "G4LogicalVolumeStore.hh"
0007 #include "G4ParticleTable.hh"
0008 #include "G4PhysicalVolumeStore.hh"
0009 #include "G4RegionStore.hh"
0010 #include "G4UnitsTable.hh"
0011 #include "G4SystemOfUnits.hh"
0012
0013 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0014 #include "FWCore/Utilities/interface/isFinite.h"
0015
0016
0017
0018 SteppingAction::SteppingAction(const CMSSteppingVerbose* sv, const edm::ParameterSet& p, bool hasW)
0019 : steppingVerbose(sv), hasWatcher(hasW) {
0020 theCriticalEnergyForVacuum = (p.getParameter<double>("CriticalEnergyForVacuum") * CLHEP::MeV);
0021 if (0.0 < theCriticalEnergyForVacuum) {
0022 killBeamPipe = true;
0023 }
0024 theCriticalDensity = (p.getParameter<double>("CriticalDensity") * CLHEP::g / CLHEP::cm3);
0025 maxZCentralCMS = p.getParameter<double>("MaxZCentralCMS") * CLHEP::m;
0026 maxTrackTime = p.getParameter<double>("MaxTrackTime") * CLHEP::ns;
0027 maxTrackTimeForward = p.getParameter<double>("MaxTrackTimeForward") * CLHEP::ns;
0028 maxTrackTimes = p.getParameter<std::vector<double> >("MaxTrackTimes");
0029 maxTimeNames = p.getParameter<std::vector<std::string> >("MaxTimeNames");
0030 deadRegionNames = p.getParameter<std::vector<std::string> >("DeadRegions");
0031 maxNumberOfSteps = p.getParameter<int>("MaxNumberOfSteps");
0032 ekinMins = p.getParameter<std::vector<double> >("EkinThresholds");
0033 ekinNames = p.getParameter<std::vector<std::string> >("EkinNames");
0034 ekinParticles = p.getParameter<std::vector<std::string> >("EkinParticles");
0035
0036 edm::LogVerbatim("SimG4CoreApplication")
0037 << "SteppingAction:: KillBeamPipe = " << killBeamPipe
0038 << " CriticalDensity = " << theCriticalDensity * CLHEP::cm3 / CLHEP::g << " g/cm3\n"
0039 << " CriticalEnergyForVacuum = " << theCriticalEnergyForVacuum / CLHEP::MeV << " Mev;"
0040 << " MaxTrackTime = " << maxTrackTime / CLHEP::ns << " ns;"
0041 << " MaxZCentralCMS = " << maxZCentralCMS / CLHEP::m << " m"
0042 << " MaxTrackTimeForward = " << maxTrackTimeForward / CLHEP::ns << " ns"
0043 << " MaxNumberOfSteps = " << maxNumberOfSteps;
0044
0045 numberTimes = maxTrackTimes.size();
0046 if (numberTimes > 0) {
0047 for (unsigned int i = 0; i < numberTimes; i++) {
0048 edm::LogVerbatim("SimG4CoreApplication")
0049 << "SteppingAction::MaxTrackTime for " << maxTimeNames[i] << " is " << maxTrackTimes[i] << " ns ";
0050 maxTrackTimes[i] *= ns;
0051 }
0052 }
0053
0054 ndeadRegions = deadRegionNames.size();
0055 if (ndeadRegions > 0) {
0056 edm::LogVerbatim("SimG4CoreApplication")
0057 << "SteppingAction: Number of DeadRegions where all trackes are killed " << ndeadRegions;
0058 for (unsigned int i = 0; i < ndeadRegions; ++i) {
0059 edm::LogVerbatim("SimG4CoreApplication") << "SteppingAction: DeadRegion " << i << ". " << deadRegionNames[i];
0060 }
0061 }
0062 numberEkins = ekinNames.size();
0063 numberPart = ekinParticles.size();
0064 if (0 == numberPart) {
0065 numberEkins = 0;
0066 }
0067
0068 if (numberEkins > 0) {
0069 edm::LogVerbatim("SimG4CoreApplication")
0070 << "SteppingAction::Kill following " << numberPart << " particles in " << numberEkins << " volumes";
0071 for (unsigned int i = 0; i < numberPart; ++i) {
0072 edm::LogVerbatim("SimG4CoreApplication")
0073 << "SteppingAction::Particle " << i << " " << ekinParticles[i] << " Threshold = " << ekinMins[i] << " MeV";
0074 ekinMins[i] *= CLHEP::MeV;
0075 }
0076 for (unsigned int i = 0; i < numberEkins; ++i) {
0077 edm::LogVerbatim("SimG4CoreApplication") << "SteppingAction::LogVolume[" << i << "] = " << ekinNames[i];
0078 }
0079 }
0080 }
0081
0082 void SteppingAction::UserSteppingAction(const G4Step* aStep) {
0083 if (!initialized) {
0084 initialized = initPointer();
0085 }
0086
0087 m_g4StepSignal(aStep);
0088
0089 G4Track* theTrack = aStep->GetTrack();
0090 TrackStatus tstat = (theTrack->GetTrackStatus() == fAlive) ? sAlive : sKilledByProcess;
0091
0092 if (theTrack->GetKineticEnergy() < 0.0) {
0093 if (nWarnings < 2) {
0094 ++nWarnings;
0095 edm::LogWarning("SimG4CoreApplication")
0096 << "SteppingAction::UserSteppingAction: Track #" << theTrack->GetTrackID() << " "
0097 << theTrack->GetDefinition()->GetParticleName() << " Ekin(MeV)= " << theTrack->GetKineticEnergy() / MeV;
0098 }
0099 theTrack->SetKineticEnergy(0.0);
0100 }
0101
0102 const G4StepPoint* preStep = aStep->GetPreStepPoint();
0103 const G4StepPoint* postStep = aStep->GetPostStepPoint();
0104
0105
0106 if (tstat == sKilledByProcess) {
0107 if (nullptr != steppingVerbose) {
0108 steppingVerbose->nextStep(aStep, fpSteppingManager, false);
0109 }
0110 return;
0111 }
0112
0113 if (sAlive == tstat && theTrack->GetCurrentStepNumber() > maxNumberOfSteps) {
0114 tstat = sNumberOfSteps;
0115 if (nWarnings < 5) {
0116 ++nWarnings;
0117 edm::LogWarning("SimG4CoreApplication")
0118 << "Track #" << theTrack->GetTrackID() << " " << theTrack->GetDefinition()->GetParticleName()
0119 << " E(MeV)= " << preStep->GetKineticEnergy() / MeV << " Nstep= " << theTrack->GetCurrentStepNumber()
0120 << " is killed due to limit on number of steps;/n PV: " << preStep->GetPhysicalVolume()->GetName() << " at "
0121 << theTrack->GetPosition() << " StepLen(mm)= " << aStep->GetStepLength();
0122 }
0123 }
0124 const double time = theTrack->GetGlobalTime();
0125
0126
0127 if (sAlive == tstat && std::abs(theTrack->GetPosition().z()) >= maxZCentralCMS) {
0128 tstat = (time > maxTrackTimeForward) ? sOutOfTime : sVeryForward;
0129 }
0130
0131
0132 if (sAlive == tstat) {
0133
0134 const G4LogicalVolume* lv = postStep->GetPhysicalVolume()->GetLogicalVolume();
0135 const G4Region* theRegion = lv->GetRegion();
0136
0137
0138 if (isInsideDeadRegion(theRegion))
0139 tstat = sDeadRegion;
0140
0141
0142 if (sAlive == tstat) {
0143 if (isOutOfTimeWindow(theRegion, time))
0144 tstat = sOutOfTime;
0145 }
0146
0147
0148 if (sAlive == tstat && numberEkins > 0) {
0149 if (isLowEnergy(lv, theTrack))
0150 tstat = sLowEnergy;
0151 }
0152
0153
0154 if (sAlive == tstat && killBeamPipe) {
0155 if (theTrack->GetKineticEnergy() < theCriticalEnergyForVacuum &&
0156 theTrack->GetDefinition()->GetPDGCharge() != 0.0 && lv->GetMaterial()->GetDensity() <= theCriticalDensity) {
0157 tstat = sLowEnergyInVacuum;
0158 }
0159 }
0160 }
0161
0162 bool isKilled = false;
0163 if (sAlive == tstat || sVeryForward == tstat) {
0164 if (preStep->GetPhysicalVolume() == tracker && postStep->GetPhysicalVolume() == calo) {
0165 TrackInformation* trkinfo = static_cast<TrackInformation*>(theTrack->GetUserInformation());
0166 if (!trkinfo->crossedBoundary()) {
0167 trkinfo->setCrossedBoundary(theTrack);
0168 }
0169 }
0170 } else {
0171 theTrack->SetTrackStatus(fStopAndKill);
0172 isKilled = true;
0173 #ifdef DebugLog
0174 PrintKilledTrack(theTrack, tstat);
0175 #endif
0176 }
0177 if (nullptr != steppingVerbose) {
0178 steppingVerbose->nextStep(aStep, fpSteppingManager, isKilled);
0179 }
0180 }
0181
0182 bool SteppingAction::isLowEnergy(const G4LogicalVolume* lv, const G4Track* theTrack) const {
0183 const double ekin = theTrack->GetKineticEnergy();
0184 int pCode = theTrack->GetDefinition()->GetPDGEncoding();
0185
0186 for (auto& vol : ekinVolumes) {
0187 if (lv == vol) {
0188 for (unsigned int i = 0; i < numberPart; ++i) {
0189 if (pCode == ekinPDG[i]) {
0190 return (ekin <= ekinMins[i]);
0191 }
0192 }
0193 break;
0194 }
0195 }
0196 return false;
0197 }
0198
0199 bool SteppingAction::initPointer() {
0200 const G4PhysicalVolumeStore* pvs = G4PhysicalVolumeStore::GetInstance();
0201 for (auto const& pvcite : *pvs) {
0202 const G4String& pvname = pvcite->GetName();
0203 if (pvname == "Tracker" || pvname == "tracker:Tracker_1") {
0204 tracker = pvcite;
0205 } else if (pvname == "CALO" || pvname == "caloBase:CALO_1") {
0206 calo = pvcite;
0207 }
0208 if (tracker && calo)
0209 break;
0210 }
0211 edm::LogVerbatim("SimG4CoreApplication")
0212 << "SteppingAction: pointer for Tracker " << tracker << " and for Calo " << calo;
0213
0214 const G4LogicalVolumeStore* lvs = G4LogicalVolumeStore::GetInstance();
0215 if (numberEkins > 0) {
0216 ekinVolumes.resize(numberEkins, nullptr);
0217 for (auto const& lvcite : *lvs) {
0218 const G4String& lvname = lvcite->GetName();
0219 for (unsigned int i = 0; i < numberEkins; ++i) {
0220 if (lvname == (G4String)(ekinNames[i])) {
0221 ekinVolumes[i] = lvcite;
0222 break;
0223 }
0224 }
0225 }
0226 for (unsigned int i = 0; i < numberEkins; ++i) {
0227 edm::LogVerbatim("SimG4CoreApplication") << ekinVolumes[i]->GetName() << " with pointer " << ekinVolumes[i];
0228 }
0229 }
0230
0231 if (numberPart > 0) {
0232 G4ParticleTable* theParticleTable = G4ParticleTable::GetParticleTable();
0233 ekinPDG.resize(numberPart, 0);
0234 for (unsigned int i = 0; i < numberPart; ++i) {
0235 const G4ParticleDefinition* part = theParticleTable->FindParticle(ekinParticles[i]);
0236 if (nullptr != part)
0237 ekinPDG[i] = part->GetPDGEncoding();
0238 edm::LogVerbatim("SimG4CoreApplication") << "Particle " << ekinParticles[i] << " with PDG code " << ekinPDG[i]
0239 << " and KE cut off " << ekinMins[i] / MeV << " MeV";
0240 }
0241 }
0242
0243 const G4RegionStore* rs = G4RegionStore::GetInstance();
0244 if (numberTimes > 0) {
0245 maxTimeRegions.resize(numberTimes, nullptr);
0246 for (auto const& rcite : *rs) {
0247 const G4String& rname = rcite->GetName();
0248 for (unsigned int i = 0; i < numberTimes; ++i) {
0249 if (rname == (G4String)(maxTimeNames[i])) {
0250 maxTimeRegions[i] = rcite;
0251 break;
0252 }
0253 }
0254 }
0255 }
0256 if (ndeadRegions > 0) {
0257 deadRegions.resize(ndeadRegions, nullptr);
0258 for (auto const& rcite : *rs) {
0259 const G4String& rname = rcite->GetName();
0260 for (unsigned int i = 0; i < ndeadRegions; ++i) {
0261 if (rname == (G4String)(deadRegionNames[i])) {
0262 deadRegions[i] = rcite;
0263 break;
0264 }
0265 }
0266 }
0267 }
0268 return true;
0269 }
0270
0271 void SteppingAction::PrintKilledTrack(const G4Track* aTrack, const TrackStatus& tst) const {
0272 std::string vname = "";
0273 std::string rname = "";
0274 std::string typ = " ";
0275 switch (tst) {
0276 case sDeadRegion:
0277 typ = " in dead region ";
0278 break;
0279 case sOutOfTime:
0280 typ = " out of time window ";
0281 break;
0282 case sLowEnergy:
0283 typ = " low energy limit ";
0284 break;
0285 case sLowEnergyInVacuum:
0286 typ = " low energy limit in vacuum ";
0287 break;
0288 case sEnergyDepNaN:
0289 typ = " energy deposition is NaN ";
0290 break;
0291 case sVeryForward:
0292 typ = " very forward track ";
0293 break;
0294 case sNumberOfSteps:
0295 typ = " too many steps ";
0296 break;
0297 default:
0298 break;
0299 }
0300 G4VPhysicalVolume* pv = aTrack->GetNextVolume();
0301 vname = pv->GetLogicalVolume()->GetName();
0302 rname = pv->GetLogicalVolume()->GetRegion()->GetName();
0303
0304 const double ekin = aTrack->GetKineticEnergy();
0305 if (ekin < 2 * CLHEP::MeV) {
0306 return;
0307 }
0308 edm::LogWarning("SimG4CoreApplication")
0309 << "Track #" << aTrack->GetTrackID() << " StepN= " << aTrack->GetCurrentStepNumber() << " "
0310 << aTrack->GetDefinition()->GetParticleName() << " E(MeV)=" << ekin / CLHEP::MeV
0311 << " T(ns)=" << aTrack->GetGlobalTime() / CLHEP::ns << " is killed due to " << typ << "\n LV: " << vname << " ("
0312 << rname << ") at " << aTrack->GetPosition() << " step(cm)=" << aTrack->GetStep()->GetStepLength() / CLHEP::cm;
0313 }