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 &)
 
 GFlashHadronShowerModel (G4String modelName, G4Region *envelope, const edm::ParameterSet &parSet)
 
G4bool IsApplicable (const G4ParticleDefinition &)
 
G4bool ModelTrigger (const G4FastTrack &)
 
 ~GFlashHadronShowerModel ()
 

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
 
GflashKaonMinusShowerProfiletheKaonMinusProfile
 
GflashKaonPlusShowerProfiletheKaonPlusProfile
 
edm::ParameterSet theParSet
 
GflashPiKShowerProfilethePiKProfile
 
GflashHadronShowerProfiletheProfile
 
GflashProtonShowerProfiletheProtonProfile
 
const G4Region * theRegion
 
G4bool theWatcherOn
 

Detailed Description

Definition at line 20 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(), theAntiProtonProfile, theGflashNavigator, theGflashStep, theGflashTouchableHandle, theKaonMinusProfile, theKaonPlusProfile, thePiKProfile, theProfile, theProtonProfile, theRegion, and theWatcherOn.

37  : G4VFastSimulationModel(modelName, envelope), theParSet(parSet)
38 {
39  theWatcherOn = parSet.getParameter<bool>("watcherOn");
40 
41  theProfile = 0;
47  //theHisto = GflashHistogram::instance();
48 
49  theRegion = const_cast<const G4Region*>(envelope);
50 
51  theGflashStep = new G4Step();
52  theGflashTouchableHandle = new G4TouchableHistory();
53  theGflashNavigator = new G4Navigator();
54 
55 }
T getParameter(std::string const &) const
GflashPiKShowerProfile * thePiKProfile
GflashKaonMinusShowerProfile * theKaonMinusProfile
G4TouchableHandle theGflashTouchableHandle
GflashKaonPlusShowerProfile * theKaonPlusProfile
GflashAntiProtonShowerProfile * theAntiProtonProfile
GflashHadronShowerProfile * theProfile
GflashProtonShowerProfile * theProtonProfile
GFlashHadronShowerModel::~GFlashHadronShowerModel ( )

Definition at line 57 of file GFlashHadronShowerModel.cc.

References theAntiProtonProfile, theGflashStep, theKaonMinusProfile, theKaonPlusProfile, thePiKProfile, and theProtonProfile.

58 {
59  delete thePiKProfile;
60  delete theKaonPlusProfile;
61  delete theKaonMinusProfile;
62  delete theProtonProfile;
63  delete theAntiProtonProfile;
64  delete theGflashStep;
65 }
GflashPiKShowerProfile * thePiKProfile
GflashKaonMinusShowerProfile * theKaonMinusProfile
GflashKaonPlusShowerProfile * theKaonPlusProfile
GflashAntiProtonShowerProfile * theAntiProtonProfile
GflashProtonShowerProfile * theProtonProfile

Member Function Documentation

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

Definition at line 113 of file GFlashHadronShowerModel.cc.

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

