CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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
 
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.

36  : G4VFastSimulationModel(modelName, envelope), theParSet(parSet)
37 {
38  theWatcherOn = parSet.getParameter<bool>("watcherOn");
39 
47 
48  theGflashStep = new G4Step();
49  theGflashTouchableHandle = new G4TouchableHistory();
50  theGflashNavigator = new G4Navigator();
51 
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 ( )

Definition at line 54 of file GflashHadronShowerModel.cc.

References theGflashStep, and theProfile.

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

Member Function Documentation

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

Definition at line 102 of file GflashHadronShowerModel.cc.

References RecoTauCleanerPlugins::charge, relval_parameters_module::energy, Gflash::findShowerType(), GeV, GflashHadronShowerProfile::hadronicParameterization(), GflashHadronShowerProfile::initialize(), GflashHadronShowerProfile::loadParameters(), makeHits(), HLT_FULL_cff::particleType, position, theAntiProtonProfile, theKaonMinusProfile, theKaonPlusProfile, thePiKProfile, theProfile, and theProtonProfile.

103 {
104  // kill the particle
105  fastStep.KillPrimaryTrack();
106  fastStep.ProposePrimaryTrackPathLength(0.0);
107 
108  // parameterize energy depostion by the particle type
109  G4ParticleDefinition* particleType = fastTrack.GetPrimaryTrack()->GetDefinition();
110 
112  if(particleType == G4KaonMinus::KaonMinusDefinition()) theProfile = theKaonMinusProfile;
113  else if(particleType == G4KaonPlus::KaonPlusDefinition()) theProfile = theKaonPlusProfile;
114  else if(particleType == G4AntiProton::AntiProtonDefinition()) theProfile = theAntiProtonProfile;
115  else if(particleType == G4Proton::ProtonDefinition()) theProfile = theProtonProfile;
116 
117  //input variables for GflashHadronShowerProfile
118  G4double energy = fastTrack.GetPrimaryTrack()->GetKineticEnergy()/GeV;
119  G4double globalTime = fastTrack.GetPrimaryTrack()->GetStep()->GetPostStepPoint()->GetGlobalTime();
120  G4double charge = fastTrack.GetPrimaryTrack()->GetStep()->GetPreStepPoint()->GetCharge();
121  G4ThreeVector position = fastTrack.GetPrimaryTrack()->GetPosition()/cm;
122  G4ThreeVector momentum = fastTrack.GetPrimaryTrack()->GetMomentum()/GeV;
123  G4int showerType = Gflash::findShowerType(position);
124 
125  theProfile->initialize(showerType,energy,globalTime,charge,position,momentum);
128 
129  // make hits
130  makeHits(fastTrack);
131 
132 }
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
tuple particleType
void initialize(int showerType, double energy, double globalTime, double charge, Gflash3Vector &position, Gflash3Vector &momentum)
G4bool GflashHadronShowerModel::excludeDetectorRegion ( const G4FastTrack &  fastTrack)
private

Definition at line 251 of file GflashHadronShowerModel.cc.

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

Referenced by ModelTrigger().

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

Definition at line 60 of file GflashHadronShowerModel.cc.

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

Definition at line 194 of file GflashHadronShowerModel.cc.

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

Referenced by ModelTrigger().

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

Definition at line 134 of file GflashHadronShowerModel.cc.

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

Referenced by DoIt().

134  {
135 
136  std::vector<GflashHit>& gflashHitList = theProfile->getGflashHitList();
137 
138  theGflashStep->SetTrack(const_cast<G4Track*>(fastTrack.GetPrimaryTrack()));
139  theGflashStep->GetPostStepPoint()->SetProcessDefinedStep(const_cast<G4VProcess*>
140  (fastTrack.GetPrimaryTrack()->GetStep()->GetPostStepPoint()->GetProcessDefinedStep()));
141  theGflashNavigator->SetWorldVolume(G4TransportationManager::GetTransportationManager()->GetNavigatorForTracking()->GetWorldVolume());
142 
143  std::vector<GflashHit>::const_iterator spotIter = gflashHitList.begin();
144  std::vector<GflashHit>::const_iterator spotIterEnd = gflashHitList.end();
145 
146  for( ; spotIter != spotIterEnd; spotIter++){
147 
148  theGflashNavigator->LocateGlobalPointAndUpdateTouchableHandle(spotIter->getPosition(),G4ThreeVector(0,0,0),
149  theGflashTouchableHandle, false);
150  updateGflashStep(spotIter->getPosition(),spotIter->getTime());
151 
152  // if there is a watcher defined in a job and the flag is turned on
153  if(theWatcherOn) {
154  theGflashStep->SetTotalEnergyDeposit(spotIter->getEnergy());
155  SteppingAction* userSteppingAction = (SteppingAction*) G4EventManager::GetEventManager()->GetUserSteppingAction();
156  userSteppingAction->m_g4StepSignal(theGflashStep);
157  }
158 
159  G4VPhysicalVolume* aCurrentVolume = theGflashStep->GetPreStepPoint()->GetPhysicalVolume();
160  if( aCurrentVolume == 0 ) continue;
161 
162  G4LogicalVolume* lv = aCurrentVolume->GetLogicalVolume();
163  if(lv->GetRegion()->GetName() != "CaloRegion") continue;
164 
165  theGflashStep->GetPreStepPoint()->SetSensitiveDetector(aCurrentVolume->GetLogicalVolume()->GetSensitiveDetector());
166  G4VSensitiveDetector* aSensitive = theGflashStep->GetPreStepPoint()->GetSensitiveDetector();
167 
168  if( aSensitive == 0 ) continue;
169 
170  G4String nameCalor = aCurrentVolume->GetName();
171  nameCalor.assign(nameCalor,0,2);
172  G4double samplingWeight = 1.0;
173  if(nameCalor == "HB" ) {
174  samplingWeight = Gflash::scaleSensitiveHB;
175  }
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 }
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)

Definition at line 71 of file GflashHadronShowerModel.cc.

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

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

References theGflashStep, and theGflashTouchableHandle.

Referenced by makeHits().

187 {
188  theGflashStep->GetPostStepPoint()->SetGlobalTime(timeGlobal);
189  theGflashStep->GetPreStepPoint()->SetPosition(spotPosition);
190  theGflashStep->GetPostStepPoint()->SetPosition(spotPosition);
191  theGflashStep->GetPreStepPoint()->SetTouchableHandle(theGflashTouchableHandle);
192 }
G4TouchableHandle theGflashTouchableHandle

Member Data Documentation

GflashAntiProtonShowerProfile* GflashHadronShowerModel::theAntiProtonProfile
private

Definition at line 51 of file GflashHadronShowerModel.h.

Referenced by DoIt(), and GflashHadronShowerModel().

G4Navigator* GflashHadronShowerModel::theGflashNavigator
private

Definition at line 54 of file GflashHadronShowerModel.h.

Referenced by GflashHadronShowerModel(), and makeHits().

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

Definition at line 55 of file GflashHadronShowerModel.h.

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

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

Definition at line 49 of file GflashHadronShowerModel.h.

Referenced by DoIt(), and GflashHadronShowerModel().

GflashKaonPlusShowerProfile* GflashHadronShowerModel::theKaonPlusProfile
private

Definition at line 48 of file GflashHadronShowerModel.h.

Referenced by DoIt(), and GflashHadronShowerModel().

edm::ParameterSet GflashHadronShowerModel::theParSet
private

Definition at line 45 of file GflashHadronShowerModel.h.

Referenced by excludeDetectorRegion().

GflashPiKShowerProfile* GflashHadronShowerModel::thePiKProfile
private

Definition at line 47 of file GflashHadronShowerModel.h.

Referenced by DoIt(), and GflashHadronShowerModel().

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

Definition at line 50 of file GflashHadronShowerModel.h.

Referenced by DoIt(), and GflashHadronShowerModel().

G4bool GflashHadronShowerModel::theWatcherOn
private

Definition at line 44 of file GflashHadronShowerModel.h.

Referenced by GflashHadronShowerModel(), and makeHits().