CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
GflashHadronShowerModel Class Reference

#include <GflashHadronShowerModel.h>

Inheritance diagram for GflashHadronShowerModel:

Public Member Functions

void DoIt (const G4FastTrack &, G4FastStep &) override
 
 GflashHadronShowerModel (G4String modelName, G4Region *envelope, const edm::ParameterSet &parSet)
 
G4bool IsApplicable (const G4ParticleDefinition &) override
 
G4bool ModelTrigger (const G4FastTrack &) override
 
 ~GflashHadronShowerModel () override
 

Private Member Functions

G4bool excludeDetectorRegion (const G4FastTrack &fastTrack)
 
G4bool isFirstInelasticInteraction (const G4FastTrack &fastTrack)
 
void makeHits (const G4FastTrack &fastTrack)
 
void updateGflashStep (const G4ThreeVector &position, G4double time)
 

Private Attributes

GflashAntiProtonShowerProfiletheAntiProtonProfile
 
G4Navigator * theGflashNavigator
 
G4Step * theGflashStep
 
G4TouchableHandle theGflashTouchableHandle
 
GflashHistogramtheHisto
 
GflashKaonMinusShowerProfiletheKaonMinusProfile
 
GflashKaonPlusShowerProfiletheKaonPlusProfile
 
edm::ParameterSet theParSet
 
GflashPiKShowerProfilethePiKProfile
 
GflashHadronShowerProfiletheProfile
 
GflashProtonShowerProfiletheProtonProfile
 
G4bool theWatcherOn
 

Detailed Description

Definition at line 19 of file GflashHadronShowerModel.h.

Constructor & Destructor Documentation

◆ GflashHadronShowerModel()

GflashHadronShowerModel::GflashHadronShowerModel ( G4String  modelName,
G4Region *  envelope,
const edm::ParameterSet parSet 
)

Definition at line 34 of file GflashHadronShowerModel.cc.

References edm::ParameterSet::getParameter(), GflashHistogram::instance(), theAntiProtonProfile, theGflashNavigator, theGflashStep, theGflashTouchableHandle, theHisto, theKaonMinusProfile, theKaonPlusProfile, thePiKProfile, theProfile, theProtonProfile, and theWatcherOn.

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 }
static GflashHistogram * instance()
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
GflashKaonMinusShowerProfile * theKaonMinusProfile
GflashAntiProtonShowerProfile * theAntiProtonProfile
GflashHadronShowerProfile * theProfile
G4TouchableHandle theGflashTouchableHandle
GflashKaonPlusShowerProfile * theKaonPlusProfile
GflashPiKShowerProfile * thePiKProfile
GflashProtonShowerProfile * theProtonProfile

◆ ~GflashHadronShowerModel()

GflashHadronShowerModel::~GflashHadronShowerModel ( )
override

Definition at line 53 of file GflashHadronShowerModel.cc.

References theGflashStep, and theProfile.

53  {
54  if (theProfile)
55  delete theProfile;
56  if (theGflashStep)
57  delete theGflashStep;
58 }
GflashHadronShowerProfile * theProfile

Member Function Documentation

◆ DoIt()

void GflashHadronShowerModel::DoIt ( const G4FastTrack &  fastTrack,
G4FastStep &  fastStep 
)
override

Definition at line 96 of file GflashHadronShowerModel.cc.

References ALCARECOTkAlJpsiMuMu_cff::charge, hcalRecHitTable_cff::energy, Gflash::findShowerType(), GflashHadronShowerProfile::hadronicParameterization(), GflashHadronShowerProfile::initialize(), GflashHadronShowerProfile::loadParameters(), makeHits(), PbPb_ZMuSkimMuonDPG_cff::particleType, position, theAntiProtonProfile, theKaonMinusProfile, theKaonPlusProfile, thePiKProfile, theProfile, and theProtonProfile.

