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 ( G4String  modelName,
G4Region *  envelope,
const edm::ParameterSet parSet 
)

Definition at line 35 of file GflashHadronShowerModel.cc.

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

38  : G4VFastSimulationModel(modelName, envelope), theParSet(parSet) {
39  theWatcherOn = parSet.getParameter<bool>("watcherOn");
40 
48 
49  theGflashStep = new G4Step();
50  theGflashTouchableHandle = new G4TouchableHistory();
51  theGflashNavigator = new G4Navigator();
52 }
T getParameter(std::string const &) const
static GflashHistogram * instance()
GflashKaonMinusShowerProfile * theKaonMinusProfile
GflashAntiProtonShowerProfile * theAntiProtonProfile
GflashHadronShowerProfile * theProfile
G4TouchableHandle theGflashTouchableHandle
GflashKaonPlusShowerProfile * theKaonPlusProfile
GflashPiKShowerProfile * thePiKProfile
GflashProtonShowerProfile * theProtonProfile
GflashHadronShowerModel::~GflashHadronShowerModel ( )
override

Definition at line 54 of file GflashHadronShowerModel.cc.

References theGflashStep, and theProfile.

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

Member Function Documentation

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

Definition at line 97 of file GflashHadronShowerModel.cc.

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

97  {
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 }
const double GeV
Definition: MathUtil.h:16
GflashKaonMinusShowerProfile * theKaonMinusProfile
GflashAntiProtonShowerProfile * theAntiProtonProfile
void makeHits(const G4FastTrack &fastTrack)
GflashHadronShowerProfile * theProfile
int findShowerType(const Gflash3Vector &position)
static int position[264][3]
Definition: ReadPGInfo.cc:509
GflashKaonPlusShowerProfile * theKaonPlusProfile
GflashPiKShowerProfile * thePiKProfile
GflashProtonShowerProfile * theProtonProfile
void initialize(int showerType, double energy, double globalTime, double charge, Gflash3Vector &position, Gflash3Vector &momentum)
G4bool GflashHadronShowerModel::excludeDetectorRegion ( const G4FastTrack &  fastTrack)
private

Definition at line 247 of file GflashHadronShowerModel.cc.

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

Referenced by ModelTrigger().

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

Definition at line 61 of file GflashHadronShowerModel.cc.

61  {
62  return &particleType == G4PionMinus::PionMinusDefinition() || &particleType == G4PionPlus::PionPlusDefinition() ||
63  &particleType == G4KaonMinus::KaonMinusDefinition() || &particleType == G4KaonPlus::KaonPlusDefinition() ||
64  &particleType == G4AntiProton::AntiProtonDefinition() || &particleType == G4Proton::ProtonDefinition();
65 }
G4bool GflashHadronShowerModel::isFirstInelasticInteraction ( const G4FastTrack &  fastTrack)
private

Definition at line 190 of file GflashHadronShowerModel.cc.

References GflashHistogram::deltaStep, GflashHistogram::energyLoss, GflashHistogram::energyRatio, GflashHistogram::getStoreFlag(), GeV, cuy::isFirst, GflashHistogram::kineticEnergy, objects.autophobj::particleType, GflashHistogram::postStepPosition, GflashHistogram::preStepPosition, Gflash::QuasiElasticLike, and theHisto.

Referenced by ModelTrigger().

190  {
191  G4bool isFirst = false;
192 
193  G4StepPoint *preStep = fastTrack.GetPrimaryTrack()->GetStep()->GetPreStepPoint();
194  G4StepPoint *postStep = fastTrack.GetPrimaryTrack()->GetStep()->GetPostStepPoint();
195 
196  G4String procName = postStep->GetProcessDefinedStep()->GetProcessName();
197  G4ParticleDefinition *particleType = fastTrack.GetPrimaryTrack()->GetDefinition();
198 
199  //@@@ this part is still temporary and the cut for the variable ratio should
200  // be optimized later
201 
202  if ((particleType == G4PionPlus::PionPlusDefinition() && procName == "WrappedPionPlusInelastic") ||
203  (particleType == G4PionMinus::PionMinusDefinition() && procName == "WrappedPionMinusInelastic") ||
204  (particleType == G4KaonPlus::KaonPlusDefinition() && procName == "WrappedKaonPlusInelastic") ||
205  (particleType == G4KaonMinus::KaonMinusDefinition() && procName == "WrappedKaonMinusInelastic") ||
206  (particleType == G4AntiProton::AntiProtonDefinition() && procName == "WrappedAntiProtonInelastic") ||
207  (particleType == G4Proton::ProtonDefinition() && procName == "WrappedProtonInelastic")) {
208  // skip to the second interaction if the first inelastic is a quasi-elastic
209  // 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
216  // processes.
217  //@@@may require an additional condition only for hadron interaction with
218  // the process name, 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  // Fill debugging histograms and check information on secondaries -
233  // remove after final implimentation
234 
235  if (theHisto->getStoreFlag()) {
236  theHisto->preStepPosition->Fill(preStep->GetPosition().getRho() / cm);
237  theHisto->postStepPosition->Fill(postStep->GetPosition().getRho() / cm);
238  theHisto->deltaStep->Fill((postStep->GetPosition() - preStep->GetPosition()).getRho() / cm);
239  theHisto->kineticEnergy->Fill(fastTrack.GetPrimaryTrack()->GetKineticEnergy() / GeV);
240  theHisto->energyLoss->Fill(fabs(fastTrack.GetPrimaryTrack()->GetStep()->GetDeltaEnergy() / GeV));
241  theHisto->energyRatio->Fill(leadingEnergy / preStep->GetTotalEnergy());
242  }
243  }
244  return isFirst;
245 }
const double GeV
Definition: MathUtil.h:16
isFirst
Definition: cuy.py:419
const double QuasiElasticLike
void GflashHadronShowerModel::makeHits ( const G4FastTrack &  fastTrack)
private

