CMS 3D CMS Logo

GflashHadronShowerModel.cc
Go to the documentation of this file.
2 
4 
14 
15 #include "G4EventManager.hh"
16 #include "G4FastSimulationManager.hh"
17 #include "G4TouchableHandle.hh"
18 #include "G4TransportationManager.hh"
19 #include "G4VPhysicalVolume.hh"
20 #include "G4VSensitiveDetector.hh"
21 
22 #include "G4AntiProton.hh"
23 #include "G4KaonMinus.hh"
24 #include "G4KaonPlus.hh"
25 #include "G4PionMinus.hh"
26 #include "G4PionPlus.hh"
27 #include "G4Proton.hh"
28 #include "G4VProcess.hh"
29 
30 #include <vector>
31 
32 using namespace CLHEP;
33 
35  G4Region *envelope,
36  const edm::ParameterSet &parSet)
37  : G4VFastSimulationModel(modelName, envelope), theParSet(parSet) {
38  theWatcherOn = parSet.getParameter<bool>("watcherOn");
39 
47 
48  theGflashStep = new G4Step();
49  theGflashTouchableHandle = new G4TouchableHistory();
50  theGflashNavigator = new G4Navigator();
51 }
52 
54  delete theProfile;
55  delete theGflashStep;
56 }
57 
58 G4bool GflashHadronShowerModel::IsApplicable(const G4ParticleDefinition &particleType) {
59  return &particleType == G4PionMinus::PionMinusDefinition() || &particleType == G4PionPlus::PionPlusDefinition() ||
60  &particleType == G4KaonMinus::KaonMinusDefinition() || &particleType == G4KaonPlus::KaonPlusDefinition() ||
61  &particleType == G4AntiProton::AntiProtonDefinition() || &particleType == G4Proton::ProtonDefinition();
62 }
63 
64 G4bool GflashHadronShowerModel::ModelTrigger(const G4FastTrack &fastTrack) {
65  // ModelTrigger returns false for Gflash Hadronic Shower Model if it is not
66  // tested from the corresponding wrapper process, GflashHadronWrapperProcess.
67  // Temporarily the track status is set to fPostponeToNextEvent at the wrapper
68  // process before ModelTrigger is really tested for the second time through
69  // PostStepGPIL of enviaged parameterization processes. The better
70  // implmentation may be using via G4VUserTrackInformation of each track, which
71  // requires to modify a geant source code of stepping (G4SteppingManager2)
72 
73  G4bool trigger = false;
74 
75  // mininum energy cutoff to parameterize
76  if (fastTrack.GetPrimaryTrack()->GetKineticEnergy() < Gflash::energyCutOff * CLHEP::GeV)
77  return trigger;
78 
79  // check whether this is called from the normal GPIL or the wrapper process
80  // GPIL
81  if (fastTrack.GetPrimaryTrack()->GetTrackStatus() == fPostponeToNextEvent) {
82  // Shower pameterization start at the first inelastic interaction point
83  G4bool isInelastic = isFirstInelasticInteraction(fastTrack);
84 
85  // Other conditions
86  if (isInelastic) {
87  trigger = (!excludeDetectorRegion(fastTrack));
88  }
89  }
90 
91  return trigger;
92 }
93 
94 void GflashHadronShowerModel::DoIt(const G4FastTrack &fastTrack, G4FastStep &fastStep) {
95  // kill the particle
96  fastStep.KillPrimaryTrack();
97  fastStep.ProposePrimaryTrackPathLength(0.0);
98 
99  // parameterize energy depostion by the particle type
100  G4ParticleDefinition *particleType = fastTrack.GetPrimaryTrack()->GetDefinition();
101 
103  if (particleType == G4KaonMinus::KaonMinusDefinition())
105  else if (particleType == G4KaonPlus::KaonPlusDefinition())
107  else if (particleType == G4AntiProton::AntiProtonDefinition())
109  else if (particleType == G4Proton::ProtonDefinition())
111 
112  // input variables for GflashHadronShowerProfile
113  G4double energy = fastTrack.GetPrimaryTrack()->GetKineticEnergy() / GeV;
114  G4double globalTime = fastTrack.GetPrimaryTrack()->GetStep()->GetPostStepPoint()->GetGlobalTime();
115  G4double charge = fastTrack.GetPrimaryTrack()->GetStep()->GetPreStepPoint()->GetCharge();
116  G4ThreeVector position = fastTrack.GetPrimaryTrack()->GetPosition() / cm;
117  G4ThreeVector momentum = fastTrack.GetPrimaryTrack()->GetMomentum() / GeV;
118  G4int showerType = Gflash::findShowerType(position);
119 
120  theProfile->initialize(showerType, energy, globalTime, charge, position, momentum);
123 
124  // make hits
125  makeHits(fastTrack);
126 }
127 
128 void GflashHadronShowerModel::makeHits(const G4FastTrack &fastTrack) {
129  std::vector<GflashHit> &gflashHitList = theProfile->getGflashHitList();
130 
131  theGflashStep->SetTrack(const_cast<G4Track *>(fastTrack.GetPrimaryTrack()));
132  theGflashStep->GetPostStepPoint()->SetProcessDefinedStep(
133  const_cast<G4VProcess *>(fastTrack.GetPrimaryTrack()->GetStep()->GetPostStepPoint()->GetProcessDefinedStep()));
134  theGflashNavigator->SetWorldVolume(
135  G4TransportationManager::GetTransportationManager()->GetNavigatorForTracking()->GetWorldVolume());
136 
137  for (auto &hit : gflashHitList) {
138  theGflashNavigator->LocateGlobalPointAndUpdateTouchableHandle(
139  hit.getPosition(), G4ThreeVector(0, 0, 0), theGflashTouchableHandle, false);
140  updateGflashStep(hit.getPosition(), hit.getTime());
141 
142  const G4VPhysicalVolume *aCurrentVolume = theGflashStep->GetPreStepPoint()->GetPhysicalVolume();
143  if (aCurrentVolume == nullptr)
144  continue;
145 
146  const G4LogicalVolume *lv = aCurrentVolume->GetLogicalVolume();
147  if (lv->GetRegion()->GetName() != "CaloRegion")
148  continue;
149 
150  theGflashStep->GetPreStepPoint()->SetSensitiveDetector(aCurrentVolume->GetLogicalVolume()->GetSensitiveDetector());
151  G4VSensitiveDetector *aSensitive = theGflashStep->GetPreStepPoint()->GetSensitiveDetector();
152 
153  if (aSensitive == nullptr)
154  continue;
155 
156  G4String nameCalor = aCurrentVolume->GetName();
157  nameCalor.assign(nameCalor, 0, 2);
158  G4double samplingWeight = 1.0;
159  if (nameCalor == "HB") {
160  samplingWeight = Gflash::scaleSensitiveHB;
161  } else if (nameCalor == "HE" || nameCalor == "HT") {
162  samplingWeight = Gflash::scaleSensitiveHE;
163  }
164  theGflashStep->SetTotalEnergyDeposit(hit.getEnergy() * samplingWeight);
165 
166  aSensitive->Hit(theGflashStep);
167  }
168 }
169 
170 void GflashHadronShowerModel::updateGflashStep(const G4ThreeVector &spotPosition, G4double timeGlobal) {
171  theGflashStep->GetPostStepPoint()->SetGlobalTime(timeGlobal);
172  theGflashStep->GetPreStepPoint()->SetPosition(spotPosition);
173  theGflashStep->GetPostStepPoint()->SetPosition(spotPosition);
174  theGflashStep->GetPreStepPoint()->SetTouchableHandle(theGflashTouchableHandle);
175 }
176 
177 G4bool GflashHadronShowerModel::isFirstInelasticInteraction(const G4FastTrack &fastTrack) {
178  G4bool isFirst = false;
179 
180  G4StepPoint *preStep = fastTrack.GetPrimaryTrack()->GetStep()->GetPreStepPoint();
181  G4StepPoint *postStep = fastTrack.GetPrimaryTrack()->GetStep()->GetPostStepPoint();
182 
183  G4String procName = postStep->GetProcessDefinedStep()->GetProcessName();
184  G4ParticleDefinition *particleType = fastTrack.GetPrimaryTrack()->GetDefinition();
185 
186  //@@@ this part is still temporary and the cut for the variable ratio should
187  // be optimized later
188 
189  if ((particleType == G4PionPlus::PionPlusDefinition() && procName == "WrappedPionPlusInelastic") ||
190  (particleType == G4PionMinus::PionMinusDefinition() && procName == "WrappedPionMinusInelastic") ||
191  (particleType == G4KaonPlus::KaonPlusDefinition() && procName == "WrappedKaonPlusInelastic") ||
192  (particleType == G4KaonMinus::KaonMinusDefinition() && procName == "WrappedKaonMinusInelastic") ||
193  (particleType == G4AntiProton::AntiProtonDefinition() && procName == "WrappedAntiProtonInelastic") ||
194  (particleType == G4Proton::ProtonDefinition() && procName == "WrappedProtonInelastic")) {
195  // skip to the second interaction if the first inelastic is a quasi-elastic
196  // like interaction
197  //@@@ the cut may be optimized later
198 
199  const G4TrackVector *fSecondaryVector = fastTrack.GetPrimaryTrack()->GetStep()->GetSecondary();
200  G4double leadingEnergy = 0.0;
201 
202  // loop over 'all' secondaries including those produced by continuous
203  // processes.
204  //@@@may require an additional condition only for hadron interaction with
205  // the process name, but it will not change the result anyway
206 
207  for (unsigned int isec = 0; isec < fSecondaryVector->size(); isec++) {
208  G4Track *fSecondaryTrack = (*fSecondaryVector)[isec];
209  G4double secondaryEnergy = fSecondaryTrack->GetKineticEnergy();
210 
211  if (secondaryEnergy > leadingEnergy) {
212  leadingEnergy = secondaryEnergy;
213  }
214  }
215 
216  if ((preStep->GetTotalEnergy() != 0) && (leadingEnergy / preStep->GetTotalEnergy() < Gflash::QuasiElasticLike))
217  isFirst = true;
218 
219  // Fill debugging histograms and check information on secondaries -
220  // remove after final implimentation
221 
222  if (theHisto->getStoreFlag()) {
223  theHisto->preStepPosition->Fill(preStep->GetPosition().getRho() / cm);
224  theHisto->postStepPosition->Fill(postStep->GetPosition().getRho() / cm);
225  theHisto->deltaStep->Fill((postStep->GetPosition() - preStep->GetPosition()).getRho() / cm);
226  theHisto->kineticEnergy->Fill(fastTrack.GetPrimaryTrack()->GetKineticEnergy() / GeV);
227  theHisto->energyLoss->Fill(fabs(fastTrack.GetPrimaryTrack()->GetStep()->GetDeltaEnergy() / GeV));
228  theHisto->energyRatio->Fill(leadingEnergy / preStep->GetTotalEnergy());
229  }
230  }
231  return isFirst;
232 }
233 
234 G4bool GflashHadronShowerModel::excludeDetectorRegion(const G4FastTrack &fastTrack) {
235  const double invcm = 1.0 / CLHEP::cm;
236  G4bool isExcluded = false;
237  int verbosity = theParSet.getUntrackedParameter<int>("Verbosity");
238 
239  // exclude regions where geometry are complicated
240  //+- one supermodule around the EB/EE boundary: 1.479 +- 0.0174*5
241  G4double eta = fastTrack.GetPrimaryTrack()->GetPosition().pseudoRapidity();
242  if (std::fabs(eta) > 1.392 && std::fabs(eta) < 1.566) {
243  if (verbosity > 0) {
244  edm::LogVerbatim("SimGeneralGFlash") << "GflashHadronShowerModel: excluding region of eta = " << eta;
245  }
246  return true;
247  } else {
248  const G4StepPoint *postStep = fastTrack.GetPrimaryTrack()->GetStep()->GetPostStepPoint();
249 
250  Gflash::CalorimeterNumber kCalor = Gflash::getCalorimeterNumber(postStep->GetPosition() * invcm);
251  G4double distOut = 9999.0;
252 
253  // exclude the region where the shower starting point is inside the
254  // preshower
255  if (std::fabs(eta) > Gflash::EtaMin[Gflash::kENCA] &&
256  std::fabs((postStep->GetPosition()).getZ() / CLHEP::cm) < Gflash::Zmin[Gflash::kENCA]) {
257  return true;
258  }
259 
260  //<---the shower starting point is always inside envelopes
261  //@@@exclude the region where the shower starting point is too close to the
262  // end of the hadronic envelopes (may need to be optimized further!)
263  //@@@if we extend parameterization including Magnet/HO, we need to change
264  // this strategy
265  if (kCalor == Gflash::kHB) {
266  distOut = Gflash::Rmax[Gflash::kHB] - postStep->GetPosition().getRho() * invcm;
267  if (distOut < Gflash::MinDistanceToOut)
268  isExcluded = true;
269  } else if (kCalor == Gflash::kHE) {
270  distOut = Gflash::Zmax[Gflash::kHE] - std::fabs(postStep->GetPosition().getZ() * invcm);
271  if (distOut < Gflash::MinDistanceToOut)
272  isExcluded = true;
273  }
274 
275  //@@@remove this print statement later
276  if (isExcluded && verbosity > 0) {
277  G4cout << "GflashHadronShowerModel: skipping kCalor = " << kCalor << " DistanceToOut " << distOut << " from ("
278  << (postStep->GetPosition()).getRho() * invcm << ":" << (postStep->GetPosition()).getZ() * invcm
279  << ") of KE = " << fastTrack.GetPrimaryTrack()->GetKineticEnergy() / CLHEP::GeV << G4endl;
280  }
281  }
282 
283  return isExcluded;
284 }
Log< level::Info, true > LogVerbatim
const double MinDistanceToOut
static GflashHistogram * instance()
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
const double Zmax[kNumberCalorimeter]
GflashKaonMinusShowerProfile * theKaonMinusProfile
GflashAntiProtonShowerProfile * theAntiProtonProfile
G4bool ModelTrigger(const G4FastTrack &) override
void updateGflashStep(const G4ThreeVector &position, G4double time)
void makeHits(const G4FastTrack &fastTrack)
G4bool IsApplicable(const G4ParticleDefinition &) override
T getUntrackedParameter(std::string const &, T const &) const
GflashHadronShowerProfile * theProfile
CalorimeterNumber getCalorimeterNumber(const Gflash3Vector &position)
const double energyCutOff
void DoIt(const G4FastTrack &, G4FastStep &) override
const double EtaMin[kNumberCalorimeter]
isFirst
Definition: cuy.py:418
const int verbosity
int findShowerType(const Gflash3Vector &position)
G4TouchableHandle theGflashTouchableHandle
G4bool isFirstInelasticInteraction(const G4FastTrack &fastTrack)
const double Zmin[kNumberCalorimeter]
const double QuasiElasticLike
const double scaleSensitiveHE
G4bool excludeDetectorRegion(const G4FastTrack &fastTrack)
static int position[264][3]
Definition: ReadPGInfo.cc:289
GflashKaonPlusShowerProfile * theKaonPlusProfile
GflashPiKShowerProfile * thePiKProfile
std::vector< GflashHit > & getGflashHitList()
GflashProtonShowerProfile * theProtonProfile
const double Rmax[kNumberCalorimeter]
const double scaleSensitiveHB
void initialize(int showerType, double energy, double globalTime, double charge, Gflash3Vector &position, Gflash3Vector &momentum)
GflashHadronShowerModel(G4String modelName, G4Region *envelope, const edm::ParameterSet &parSet)