96  {
97  // kill the particle
98  fastStep.KillPrimaryTrack();
99  fastStep.ProposePrimaryTrackPathLength(0.0);
100 
101  // parameterize energy depostion by the particle type
102  G4ParticleDefinition *particleType = fastTrack.GetPrimaryTrack()->GetDefinition();
103 
105  if (particleType == G4KaonMinus::KaonMinusDefinition())
107  else if (particleType == G4KaonPlus::KaonPlusDefinition())
109  else if (particleType == G4AntiProton::AntiProtonDefinition())
111  else if (particleType == G4Proton::ProtonDefinition())
113 
114  // input variables for GflashHadronShowerProfile
115  G4double energy = fastTrack.GetPrimaryTrack()->GetKineticEnergy() / GeV;
116  G4double globalTime = fastTrack.GetPrimaryTrack()->GetStep()->GetPostStepPoint()->GetGlobalTime();
117  G4double charge = fastTrack.GetPrimaryTrack()->GetStep()->GetPreStepPoint()->GetCharge();
118  G4ThreeVector position = fastTrack.GetPrimaryTrack()->GetPosition() / cm;
119  G4ThreeVector momentum = fastTrack.GetPrimaryTrack()->GetMomentum() / GeV;
120  G4int showerType = Gflash::findShowerType(position);
121 
122  theProfile->initialize(showerType, energy, globalTime, charge, position, momentum);
125 
126  // make hits
127  makeHits(fastTrack);
128 }
GflashKaonMinusShowerProfile * theKaonMinusProfile
GflashAntiProtonShowerProfile * theAntiProtonProfile
void makeHits(const G4FastTrack &fastTrack)
GflashHadronShowerProfile * theProfile
int findShowerType(const Gflash3Vector &position)
static int position[264][3]
Definition: ReadPGInfo.cc:289
GflashKaonPlusShowerProfile * theKaonPlusProfile
GflashPiKShowerProfile * thePiKProfile
GflashProtonShowerProfile * theProtonProfile
void initialize(int showerType, double energy, double globalTime, double charge, Gflash3Vector &position, Gflash3Vector &momentum)

◆ excludeDetectorRegion()

G4bool GflashHadronShowerModel::excludeDetectorRegion ( const G4FastTrack &  fastTrack)
private

Definition at line 239 of file GflashHadronShowerModel.cc.

References PVValHelper::eta, Gflash::EtaMin, ecalTB2006H4_GenSimDigiReco_cfg::G4cout, Gflash::getCalorimeterNumber(), edm::ParameterSet::getUntrackedParameter(), Gflash::kENCA, Gflash::kHB, Gflash::kHE, Gflash::MinDistanceToOut, Gflash::Rmax, theParSet, verbosity, Gflash::Zmax, and Gflash::Zmin.

Referenced by ModelTrigger().

239  {
240  const double invcm = 1.0 / CLHEP::cm;
241  G4bool isExcluded = false;
242  int verbosity = theParSet.getUntrackedParameter<int>("Verbosity");
243 
244  // exclude regions where geometry are complicated
245  //+- one supermodule around the EB/EE boundary: 1.479 +- 0.0174*5
246  G4double eta = fastTrack.GetPrimaryTrack()->GetPosition().pseudoRapidity();
247  if (std::fabs(eta) > 1.392 && std::fabs(eta) < 1.566) {
248  if (verbosity > 0) {
249  edm::LogVerbatim("SimGeneralGFlash") << "GflashHadronShowerModel: excluding region of eta = " << eta;
250  }
251  return true;
252  } else {
253  const G4StepPoint *postStep = fastTrack.GetPrimaryTrack()->GetStep()->GetPostStepPoint();
254 
255  Gflash::CalorimeterNumber kCalor = Gflash::getCalorimeterNumber(postStep->GetPosition() * invcm);
256  G4double distOut = 9999.0;
257 
258  // exclude the region where the shower starting point is inside the
259  // preshower
260  if (std::fabs(eta) > Gflash::EtaMin[Gflash::kENCA] &&
261  std::fabs((postStep->GetPosition()).getZ() / CLHEP::cm) < Gflash::Zmin[Gflash::kENCA]) {
262  return true;
263  }
264 
265  //<---the shower starting point is always inside envelopes
266  //@@@exclude the region where the shower starting point is too close to the
267  // end of the hadronic envelopes (may need to be optimized further!)
268  //@@@if we extend parameterization including Magnet/HO, we need to change
269  // this strategy
270  if (kCalor == Gflash::kHB) {
271  distOut = Gflash::Rmax[Gflash::kHB] - postStep->GetPosition().getRho() * invcm;
272  if (distOut < Gflash::MinDistanceToOut)
273  isExcluded = true;
274  } else if (kCalor == Gflash::kHE) {
275  distOut = Gflash::Zmax[Gflash::kHE] - std::fabs(postStep->GetPosition().getZ() * invcm);
276  if (distOut < Gflash::MinDistanceToOut)
277  isExcluded = true;
278  }
279 
280  //@@@remove this print statement later
281  if (isExcluded && verbosity > 0) {
282  G4cout << "GflashHadronShowerModel: skipping kCalor = " << kCalor << " DistanceToOut " << distOut << " from ("
283  << (postStep->GetPosition()).getRho() * invcm << ":" << (postStep->GetPosition()).getZ() * invcm
284  << ") of KE = " << fastTrack.GetPrimaryTrack()->GetKineticEnergy() / CLHEP::GeV << G4endl;
285  }
286  }
287 
288  return isExcluded;
289 }
Log< level::Info, true > LogVerbatim
const double MinDistanceToOut
const double Zmax[kNumberCalorimeter]
T getUntrackedParameter(std::string const &, T const &) const
CalorimeterNumber getCalorimeterNumber(const Gflash3Vector &position)
const double EtaMin[kNumberCalorimeter]
const int verbosity
const double Zmin[kNumberCalorimeter]
const double Rmax[kNumberCalorimeter]

