CMS 3D CMS Logo

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