114 {
115  // kill the particle
116  fastStep.KillPrimaryTrack();
117  fastStep.ProposePrimaryTrackPathLength(0.0);
118 
119  // parameterize energy depostion by the particle type
120  G4ParticleDefinition* particleType = fastTrack.GetPrimaryTrack()->GetDefinition();
121 
123  if(particleType == G4KaonMinus::KaonMinusDefinition()) theProfile = theKaonMinusProfile;
124  else if(particleType == G4KaonPlus::KaonPlusDefinition()) theProfile = theKaonPlusProfile;
125  else if(particleType == G4AntiProton::AntiProtonDefinition()) theProfile = theAntiProtonProfile;
126  else if(particleType == G4Proton::ProtonDefinition()) theProfile = theProtonProfile;
127 
128  //input variables for GflashHadronShowerProfile
129  G4double energy = fastTrack.GetPrimaryTrack()->GetKineticEnergy()/GeV;
130  G4double globalTime = fastTrack.GetPrimaryTrack()->GetStep()->GetPostStepPoint()->GetGlobalTime();
131  G4double charge = fastTrack.GetPrimaryTrack()->GetStep()->GetPreStepPoint()->GetCharge();
132  G4ThreeVector position = fastTrack.GetPrimaryTrack()->GetPosition()/cm;
133  G4ThreeVector momentum = fastTrack.GetPrimaryTrack()->GetMomentum()/GeV;
134  G4int showerType = Gflash::findShowerType(position);
135 
136  theProfile->initialize(showerType,energy,globalTime,charge,position,momentum);
139 
140  // make hits
141  makeHits(fastTrack);
142 
143 }
const double GeV
Definition: MathUtil.h:16
GflashPiKShowerProfile * thePiKProfile
GflashKaonMinusShowerProfile * theKaonMinusProfile
GflashKaonPlusShowerProfile * theKaonPlusProfile
GflashAntiProtonShowerProfile * theAntiProtonProfile
GflashHadronShowerProfile * theProfile
void makeHits(const G4FastTrack &fastTrack)
GflashProtonShowerProfile * theProtonProfile
int findShowerType(const Gflash3Vector &position)
static int position[264][3]
Definition: ReadPGInfo.cc:509
void initialize(int showerType, double energy, double globalTime, double charge, Gflash3Vector &position, Gflash3Vector &momentum)
G4bool GFlashHadronShowerModel::excludeDetectorRegion ( const G4FastTrack &  fastTrack)
private

Definition at line 264 of file GFlashHadronShowerModel.cc.

References 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().

265 {
266  G4bool isExcluded=false;
267  int verbosity = theParSet.getUntrackedParameter<int>("Verbosity");
268 
269  //exclude regions where geometry are complicated
270  //+- one supermodule around the EB/EE boundary: 1.479 +- 0.0174*5
271  G4double eta = fastTrack.GetPrimaryTrack()->GetPosition().pseudoRapidity() ;
272  if(std::fabs(eta) > 1.392 && std::fabs(eta) < 1.566) {
273  if(verbosity>0) {
274  edm::LogInfo("SimG4CoreApplication")
275  << "GFlashHadronShowerModel: excluding region of eta = " << eta;
276  }
277  return true;
278  }
279  else {
280  G4StepPoint* postStep = fastTrack.GetPrimaryTrack()->GetStep()->GetPostStepPoint();
281 
282  Gflash::CalorimeterNumber kCalor = Gflash::getCalorimeterNumber(postStep->GetPosition()/cm);
283  G4double distOut = 9999.0;
284 
285  //exclude the region where the shower starting point is inside the preshower
286  if( std::fabs(eta) > Gflash::EtaMin[Gflash::kENCA] &&
287  std::fabs((postStep->GetPosition()).getZ()/cm) < Gflash::Zmin[Gflash::kENCA]) {
288  return true;
289  }
290 
291  //<---the shower starting point is always inside envelopes
292  //@@@exclude the region where the shower starting point is too close to the end of
293  //the hadronic envelopes (may need to be optimized further!)
294  //@@@if we extend parameterization including Magnet/HO, we need to change this strategy
295  if(kCalor == Gflash::kHB) {
296  distOut = Gflash::Rmax[Gflash::kHB] - postStep->GetPosition().getRho()/cm;
297  if (distOut < Gflash::MinDistanceToOut ) isExcluded = true;
298  }
299  else if(kCalor == Gflash::kHE) {
300  distOut = Gflash::Zmax[Gflash::kHE] - std::fabs(postStep->GetPosition().getZ()/cm);
301  if (distOut < Gflash::MinDistanceToOut ) isExcluded = true;
302  }
303 
304  //@@@remove this print statement later
305  if(isExcluded && verbosity > 0) {
306  edm::LogInfo("SimG4CoreApplication")
307  << "GFlashHadronShowerModel: skipping kCalor = " << kCalor <<
308  " DistanceToOut " << distOut << " from (" << (postStep->GetPosition()).getRho()/cm <<
309  ":" << (postStep->GetPosition()).getZ()/cm << ") of KE = "
310  << fastTrack.GetPrimaryTrack()->GetKineticEnergy()/GeV;
311  }
312  }
313 
314  return isExcluded;
315 }
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)