◆ IsApplicable()

G4bool GflashHadronShowerModel::IsApplicable ( const G4ParticleDefinition &  particleType)
override

Definition at line 60 of file GflashHadronShowerModel.cc.

References PbPb_ZMuSkimMuonDPG_cff::particleType.

60  {
61  return &particleType == G4PionMinus::PionMinusDefinition() || &particleType == G4PionPlus::PionPlusDefinition() ||
62  &particleType == G4KaonMinus::KaonMinusDefinition() || &particleType == G4KaonPlus::KaonPlusDefinition() ||
63  &particleType == G4AntiProton::AntiProtonDefinition() || &particleType == G4Proton::ProtonDefinition();
64 }

◆ isFirstInelasticInteraction()

G4bool GflashHadronShowerModel::isFirstInelasticInteraction ( const G4FastTrack &  fastTrack)
private

Definition at line 182 of file GflashHadronShowerModel.cc.

References GflashHistogram::deltaStep, GflashHistogram::energyLoss, GflashHistogram::energyRatio, GflashHistogram::getStoreFlag(), l1ctLayer1_patternWriters_cff::isec, cuy::isFirst, GflashHistogram::kineticEnergy, PbPb_ZMuSkimMuonDPG_cff::particleType, GflashHistogram::postStepPosition, GflashHistogram::preStepPosition, Gflash::QuasiElasticLike, and theHisto.

Referenced by ModelTrigger().

182  {
183  G4bool isFirst = false;
184 
185  G4StepPoint *preStep = fastTrack.GetPrimaryTrack()->GetStep()->GetPreStepPoint();
186  G4StepPoint *postStep = fastTrack.GetPrimaryTrack()->GetStep()->GetPostStepPoint();
187 
188  G4String procName = postStep->GetProcessDefinedStep()->GetProcessName();
189  G4ParticleDefinition *particleType = fastTrack.GetPrimaryTrack()->GetDefinition();
190 
191  //@@@ this part is still temporary and the cut for the variable ratio should
192  // be optimized later
193 
194  if ((particleType == G4PionPlus::PionPlusDefinition() && procName == "WrappedPionPlusInelastic") ||
195  (particleType == G4PionMinus::PionMinusDefinition() && procName == "WrappedPionMinusInelastic") ||
196  (particleType == G4KaonPlus::KaonPlusDefinition() && procName == "WrappedKaonPlusInelastic") ||
197  (particleType == G4KaonMinus::KaonMinusDefinition() && procName == "WrappedKaonMinusInelastic") ||
198  (particleType == G4AntiProton::AntiProtonDefinition() && procName == "WrappedAntiProtonInelastic") ||
199  (particleType == G4Proton::ProtonDefinition() && procName == "WrappedProtonInelastic")) {
200  // skip to the second interaction if the first inelastic is a quasi-elastic
201  // like interaction
202  //@@@ the cut may be optimized later
203 
204  const G4TrackVector *fSecondaryVector = fastTrack.GetPrimaryTrack()->GetStep()->GetSecondary();
205  G4double leadingEnergy = 0.0;
206 
207  // loop over 'all' secondaries including those produced by continuous
208  // processes.
209  //@@@may require an additional condition only for hadron interaction with
210  // the process name, but it will not change the result anyway
211 
212  for (unsigned int isec = 0; isec < fSecondaryVector->size(); isec++) {
213  G4Track *fSecondaryTrack = (*fSecondaryVector)[isec];
214  G4double secondaryEnergy = fSecondaryTrack->GetKineticEnergy();
215 
216  if (secondaryEnergy > leadingEnergy) {
217  leadingEnergy = secondaryEnergy;
218  }
219  }
220 
221  if ((preStep->GetTotalEnergy() != 0) && (leadingEnergy / preStep->GetTotalEnergy() < Gflash::QuasiElasticLike))
222  isFirst = true;
223 
224  // Fill debugging histograms and check information on secondaries -
225  // remove after final implimentation
226 
227  if (theHisto->getStoreFlag()) {
228  theHisto->preStepPosition->Fill(preStep->GetPosition().getRho() / cm);
229  theHisto->postStepPosition->Fill(postStep->GetPosition().getRho() / cm);
230  theHisto->deltaStep->Fill((postStep->GetPosition() - preStep->GetPosition()).getRho() / cm);
231  theHisto->kineticEnergy->Fill(fastTrack.GetPrimaryTrack()->GetKineticEnergy() / GeV);
232  theHisto->energyLoss->Fill(fabs(fastTrack.GetPrimaryTrack()->GetStep()->GetDeltaEnergy() / GeV));
233  theHisto->energyRatio->Fill(leadingEnergy / preStep->GetTotalEnergy());
234  }
235  }
236  return isFirst;
237 }
isFirst
Definition: cuy.py:418
const double QuasiElasticLike

