CMS 3D CMS Logo

ElectronLimiter.cc
Go to the documentation of this file.
1 //
2 // V.Ivanchenko 2013/10/19
3 // step limiter and killer for e+,e- and other charged particles
4 //
8 
9 #include "G4ParticleDefinition.hh"
10 #include "G4VEnergyLossProcess.hh"
11 #include "G4LossTableManager.hh"
12 #include "G4DummyModel.hh"
13 #include "G4Step.hh"
14 #include "G4Track.hh"
15 #include "G4Region.hh"
16 #include <CLHEP/Units/SystemOfUnits.h>
17 #include "G4TransportationProcessType.hh"
18 
19 ElectronLimiter::ElectronLimiter(const edm::ParameterSet &p, const G4ParticleDefinition *part)
20  : G4VEmProcess("eLimiter", fGeneral),
21  ionisation_(nullptr),
22  particle_(part),
23  nRegions_(0),
24  rangeCheckFlag_(false),
25  fieldCheckFlag_(false),
26  killTrack_(false) {
27  // set Process Sub Type
28  SetProcessSubType(static_cast<int>(STEP_LIMITER));
29  minStepLimit_ = p.getParameter<double>("MinStepLimit") * CLHEP::mm;
31 }
32 
34 
35 void ElectronLimiter::InitialiseProcess(const G4ParticleDefinition *) {
36  G4LossTableManager *lManager = G4LossTableManager::Instance();
37  if (rangeCheckFlag_) {
38  ionisation_ = lManager->GetEnergyLossProcess(particle_);
39  if (!ionisation_) {
40  rangeCheckFlag_ = false;
41  }
42  }
43  AddEmModel(0, new G4DummyModel());
44 
45  if (lManager->IsMaster()) {
46  edm::LogVerbatim("SimG4CoreApplication")
47  << "ElectronLimiter::BuildPhysicsTable for " << particle_->GetParticleName() << " ioni: " << ionisation_
48  << " rangeCheckFlag: " << rangeCheckFlag_ << " fieldCheckFlag: " << fieldCheckFlag_ << " " << nRegions_
49  << " regions for tracking cuts\n";
50  }
51 }
52 
54  G4double previousLimit,
55  G4ForceCondition *condition) {
56  *condition = NotForced;
57  G4double limit = DBL_MAX;
58  killTrack_ = false;
59 
60  G4double kinEnergy = aTrack.GetKineticEnergy();
61  if (0 < nRegions_) {
62  if (regions_[0] == nullptr) {
63  if (kinEnergy < energyLimits_[0]) {
64  killTrack_ = true;
66  limit = 0.0;
67  }
68  } else {
69  const G4Region *reg = aTrack.GetVolume()->GetLogicalVolume()->GetRegion();
70  for (G4int i = 0; i < nRegions_; ++i) {
71  if (reg == regions_[i] && kinEnergy < energyLimits_[i]) {
72  killTrack_ = true;
74  limit = 0.0;
75  break;
76  }
77  }
78  }
79  }
80  if (!killTrack_ && rangeCheckFlag_) {
81  G4double safety = aTrack.GetStep()->GetPreStepPoint()->GetSafety();
82  if (safety > std::min(minStepLimit_, previousLimit)) {
83  G4double range = ionisation_->GetRange(kinEnergy, aTrack.GetMaterialCutsCouple());
84  if (safety >= range) {
85  killTrack_ = true;
86  limit = 0.0;
87  }
88  }
89  }
90  if (!killTrack_ && fieldCheckFlag_) {
92  }
93  return limit;
94 }
95 
96 G4VParticleChange *ElectronLimiter::PostStepDoIt(const G4Track &aTrack, const G4Step &) {
97  fParticleChange.Initialize(aTrack);
98  if (killTrack_) {
99  fParticleChange.ProposeTrackStatus(fStopAndKill);
100  G4double edep = trcut_->SampleEnergyDepositEcal(aTrack.GetKineticEnergy());
101  fParticleChange.ProposeLocalEnergyDeposit(edep);
102  fParticleChange.SetProposedKineticEnergy(0.0);
103  }
104  return &fParticleChange;
105 }
106 
107 G4bool ElectronLimiter::IsApplicable(const G4ParticleDefinition &) { return true; }
108 
110 
111 void ElectronLimiter::SetTrackingCutPerRegion(std::vector<const G4Region *> &reg,
112  std::vector<G4double> &cut,
113  std::vector<G4double> &fac,
114  std::vector<G4double> &rms) {
115  nRegions_ = reg.size();
116  if (nRegions_ > 0) {
117  regions_ = reg;
118  energyLimits_ = cut;
119  factors_ = fac;
120  rms_ = rms;
121  }
122 }
Log< level::Info, true > LogVerbatim
ElectronLimiter(const edm::ParameterSet &, const G4ParticleDefinition *)
void SetTrackingCutPerRegion(std::vector< const G4Region *> &, std::vector< G4double > &, std::vector< G4double > &, std::vector< G4double > &)
const G4ParticleDefinition * particle_
void StartTracking(G4Track *) override
CMSTrackingCutModel * trcut_
std::vector< G4double > rms_
G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
G4double minStepLimit_
void InitialiseProcess(const G4ParticleDefinition *) override
~ElectronLimiter() override
void InitialiseForStep(G4double fac, G4double rms)
G4VEnergyLossProcess * ionisation_
G4bool IsApplicable(const G4ParticleDefinition &) override
virtual G4double SampleEnergyDepositEcal(G4double kinEnergy)
std::vector< const G4Region * > regions_
part
Definition: HCALResponse.h:20
std::vector< G4double > energyLimits_
G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &) override
std::vector< G4double > factors_