Definition at line 67 of file GFlashHadronShowerModel.cc.

68 {
69  return
70  &particleType == G4PionMinus::PionMinusDefinition() ||
71  &particleType == G4PionPlus::PionPlusDefinition() ||
72  &particleType == G4KaonMinus::KaonMinusDefinition() ||
73  &particleType == G4KaonPlus::KaonPlusDefinition() ||
74  &particleType == G4AntiProton::AntiProtonDefinition() ||
75  &particleType == G4Proton::ProtonDefinition() ;
76 }
G4bool GFlashHadronShowerModel::isFirstInelasticInteraction ( const G4FastTrack &  fastTrack)
private

Definition at line 207 of file GFlashHadronShowerModel.cc.

References cuy::isFirst, objects.autophobj::particleType, and Gflash::QuasiElasticLike.

208 {
209  G4bool isFirst=false;
210 
211  G4StepPoint* preStep = fastTrack.GetPrimaryTrack()->GetStep()->GetPreStepPoint();
212  G4StepPoint* postStep = fastTrack.GetPrimaryTrack()->GetStep()->GetPostStepPoint();
213 
214  G4String procName = postStep->GetProcessDefinedStep()->GetProcessName();
215  G4ParticleDefinition* particleType = fastTrack.GetPrimaryTrack()->GetDefinition();
216 
217  //@@@ this part is still temporary and the cut for the variable ratio should be optimized later
218 
219  if((particleType == G4PionPlus::PionPlusDefinition() && procName == "WrappedPionPlusInelastic") ||
220  (particleType == G4PionMinus::PionMinusDefinition() && procName == "WrappedPionMinusInelastic") ||
221  (particleType == G4KaonPlus::KaonPlusDefinition() && procName == "WrappedKaonPlusInelastic") ||
222  (particleType == G4KaonMinus::KaonMinusDefinition() && procName == "WrappedKaonMinusInelastic") ||
223  (particleType == G4AntiProton::AntiProtonDefinition() && procName == "WrappedAntiProtonInelastic") ||
224  (particleType == G4Proton::ProtonDefinition() && procName == "WrappedProtonInelastic")) {
225 
226  //skip to the second interaction if the first inelastic is a quasi-elastic like interaction
227  //@@@ the cut may be optimized later
228 
229  const G4TrackVector* fSecondaryVector = fastTrack.GetPrimaryTrack()->GetStep()->GetSecondary();
230  G4double leadingEnergy = 0.0;
231 
232  //loop over 'all' secondaries including those produced by continuous processes.
233  //@@@may require an additional condition only for hadron interaction with the process name,
234  //but it will not change the result anyway
235 
236  for (unsigned int isec = 0 ; isec < fSecondaryVector->size() ; isec++) {
237  G4Track* fSecondaryTrack = (*fSecondaryVector)[isec];
238  G4double secondaryEnergy = fSecondaryTrack->GetKineticEnergy();
239 
240  if(secondaryEnergy > leadingEnergy ) {
241  leadingEnergy = secondaryEnergy;
242  }
243  }
244 
245  if((preStep->GetTotalEnergy()!=0) &&
246  (leadingEnergy/preStep->GetTotalEnergy() < Gflash::QuasiElasticLike)) isFirst = true;
247 
248  //Fill debugging histograms and check information on secondaries -
249  //remove after final implimentation
250  /*
251  if(theHisto->getStoreFlag()) {
252  theHisto->preStepPosition->Fill(preStep->GetPosition().getRho()/cm);
253  theHisto->postStepPosition->Fill(postStep->GetPosition().getRho()/cm);
254  theHisto->deltaStep->Fill((postStep->GetPosition() - preStep->GetPosition()).getRho()/cm);
255  theHisto->kineticEnergy->Fill(fastTrack.GetPrimaryTrack()->GetKineticEnergy()/GeV);
256  theHisto->energyLoss->Fill(fabs(fastTrack.GetPrimaryTrack()->GetStep()->GetDeltaEnergy()/GeV));
257  theHisto->energyRatio->Fill(leadingEnergy/preStep->GetTotalEnergy());
258  }
259  */
260  }
261  return isFirst;
262 }
isFirst
Definition: cuy.py:417
const double QuasiElasticLike
void GFlashHadronShowerModel::makeHits ( const G4FastTrack &  fastTrack)
private