◆ makeHits()

void GflashHadronShowerModel::makeHits ( const G4FastTrack &  fastTrack)
private

Definition at line 130 of file GflashHadronShowerModel.cc.

References GflashHadronShowerProfile::getGflashHitList(), Gflash::scaleSensitiveHB, Gflash::scaleSensitiveHE, theGflashNavigator, theGflashStep, theGflashTouchableHandle, theProfile, and updateGflashStep().

Referenced by DoIt().

130  {
131  std::vector<GflashHit> &gflashHitList = theProfile->getGflashHitList();
132 
133  theGflashStep->SetTrack(const_cast<G4Track *>(fastTrack.GetPrimaryTrack()));
134  theGflashStep->GetPostStepPoint()->SetProcessDefinedStep(
135  const_cast<G4VProcess *>(fastTrack.GetPrimaryTrack()->GetStep()->GetPostStepPoint()->GetProcessDefinedStep()));
136  theGflashNavigator->SetWorldVolume(
137  G4TransportationManager::GetTransportationManager()->GetNavigatorForTracking()->GetWorldVolume());
138 
139  std::vector<GflashHit>::const_iterator spotIter = gflashHitList.begin();
140  std::vector<GflashHit>::const_iterator spotIterEnd = gflashHitList.end();
141 
142  for (; spotIter != spotIterEnd; spotIter++) {
143  theGflashNavigator->LocateGlobalPointAndUpdateTouchableHandle(
144  spotIter->getPosition(), G4ThreeVector(0, 0, 0), theGflashTouchableHandle, false);
145  updateGflashStep(spotIter->getPosition(), spotIter->getTime());
146 
147  const G4VPhysicalVolume *aCurrentVolume = theGflashStep->GetPreStepPoint()->GetPhysicalVolume();
148  if (aCurrentVolume == nullptr)
149  continue;
150 
151  const G4LogicalVolume *lv = aCurrentVolume->GetLogicalVolume();
152  if (lv->GetRegion()->GetName() != "CaloRegion")
153  continue;
154 
155  theGflashStep->GetPreStepPoint()->SetSensitiveDetector(aCurrentVolume->GetLogicalVolume()->GetSensitiveDetector());
156  G4VSensitiveDetector *aSensitive = theGflashStep->GetPreStepPoint()->GetSensitiveDetector();
157 
158  if (aSensitive == nullptr)
159  continue;
160 
161  G4String nameCalor = aCurrentVolume->GetName();
162  nameCalor.assign(nameCalor, 0, 2);
163  G4double samplingWeight = 1.0;
164  if (nameCalor == "HB") {
165  samplingWeight = Gflash::scaleSensitiveHB;
166  } else if (nameCalor == "HE" || nameCalor == "HT") {
167  samplingWeight = Gflash::scaleSensitiveHE;
168  }
169  theGflashStep->SetTotalEnergyDeposit(spotIter->getEnergy() * samplingWeight);
170 
171  aSensitive->Hit(theGflashStep);
172  }
173 }
void updateGflashStep(const G4ThreeVector &position, G4double time)
GflashHadronShowerProfile * theProfile
G4TouchableHandle theGflashTouchableHandle
const double scaleSensitiveHE
std::vector< GflashHit > & getGflashHitList()
const double scaleSensitiveHB

◆ ModelTrigger()

G4bool GflashHadronShowerModel::ModelTrigger ( const G4FastTrack &  fastTrack)
override

Definition at line 66 of file GflashHadronShowerModel.cc.

References Gflash::energyCutOff, excludeDetectorRegion(), and isFirstInelasticInteraction().

