CMS 3D CMS Logo

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