Definition at line 145 of file GFlashHadronShowerModel.cc.

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

Referenced by DoIt().

145  {
146 
147  std::vector<GflashHit>& gflashHitList = theProfile->getGflashHitList();
148 
149  theGflashStep->SetTrack(const_cast<G4Track*>(fastTrack.GetPrimaryTrack()));
150  theGflashStep->GetPostStepPoint()->SetProcessDefinedStep(const_cast<G4VProcess*>
151  (fastTrack.GetPrimaryTrack()->GetStep()->GetPostStepPoint()->GetProcessDefinedStep()));
152  theGflashNavigator->SetWorldVolume(G4TransportationManager::GetTransportationManager()->GetNavigatorForTracking()->GetWorldVolume());
153 
154  std::vector<GflashHit>::const_iterator spotIter = gflashHitList.begin();
155  std::vector<GflashHit>::const_iterator spotIterEnd = gflashHitList.end();
156 
157  for( ; spotIter != spotIterEnd; spotIter++){
158 
159  theGflashNavigator->LocateGlobalPointAndUpdateTouchableHandle(spotIter->getPosition(),
160  G4ThreeVector(0,0,0),
161  theGflashTouchableHandle, false);
162  updateGflashStep(spotIter->getPosition(),spotIter->getTime());
163 
164  // if there is a watcher defined in a job and the flag is turned on
165  if(theWatcherOn) {
166  theGflashStep->SetTotalEnergyDeposit(spotIter->getEnergy());
167  SteppingAction* userSteppingAction =
168  (SteppingAction*) G4EventManager::GetEventManager()->GetUserSteppingAction();
169  userSteppingAction->m_g4StepSignal(theGflashStep);
170  }
171 
172  G4VPhysicalVolume* aCurrentVolume = theGflashStep->GetPreStepPoint()->GetPhysicalVolume();
173  if( aCurrentVolume == 0 ) { continue; }
174 
175  G4LogicalVolume* lv = aCurrentVolume->GetLogicalVolume();
176  if(lv->GetRegion() != theRegion) { continue; }
177 
178  theGflashStep->GetPreStepPoint()->SetSensitiveDetector(aCurrentVolume->GetLogicalVolume()->GetSensitiveDetector());
179  G4VSensitiveDetector* aSensitive = theGflashStep->GetPreStepPoint()->GetSensitiveDetector();
180 
181  if( aSensitive == 0 ) continue;
182 
183  G4String nameCalor = aCurrentVolume->GetName();
184  nameCalor.assign(nameCalor,0,2);
185  G4double samplingWeight = 1.0;
186  if(nameCalor == "HB" ) {
187  samplingWeight = Gflash::scaleSensitiveHB;
188  }
189  else if(nameCalor=="HE" || nameCalor == "HT") {
190  samplingWeight = Gflash::scaleSensitiveHE;
191  }
192  theGflashStep->SetTotalEnergyDeposit(spotIter->getEnergy()*samplingWeight);
193 
194  aSensitive->Hit(theGflashStep);
195 
196  }
197 }
SimActivityRegistry::G4StepSignal m_g4StepSignal
G4TouchableHandle theGflashTouchableHandle
GflashHadronShowerProfile * theProfile
void updateGflashStep(const G4ThreeVector &position, G4double time)
const double scaleSensitiveHE
std::vector< GflashHit > & getGflashHitList()
const double scaleSensitiveHB
G4bool GFlashHadronShowerModel::ModelTrigger ( const G4FastTrack &  fastTrack)

Definition at line 78 of file GFlashHadronShowerModel.cc.

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

