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 
4 
14 
15 #include "G4FastSimulationManager.hh"
16 #include "G4TransportationManager.hh"
17 #include "G4TouchableHandle.hh"
18 #include "G4VSensitiveDetector.hh"
19 #include "G4VPhysicalVolume.hh"
20 
21 #include "G4PionMinus.hh"
22 #include "G4PionPlus.hh"
23 #include "G4KaonMinus.hh"
24 #include "G4KaonPlus.hh"
25 #include "G4Proton.hh"
26 #include "G4AntiProton.hh"
27 #include "G4VProcess.hh"
28 
29 #include <vector>
30 
31 GflashHadronShowerModel::GflashHadronShowerModel(G4String modelName, G4Region* envelope, edm::ParameterSet parSet)
32  : G4VFastSimulationModel(modelName, envelope), theParSet(parSet)
33 {
41 
42  theGflashStep = new G4Step();
43  theGflashTouchableHandle = new G4TouchableHistory();
44  theGflashNavigator = new G4Navigator();
45 
46 }
47 
49 {
50  if(theProfile) delete theProfile;
51  if(theGflashStep) delete theGflashStep;
52 }
53 
54 G4bool GflashHadronShowerModel::IsApplicable(const G4ParticleDefinition& particleType)
55 {
56  return
57  &particleType == G4PionMinus::PionMinusDefinition() ||
58  &particleType == G4PionPlus::PionPlusDefinition() ||
59  &particleType == G4KaonMinus::KaonMinusDefinition() ||
60  &particleType == G4KaonPlus::KaonPlusDefinition() ||
61  &particleType == G4AntiProton::AntiProtonDefinition() ||
62  &particleType == G4Proton::ProtonDefinition() ;
63 }
64 
65 G4bool GflashHadronShowerModel::ModelTrigger(const G4FastTrack& fastTrack)
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) return trigger;
79 
80  // check whether this is called from the normal GPIL or the wrapper process GPIL
81  if(fastTrack.GetPrimaryTrack()->GetTrackStatus() == fPostponeToNextEvent ) {
82 
83  // Shower pameterization start at the first inelastic interaction point
84  G4bool isInelastic = isFirstInelasticInteraction(fastTrack);
85 
86  // Other conditions
87  if(isInelastic) {
88  trigger = (!excludeDetectorRegion(fastTrack));
89  }
90  }
91 
92  return trigger;
93 
94 }
95 
96 void GflashHadronShowerModel::DoIt(const G4FastTrack& fastTrack, G4FastStep& fastStep)
97 {
98  // kill the particle
99  fastStep.KillPrimaryTrack();
100  fastStep.ProposePrimaryTrackPathLength(0.0);
101  fastStep.ProposeTotalEnergyDeposited(fastTrack.GetPrimaryTrack()->GetKineticEnergy());
102 
103  // parameterize energy depostion by the particle type
104  G4ParticleDefinition* particleType = fastTrack.GetPrimaryTrack()->GetDefinition();
105 
107  if(particleType == G4KaonMinus::KaonMinusDefinition()) theProfile = theKaonMinusProfile;
108  else if(particleType == G4KaonPlus::KaonPlusDefinition()) theProfile = theKaonPlusProfile;
109  else if(particleType == G4AntiProton::AntiProtonDefinition()) theProfile = theAntiProtonProfile;
110  else if(particleType == G4Proton::ProtonDefinition()) theProfile = theProtonProfile;
111 
112  //input variables for GflashHadronShowerProfile
113  G4double energy = fastTrack.GetPrimaryTrack()->GetKineticEnergy()/GeV;
114  G4double globalTime = fastTrack.GetPrimaryTrack()->GetStep()->GetPostStepPoint()->GetGlobalTime();
115  G4double charge = fastTrack.GetPrimaryTrack()->GetStep()->GetPreStepPoint()->GetCharge();
116  G4ThreeVector position = fastTrack.GetPrimaryTrack()->GetPosition()/cm;
117  G4ThreeVector momentum = fastTrack.GetPrimaryTrack()->GetMomentum()/GeV;
118  G4int showerType = Gflash::findShowerType(position);
119 
120  theProfile->initialize(showerType,energy,globalTime,charge,position,momentum);
123 
124  // make hits
125  makeHits(fastTrack);
126 
127 }
128 
129 void GflashHadronShowerModel::makeHits(const G4FastTrack& fastTrack) {
130 
131  std::vector<GflashHit>& gflashHitList = theProfile->getGflashHitList();
132 
133  theGflashStep->SetTrack(const_cast<G4Track*>(fastTrack.GetPrimaryTrack()));
134  theGflashStep->GetPostStepPoint()->SetProcessDefinedStep(const_cast<G4VProcess*>
135  (fastTrack.GetPrimaryTrack()->GetStep()->GetPostStepPoint()->GetProcessDefinedStep()));
136  theGflashNavigator->SetWorldVolume(G4TransportationManager::GetTransportationManager()->GetNavigatorForTracking()->GetWorldVolume());
137 
138  std::vector<GflashHit>::const_iterator spotIter = gflashHitList.begin();
139  std::vector<GflashHit>::const_iterator spotIterEnd = gflashHitList.end();
140 
141  for( ; spotIter != spotIterEnd; spotIter++){
142 
143  theGflashNavigator->LocateGlobalPointAndUpdateTouchableHandle(spotIter->getPosition(),G4ThreeVector(0,0,0),
144  theGflashTouchableHandle, false);
145  updateGflashStep(spotIter->getPosition(),spotIter->getTime());
146 
147  G4VPhysicalVolume* aCurrentVolume = theGflashStep->GetPreStepPoint()->GetPhysicalVolume();
148  if( aCurrentVolume == 0 ) continue;
149 
150  G4LogicalVolume* lv = aCurrentVolume->GetLogicalVolume();
151  if(lv->GetRegion()->GetName() != "CaloRegion") continue;
152 
153  theGflashStep->GetPreStepPoint()->SetSensitiveDetector(aCurrentVolume->GetLogicalVolume()->GetSensitiveDetector());
154  G4VSensitiveDetector* aSensitive = theGflashStep->GetPreStepPoint()->GetSensitiveDetector();
155 
156  if( aSensitive == 0 ) continue;
157 
158  G4String nameCalor = aCurrentVolume->GetName();
159  nameCalor.assign(nameCalor,0,2);
160  G4double samplingWeight = 1.0;
161  if(nameCalor == "HB" ) {
162  samplingWeight = Gflash::scaleSensitiveHB;
163  }
164  else if(nameCalor=="HE" || nameCalor == "HT") {
165  samplingWeight = Gflash::scaleSensitiveHE;
166  }
167  theGflashStep->SetTotalEnergyDeposit(spotIter->getEnergy()*samplingWeight);
168 
169  aSensitive->Hit(theGflashStep);
170 
171  }
172 }
173 
174 void GflashHadronShowerModel::updateGflashStep(G4ThreeVector spotPosition, G4double timeGlobal)
175 {
176  theGflashStep->GetPostStepPoint()->SetGlobalTime(timeGlobal);
177  theGflashStep->GetPreStepPoint()->SetPosition(spotPosition);
178  theGflashStep->GetPostStepPoint()->SetPosition(spotPosition);
179  theGflashStep->GetPreStepPoint()->SetTouchableHandle(theGflashTouchableHandle);
180 }
181 
182 G4bool GflashHadronShowerModel::isFirstInelasticInteraction(const G4FastTrack& fastTrack)
183 {
184  G4bool isFirst=false;
185 
186  G4StepPoint* preStep = fastTrack.GetPrimaryTrack()->GetStep()->GetPreStepPoint();
187  G4StepPoint* postStep = fastTrack.GetPrimaryTrack()->GetStep()->GetPostStepPoint();
188 
189  G4String procName = postStep->GetProcessDefinedStep()->GetProcessName();
190  G4ParticleDefinition* particleType = fastTrack.GetPrimaryTrack()->GetDefinition();
191 
192  //@@@ this part is still temporary and the cut for the variable ratio should 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 
201  //skip to the second interaction if the first inelastic is a quasi-elastic 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 processes.
208  //@@@may require an additional condition only for hadron interaction with the process name,
209  //but it will not change the result anyway
210 
211  for (unsigned int isec = 0 ; isec < fSecondaryVector->size() ; isec++) {
212  G4Track* fSecondaryTrack = (*fSecondaryVector)[isec];
213  G4double secondaryEnergy = fSecondaryTrack->GetKineticEnergy();
214 
215  if(secondaryEnergy > leadingEnergy ) {
216  leadingEnergy = secondaryEnergy;
217  }
218  }
219 
220  if((preStep->GetTotalEnergy()!=0) &&
221  (leadingEnergy/preStep->GetTotalEnergy() < Gflash::QuasiElasticLike)) isFirst = true;
222 
223  //Fill debugging histograms and check information on secondaries -
224  //remove after final implimentation
225 
226  if(theHisto->getStoreFlag()) {
227  theHisto->preStepPosition->Fill(preStep->GetPosition().getRho()/cm);
228  theHisto->postStepPosition->Fill(postStep->GetPosition().getRho()/cm);
229  theHisto->deltaStep->Fill((postStep->GetPosition() - preStep->GetPosition()).getRho()/cm);
230  theHisto->kineticEnergy->Fill(fastTrack.GetPrimaryTrack()->GetKineticEnergy()/GeV);
231  theHisto->energyLoss->Fill(fabs(fastTrack.GetPrimaryTrack()->GetStep()->GetDeltaEnergy()/GeV));
232  theHisto->energyRatio->Fill(leadingEnergy/preStep->GetTotalEnergy());
233  }
234 
235  }
236  return isFirst;
237 }
238 
239 G4bool GflashHadronShowerModel::excludeDetectorRegion(const G4FastTrack& fastTrack)
240 {
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::LogInfo("SimGeneralGFlash") << "GflashHadronShowerModel: excluding region of eta = " << eta;
250  }
251  return true;
252  }
253  else {
254  G4StepPoint* postStep = fastTrack.GetPrimaryTrack()->GetStep()->GetPostStepPoint();
255 
256  Gflash::CalorimeterNumber kCalor = Gflash::getCalorimeterNumber(postStep->GetPosition()/cm);
257  G4double distOut = 9999.0;
258 
259  //exclude the region where the shower starting point is inside the preshower
260  if( std::fabs(eta) > Gflash::EtaMin[Gflash::kENCA] &&
261  std::fabs((postStep->GetPosition()).getZ()/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 end of
267  //the hadronic envelopes (may need to be optimized further!)
268  //@@@if we extend parameterization including Magnet/HO, we need to change this strategy
269  if(kCalor == Gflash::kHB) {
270  distOut = Gflash::Rmax[Gflash::kHB] - postStep->GetPosition().getRho()/cm;
271  if (distOut < Gflash::MinDistanceToOut ) isExcluded = true;
272  }
273  else if(kCalor == Gflash::kHE) {
274  distOut = Gflash::Zmax[Gflash::kHE] - std::fabs(postStep->GetPosition().getZ()/cm);
275  if (distOut < Gflash::MinDistanceToOut ) isExcluded = true;
276  }
277 
278  //@@@remove this print statement later
279  if(isExcluded && verbosity > 0) {
280  std::cout << "GflashHadronShowerModel: skipping kCalor = " << kCalor <<
281  " DistanceToOut " << distOut << " from (" << (postStep->GetPosition()).getRho()/cm <<
282  ":" << (postStep->GetPosition()).getZ()/cm << ") of KE = " << fastTrack.GetPrimaryTrack()->GetKineticEnergy()/GeV << std::endl;
283  }
284  }
285 
286  return isExcluded;
287 }
T getUntrackedParameter(std::string const &, T const &) const
const double MinDistanceToOut
static GflashHistogram * instance()
void updateGflashStep(G4ThreeVector position, G4double time)
G4bool ModelTrigger(const G4FastTrack &)
const double Zmax[kNumberCalorimeter]
GflashKaonMinusShowerProfile * theKaonMinusProfile
GflashAntiProtonShowerProfile * theAntiProtonProfile
GflashHadronShowerModel(G4String modelName, G4Region *envelope, edm::ParameterSet parSet)
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]
const int verbosity
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:41
const double scaleSensitiveHB
G4bool IsApplicable(const G4ParticleDefinition &)
void initialize(int showerType, double energy, double globalTime, double charge, Gflash3Vector &position, Gflash3Vector &momentum)