66  {
67  // ModelTrigger returns false for Gflash Hadronic Shower Model if it is not
68  // tested from the corresponding wrapper process, GflashHadronWrapperProcess.
69  // Temporarily the track status is set to fPostponeToNextEvent at the wrapper
70  // process before ModelTrigger is really tested for the second time through
71  // PostStepGPIL of enviaged parameterization processes. The better
72  // implmentation may be using via G4VUserTrackInformation of each track, which
73  // requires to modify a geant source code of stepping (G4SteppingManager2)
74 
75  G4bool trigger = false;
76 
77  // mininum energy cutoff to parameterize
78  if (fastTrack.GetPrimaryTrack()->GetKineticEnergy() < Gflash::energyCutOff * CLHEP::GeV)
79  return trigger;
80 
81  // check whether this is called from the normal GPIL or the wrapper process
82  // GPIL
83  if (fastTrack.GetPrimaryTrack()->GetTrackStatus() == fPostponeToNextEvent) {
84  // Shower pameterization start at the first inelastic interaction point
85  G4bool isInelastic = isFirstInelasticInteraction(fastTrack);
86 
87  // Other conditions
88  if (isInelastic) {
89  trigger = (!excludeDetectorRegion(fastTrack));
90  }
91  }
92 
93  return trigger;
94 }
const double energyCutOff
G4bool isFirstInelasticInteraction(const G4FastTrack &fastTrack)
G4bool excludeDetectorRegion(const G4FastTrack &fastTrack)

◆ updateGflashStep()

void GflashHadronShowerModel::updateGflashStep ( const G4ThreeVector &  position,
G4double  time 
)
private

Definition at line 175 of file GflashHadronShowerModel.cc.

References theGflashStep, and theGflashTouchableHandle.

Referenced by makeHits().

175  {
176  theGflashStep->GetPostStepPoint()->SetGlobalTime(timeGlobal);
177  theGflashStep->GetPreStepPoint()->SetPosition(spotPosition);
178  theGflashStep->GetPostStepPoint()->SetPosition(spotPosition);
179  theGflashStep->GetPreStepPoint()->SetTouchableHandle(theGflashTouchableHandle);
180 }
G4TouchableHandle theGflashTouchableHandle

Member Data Documentation

◆ theAntiProtonProfile

GflashAntiProtonShowerProfile* GflashHadronShowerModel::theAntiProtonProfile
private

Definition at line 49 of file GflashHadronShowerModel.h.

Referenced by DoIt(), and GflashHadronShowerModel().

◆ theGflashNavigator

G4Navigator* GflashHadronShowerModel::theGflashNavigator
private

Definition at line 52 of file GflashHadronShowerModel.h.

Referenced by GflashHadronShowerModel(), and makeHits().

◆ theGflashStep

G4Step* GflashHadronShowerModel::theGflashStep
private

◆ theGflashTouchableHandle

G4TouchableHandle GflashHadronShowerModel::theGflashTouchableHandle
private

Definition at line 53 of file GflashHadronShowerModel.h.

Referenced by GflashHadronShowerModel(), makeHits(), and updateGflashStep().

◆ theHisto

GflashHistogram* GflashHadronShowerModel::theHisto
private

◆ theKaonMinusProfile

GflashKaonMinusShowerProfile* GflashHadronShowerModel::theKaonMinusProfile
private

Definition at line 47 of file GflashHadronShowerModel.h.

Referenced by DoIt(), and GflashHadronShowerModel().

◆ theKaonPlusProfile

GflashKaonPlusShowerProfile* GflashHadronShowerModel::theKaonPlusProfile
private

Definition at line 46 of file GflashHadronShowerModel.h.

Referenced by DoIt(), and GflashHadronShowerModel().

◆ theParSet

edm::ParameterSet GflashHadronShowerModel::theParSet
private

Definition at line 43 of file GflashHadronShowerModel.h.

Referenced by excludeDetectorRegion().

◆ thePiKProfile

GflashPiKShowerProfile* GflashHadronShowerModel::thePiKProfile
private

Definition at line 45 of file GflashHadronShowerModel.h.

Referenced by DoIt(), and GflashHadronShowerModel().

◆ theProfile

GflashHadronShowerProfile* GflashHadronShowerModel::theProfile
private

◆ theProtonProfile

GflashProtonShowerProfile* GflashHadronShowerModel::theProtonProfile
private

Definition at line 48 of file GflashHadronShowerModel.h.

Referenced by DoIt(), and GflashHadronShowerModel().

◆ theWatcherOn

G4bool GflashHadronShowerModel::theWatcherOn
private

Definition at line 42 of file GflashHadronShowerModel.h.

Referenced by GflashHadronShowerModel().