79 {
80  // ModelTrigger returns false for Gflash Hadronic Shower Model if it is not
81  // tested from the corresponding wrapper process, GflashHadronWrapperProcess.
82  // Temporarily the track status is set to fPostponeToNextEvent at the wrapper
83  // process before ModelTrigger is really tested for the second time through
84  // PostStepGPIL of enviaged parameterization processes. The better
85  // implmentation may be using via G4VUserTrackInformation of each track, which
86  // requires to modify a geant source code of stepping (G4SteppingManager2)
87 
88  // mininum energy cutoff to parameterize
89  if (fastTrack.GetPrimaryTrack()->GetKineticEnergy() < GeV*Gflash::energyCutOff )
90  { return false; }
91 
92  // This will be changed accordingly when the way
93  // dealing with CaloRegion changes later.
94  G4TouchableHistory* touch =
95  (G4TouchableHistory*)(fastTrack.GetPrimaryTrack()->GetTouchable());
96  G4VPhysicalVolume* pCurrentVolume = touch->GetVolume();
97  if( pCurrentVolume == 0) { return false; }
98 
99  G4LogicalVolume* lv = pCurrentVolume->GetLogicalVolume();
100  if(lv->GetRegion() != theRegion) { return false; }
101 
102  // check whether this is called from the normal GPIL or the wrapper process GPIL
103  //if(fastTrack.GetPrimaryTrack()->GetTrackStatus() == fPostponeToNextEvent ) {
104 
105  // Shower pameterization start at the first inelastic interaction point
106  //if(isFirstInelasticInteraction(fastTrack) && excludeDetectorRegion(fastTrack))
107  if(excludeDetectorRegion(fastTrack))
108  { return false; }
109 
110  return true;
111 }
const double GeV
Definition: MathUtil.h:16
G4bool excludeDetectorRegion(const G4FastTrack &fastTrack)
const double energyCutOff
void GFlashHadronShowerModel::updateGflashStep ( const G4ThreeVector &  position,
G4double  time 
)
private

Definition at line 199 of file GFlashHadronShowerModel.cc.

References theGflashStep, and theGflashTouchableHandle.

Referenced by makeHits().

200 {
201  theGflashStep->GetPostStepPoint()->SetGlobalTime(timeGlobal);
202  theGflashStep->GetPreStepPoint()->SetPosition(spotPosition);
203  theGflashStep->GetPostStepPoint()->SetPosition(spotPosition);
204  theGflashStep->GetPreStepPoint()->SetTouchableHandle(theGflashTouchableHandle);
205 }
G4TouchableHandle theGflashTouchableHandle

Member Data Documentation

GflashAntiProtonShowerProfile* GFlashHadronShowerModel::theAntiProtonProfile
private
G4Navigator* GFlashHadronShowerModel::theGflashNavigator
private

Definition at line 58 of file GFlashHadronShowerModel.h.

Referenced by GFlashHadronShowerModel(), and makeHits().

G4Step* GFlashHadronShowerModel::theGflashStep
private
G4TouchableHandle GFlashHadronShowerModel::theGflashTouchableHandle
private

Definition at line 59 of file GFlashHadronShowerModel.h.

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

GflashKaonMinusShowerProfile* GFlashHadronShowerModel::theKaonMinusProfile
private
GflashKaonPlusShowerProfile* GFlashHadronShowerModel::theKaonPlusProfile
private
edm::ParameterSet GFlashHadronShowerModel::theParSet
private

Definition at line 47 of file GFlashHadronShowerModel.h.

Referenced by excludeDetectorRegion().

GflashPiKShowerProfile* GFlashHadronShowerModel::thePiKProfile
private
GflashHadronShowerProfile* GFlashHadronShowerModel::theProfile
private

Definition at line 48 of file GFlashHadronShowerModel.h.

Referenced by DoIt(), GFlashHadronShowerModel(), and makeHits().

GflashProtonShowerProfile* GFlashHadronShowerModel::theProtonProfile
private
const G4Region* GFlashHadronShowerModel::theRegion
private

Definition at line 55 of file GFlashHadronShowerModel.h.

Referenced by GFlashHadronShowerModel(), makeHits(), and ModelTrigger().

G4bool GFlashHadronShowerModel::theWatcherOn
private

Definition at line 46 of file GFlashHadronShowerModel.h.

Referenced by GFlashHadronShowerModel(), and makeHits().