Definition at line 131 of file GflashHadronShowerModel.cc.

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

Referenced by DoIt().

131  {
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  G4LogicalVolume *lv = aCurrentVolume->GetLogicalVolume();
160  if (lv->GetRegion()->GetName() != "CaloRegion")
161  continue;
162 
163  theGflashStep->GetPreStepPoint()->SetSensitiveDetector(aCurrentVolume->GetLogicalVolume()->GetSensitiveDetector());
164  G4VSensitiveDetector *aSensitive = theGflashStep->GetPreStepPoint()->GetSensitiveDetector();
165 
166  if (aSensitive == nullptr)
167  continue;
168 
169  G4String nameCalor = aCurrentVolume->GetName();
170  nameCalor.assign(nameCalor, 0, 2);
171  G4double samplingWeight = 1.0;
172  if (nameCalor == "HB") {
173  samplingWeight = Gflash::scaleSensitiveHB;
174  } else if (nameCalor == "HE" || nameCalor == "HT") {
175  samplingWeight = Gflash::scaleSensitiveHE;
176  }
177  theGflashStep->SetTotalEnergyDeposit(spotIter->getEnergy() * samplingWeight);
178 
179  aSensitive->Hit(theGflashStep);
180  }
181 }
SimActivityRegistry::G4StepSignal m_g4StepSignal
void updateGflashStep(const G4ThreeVector &position, G4double time)
GflashHadronShowerProfile * theProfile
G4TouchableHandle theGflashTouchableHandle
const double scaleSensitiveHE
std::vector< GflashHit > & getGflashHitList()
const double scaleSensitiveHB
G4bool GflashHadronShowerModel::ModelTrigger ( const G4FastTrack &  fastTrack)
override

Definition at line 67 of file GflashHadronShowerModel.cc.

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

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

Definition at line 183 of file GflashHadronShowerModel.cc.

References theGflashStep, and theGflashTouchableHandle.

Referenced by makeHits().

183  {
184  theGflashStep->GetPostStepPoint()->SetGlobalTime(timeGlobal);
185  theGflashStep->GetPreStepPoint()->SetPosition(spotPosition);
186  theGflashStep->GetPostStepPoint()->SetPosition(spotPosition);
187  theGflashStep->GetPreStepPoint()->SetTouchableHandle(theGflashTouchableHandle);
188 }
G4TouchableHandle theGflashTouchableHandle

Member Data Documentation

GflashAntiProtonShowerProfile* GflashHadronShowerModel::theAntiProtonProfile
private

Definition at line 49 of file GflashHadronShowerModel.h.

Referenced by DoIt(), and GflashHadronShowerModel().

G4Navigator* GflashHadronShowerModel::theGflashNavigator
private

Definition at line 52 of file GflashHadronShowerModel.h.

Referenced by GflashHadronShowerModel(), and makeHits().

G4Step* GflashHadronShowerModel::theGflashStep
private
G4TouchableHandle GflashHadronShowerModel::theGflashTouchableHandle
private

Definition at line 53 of file GflashHadronShowerModel.h.

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

GflashHistogram* GflashHadronShowerModel::theHisto
private
GflashKaonMinusShowerProfile* GflashHadronShowerModel::theKaonMinusProfile
private

Definition at line 47 of file GflashHadronShowerModel.h.

Referenced by DoIt(), and GflashHadronShowerModel().

GflashKaonPlusShowerProfile* GflashHadronShowerModel::theKaonPlusProfile
private

Definition at line 46 of file GflashHadronShowerModel.h.

Referenced by DoIt(), and GflashHadronShowerModel().

edm::ParameterSet GflashHadronShowerModel::theParSet
private

Definition at line 43 of file GflashHadronShowerModel.h.

Referenced by excludeDetectorRegion().

GflashPiKShowerProfile* GflashHadronShowerModel::thePiKProfile
private

Definition at line 45 of file GflashHadronShowerModel.h.

Referenced by DoIt(), and GflashHadronShowerModel().

GflashHadronShowerProfile* GflashHadronShowerModel::theProfile
private
GflashProtonShowerProfile* GflashHadronShowerModel::theProtonProfile
private

Definition at line 48 of file GflashHadronShowerModel.h.

Referenced by DoIt(), and GflashHadronShowerModel().

G4bool GflashHadronShowerModel::theWatcherOn
private

Definition at line 42 of file GflashHadronShowerModel.h.

Referenced by GflashHadronShowerModel(), and makeHits().