CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
GflashHadronShowerModel.cc
Go to the documentation of this file.
2 
5 
15 
16 #include "G4FastSimulationManager.hh"
17 #include "G4TransportationManager.hh"
18 #include "G4TouchableHandle.hh"
19 #include "G4VSensitiveDetector.hh"
20 #include "G4VPhysicalVolume.hh"
21 #include "G4EventManager.hh"
22 
23 #include "G4PionMinus.hh"
24 #include "G4PionPlus.hh"
25 #include "G4KaonMinus.hh"
26 #include "G4KaonPlus.hh"
27 #include "G4Proton.hh"
28 #include "G4AntiProton.hh"
29 #include "G4VProcess.hh"
30 
31 #include <vector>
32 
33 GflashHadronShowerModel::GflashHadronShowerModel(G4String modelName, G4Region* envelope, const edm::ParameterSet& parSet)
34  : G4VFastSimulationModel(modelName, envelope), theParSet(parSet)
35 {
36  theWatcherOn = parSet.getParameter<bool>("watcherOn");
37 
45 
46  theGflashStep = new G4Step();
47  theGflashTouchableHandle = new G4TouchableHistory();
48  theGflashNavigator = new G4Navigator();
49 
50 }
51 
53 {
54  if(theProfile) delete theProfile;
55  if(theGflashStep) delete theGflashStep;
56 }
57 
58 G4bool GflashHadronShowerModel::IsApplicable(const G4ParticleDefinition& particleType)
59 {
60  return
61  &particleType == G4PionMinus::PionMinusDefinition() ||
62  &particleType == G4PionPlus::PionPlusDefinition() ||
63  &particleType == G4KaonMinus::KaonMinusDefinition() ||
64  &particleType == G4KaonPlus::KaonPlusDefinition() ||
65  &particleType == G4AntiProton::AntiProtonDefinition() ||
66  &particleType == G4Proton::ProtonDefinition() ;
67 }
68 
69 G4bool GflashHadronShowerModel::ModelTrigger(const G4FastTrack& fastTrack)
70 {
71  // ModelTrigger returns false for Gflash Hadronic Shower Model if it is not
72  // tested from the corresponding wrapper process, GflashHadronWrapperProcess.
73  // Temporarily the track status is set to fPostponeToNextEvent at the wrapper
74  // process before ModelTrigger is really tested for the second time through
75  // PostStepGPIL of enviaged parameterization processes. The better
76  // implmentation may be using via G4VUserTrackInformation of each track, which
77  // requires to modify a geant source code of stepping (G4SteppingManager2)
78 
79  G4bool trigger = false;
80 
81  // mininum energy cutoff to parameterize
82  if (fastTrack.GetPrimaryTrack()->GetKineticEnergy()/GeV < Gflash::energyCutOff) return trigger;
83 
84  // check whether this is called from the normal GPIL or the wrapper process GPIL
85  if(fastTrack.GetPrimaryTrack()->GetTrackStatus() == fPostponeToNextEvent ) {
86 
87  // Shower pameterization start at the first inelastic interaction point
88  G4bool isInelastic = isFirstInelasticInteraction(fastTrack);
89 
90  // Other conditions
91  if(isInelastic) {
92  trigger = (!excludeDetectorRegion(fastTrack));
93  }
94  }
95 
96  return trigger;
97 
98 }
99 
100 void GflashHadronShowerModel::DoIt(const G4FastTrack& fastTrack, G4FastStep& fastStep)
101 {
102  // kill the particle
103  fastStep.KillPrimaryTrack();
104  fastStep.ProposePrimaryTrackPathLength(0.0);
105 
106  // parameterize energy depostion by the particle type
107  G4ParticleDefinition* particleType = fastTrack.GetPrimaryTrack()->GetDefinition();
108 
110  if(particleType == G4KaonMinus::KaonMinusDefinition()) theProfile = theKaonMinusProfile;
111  else if(particleType == G4KaonPlus::KaonPlusDefinition()) theProfile = theKaonPlusProfile;
112  else if(particleType == G4AntiProton::AntiProtonDefinition()) theProfile = theAntiProtonProfile;
113  else if(particleType == G4Proton::ProtonDefinition()) theProfile = theProtonProfile;
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 
132 void GflashHadronShowerModel::makeHits(const G4FastTrack& fastTrack) {
133 
134  std::vector<GflashHit>& gflashHitList = theProfile->getGflashHitList();
135 
136  theGflashStep->SetTrack(const_cast<G4Track*>(fastTrack.GetPrimaryTrack()));
137  theGflashStep->GetPostStepPoint()->SetProcessDefinedStep(const_cast<G4VProcess*>
138  (fastTrack.GetPrimaryTrack()->GetStep()->GetPostStepPoint()->GetProcessDefinedStep()));
139  theGflashNavigator->SetWorldVolume(G4TransportationManager::GetTransportationManager()->GetNavigatorForTracking()->GetWorldVolume());
140 
141  std::vector<GflashHit>::const_iterator spotIter = gflashHitList.begin();
142  std::vector<GflashHit>::const_iterator spotIterEnd = gflashHitList.end();
143 
144  for( ; spotIter != spotIterEnd; spotIter++){
145 
146  theGflashNavigator->LocateGlobalPointAndUpdateTouchableHandle(spotIter->getPosition(),G4ThreeVector(0,0,0),
147  theGflashTouchableHandle, false);
148  updateGflashStep(spotIter->getPosition(),spotIter->getTime());
149 
150  // if there is a watcher defined in a job and the flag is turned on
151  if(theWatcherOn) {
152  theGflashStep->SetTotalEnergyDeposit(spotIter->getEnergy());
153  SteppingAction* userSteppingAction = (SteppingAction*) G4EventManager::GetEventManager()->GetUserSteppingAction();
154  userSteppingAction->m_g4StepSignal(theGflashStep);
155  }
156 
157  G4VPhysicalVolume* aCurrentVolume = theGflashStep->GetPreStepPoint()->GetPhysicalVolume();
158  if( aCurrentVolume == 0 ) continue;
159 
160  G4LogicalVolume* lv = aCurrentVolume->GetLogicalVolume();
161  if(lv->GetRegion()->GetName() != "CaloRegion") continue;
162 
163  theGflashStep->GetPreStepPoint()->SetSensitiveDetector(aCurrentVolume->GetLogicalVolume()->GetSensitiveDetector());
164  G4VSensitiveDetector* aSensitive = theGflashStep->GetPreStepPoint()->GetSensitiveDetector();
165 
166  if( aSensitive == 0 ) continue;
167 
168  G4String nameCalor = aCurrentVolume->GetName();
169  nameCalor.assign(nameCalor,0,2);
170  G4double samplingWeight = 1.0;
171  if(nameCalor == "HB" ) {
172  samplingWeight = Gflash::scaleSensitiveHB;
173  }
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  }
182 }
183 
184 void GflashHadronShowerModel::updateGflashStep(const G4ThreeVector& spotPosition, G4double timeGlobal)
185 {
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 {
194  G4bool isFirst=false;
195 
196  G4StepPoint* preStep = fastTrack.GetPrimaryTrack()->GetStep()->GetPreStepPoint();
197  G4StepPoint* postStep = fastTrack.GetPrimaryTrack()->GetStep()->GetPostStepPoint();
198 
199  G4String procName = postStep->GetProcessDefinedStep()->GetProcessName();
200  G4ParticleDefinition* particleType = fastTrack.GetPrimaryTrack()->GetDefinition();
201 
202  //@@@ this part is still temporary and the cut for the variable ratio should be optimized later
203 
204  if((particleType == G4PionPlus::PionPlusDefinition() && procName == "WrappedPionPlusInelastic") ||
205  (particleType == G4PionMinus::PionMinusDefinition() && procName == "WrappedPionMinusInelastic") ||
206  (particleType == G4KaonPlus::KaonPlusDefinition() && procName == "WrappedKaonPlusInelastic") ||
207  (particleType == G4KaonMinus::KaonMinusDefinition() && procName == "WrappedKaonMinusInelastic") ||
208  (particleType == G4AntiProton::AntiProtonDefinition() && procName == "WrappedAntiProtonInelastic") ||
209  (particleType == G4Proton::ProtonDefinition() && procName == "WrappedProtonInelastic")) {
210 
211  //skip to the second interaction if the first inelastic is a quasi-elastic like interaction
212  //@@@ the cut may be optimized later
213 
214  const G4TrackVector* fSecondaryVector = fastTrack.GetPrimaryTrack()->GetStep()->GetSecondary();
215  G4double leadingEnergy = 0.0;
216 
217  //loop over 'all' secondaries including those produced by continuous processes.
218  //@@@may require an additional condition only for hadron interaction with the process name,
219  //but it will not change the result anyway
220 
221  for (unsigned int isec = 0 ; isec < fSecondaryVector->size() ; isec++) {
222  G4Track* fSecondaryTrack = (*fSecondaryVector)[isec];
223  G4double secondaryEnergy = fSecondaryTrack->GetKineticEnergy();
224 
225  if(secondaryEnergy > leadingEnergy ) {
226  leadingEnergy = secondaryEnergy;
227  }
228  }
229 
230  if((preStep->GetTotalEnergy()!=0) &&
231  (leadingEnergy/preStep->GetTotalEnergy() < Gflash::QuasiElasticLike)) isFirst = true;
232 
233  //Fill debugging histograms and check information on secondaries -
234  //remove after final implimentation
235 
236  if(theHisto->getStoreFlag()) {
237  theHisto->preStepPosition->Fill(preStep->GetPosition().getRho()/cm);
238  theHisto->postStepPosition->Fill(postStep->GetPosition().getRho()/cm);
239  theHisto->deltaStep->Fill((postStep->GetPosition() - preStep->GetPosition()).getRho()/cm);
240  theHisto->kineticEnergy->Fill(fastTrack.GetPrimaryTrack()->GetKineticEnergy()/GeV);
241  theHisto->energyLoss->Fill(fabs(fastTrack.GetPrimaryTrack()->GetStep()->GetDeltaEnergy()/GeV));
242  theHisto->energyRatio->Fill(leadingEnergy/preStep->GetTotalEnergy());
243  }
244 
245  }
246  return isFirst;
247 }
248 
249 G4bool GflashHadronShowerModel::excludeDetectorRegion(const G4FastTrack& fastTrack)
250 {
251  G4bool isExcluded=false;
252  int verbosity = theParSet.getUntrackedParameter<int>("Verbosity");
253 
254  //exclude regions where geometry are complicated
255  //+- one supermodule around the EB/EE boundary: 1.479 +- 0.0174*5
256  G4double eta = fastTrack.GetPrimaryTrack()->GetPosition().pseudoRapidity() ;
257  if(std::fabs(eta) > 1.392 && std::fabs(eta) < 1.566) {
258  if(verbosity>0) {
259  edm::LogInfo("SimGeneralGFlash") << "GflashHadronShowerModel: excluding region of eta = " << eta;
260  }
261  return true;
262  }
263  else {
264  G4StepPoint* postStep = fastTrack.GetPrimaryTrack()->GetStep()->GetPostStepPoint();
265 
266  Gflash::CalorimeterNumber kCalor = Gflash::getCalorimeterNumber(postStep->GetPosition()/cm);
267  G4double distOut = 9999.0;
268 
269  //exclude the region where the shower starting point is inside the preshower
270  if( std::fabs(eta) > Gflash::EtaMin[Gflash::kENCA] &&
271  std::fabs((postStep->GetPosition()).getZ()/cm) < Gflash::Zmin[Gflash::kENCA]) {
272  return true;
273  }
274 
275  //<---the shower starting point is always inside envelopes
276  //@@@exclude the region where the shower starting point is too close to the end of
277  //the hadronic envelopes (may need to be optimized further!)
278  //@@@if we extend parameterization including Magnet/HO, we need to change this strategy
279  if(kCalor == Gflash::kHB) {
280  distOut = Gflash::Rmax[Gflash::kHB] - postStep->GetPosition().getRho()/cm;
281  if (distOut < Gflash::MinDistanceToOut ) isExcluded = true;
282  }
283  else if(kCalor == Gflash::kHE) {
284  distOut = Gflash::Zmax[Gflash::kHE] - std::fabs(postStep->GetPosition().getZ()/cm);
285  if (distOut < Gflash::MinDistanceToOut ) isExcluded = true;
286  }
287 
288  //@@@remove this print statement later
289  if(isExcluded && verbosity > 0) {
290  std::cout << "GflashHadronShowerModel: skipping kCalor = " << kCalor <<
291  " DistanceToOut " << distOut << " from (" << (postStep->GetPosition()).getRho()/cm <<
292  ":" << (postStep->GetPosition()).getZ()/cm << ") of KE = " << fastTrack.GetPrimaryTrack()->GetKineticEnergy()/GeV << std::endl;
293  }
294  }
295 
296  return isExcluded;
297 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
const double MinDistanceToOut
static GflashHistogram * instance()
G4bool ModelTrigger(const G4FastTrack &)
const double Zmax[kNumberCalorimeter]
SimActivityRegistry::G4StepSignal m_g4StepSignal
GflashKaonMinusShowerProfile * theKaonMinusProfile
GflashAntiProtonShowerProfile * theAntiProtonProfile
isFirst
Definition: cuy.py:417
void updateGflashStep(const G4ThreeVector &position, G4double time)
void makeHits(const G4FastTrack &fastTrack)
T eta() const
double charge(const std::vector< uint8_t > &Ampls)
static int position[TOTALCHAMBERS][3]
Definition: ReadPGInfo.cc:509
void DoIt(const G4FastTrack &, G4FastStep &)
GflashHadronShowerProfile * theProfile
int findShowerType(const Gflash3Vector position)
const double energyCutOff
const double EtaMin[kNumberCalorimeter]
G4TouchableHandle theGflashTouchableHandle
G4bool isFirstInelasticInteraction(const G4FastTrack &fastTrack)
const double Zmin[kNumberCalorimeter]
const double QuasiElasticLike
const double scaleSensitiveHE
G4bool excludeDetectorRegion(const G4FastTrack &fastTrack)
CalorimeterNumber getCalorimeterNumber(const Gflash3Vector position)
GflashKaonPlusShowerProfile * theKaonPlusProfile
GflashPiKShowerProfile * thePiKProfile
std::vector< GflashHit > & getGflashHitList()
GflashProtonShowerProfile * theProtonProfile
const double Rmax[kNumberCalorimeter]
tuple cout
Definition: gather_cfg.py:121
const double scaleSensitiveHB
G4bool IsApplicable(const G4ParticleDefinition &)
void initialize(int showerType, double energy, double globalTime, double charge, Gflash3Vector &position, Gflash3Vector &momentum)
GflashHadronShowerModel(G4String modelName, G4Region *envelope, const edm::ParameterSet &parSet)