CMS 3D CMS Logo

GflashEMShowerModel.cc
Go to the documentation of this file.
1 //
2 // initial setup : E.Barberio & Joanna Weng
3 // big changes : Soon Jun & Dongwook Jang
4 //
7 
10 
11 #include "G4Electron.hh"
12 #include "G4Positron.hh"
13 #include "G4VProcess.hh"
14 #include "G4VPhysicalVolume.hh"
15 #include "G4LogicalVolume.hh"
16 #include "G4TransportationManager.hh"
17 #include "G4EventManager.hh"
18 #include "G4FastSimulationManager.hh"
19 #include "G4TouchableHandle.hh"
20 #include "G4VSensitiveDetector.hh"
21 
22 using namespace CLHEP;
23 
25  G4Envelope* envelope,
26  const edm::ParameterSet& parSet)
27  : G4VFastSimulationModel(modelName, envelope), theParSet(parSet) {
28 
29  theProfile = new GflashEMShowerProfile(parSet);
30  theRegion = const_cast<const G4Region*>(envelope);
31 
32  theGflashStep = new G4Step();
33  theGflashTouchableHandle = new G4TouchableHistory();
34  theGflashNavigator = new G4Navigator();
35 
36 }
37 
38 // -----------------------------------------------------------------------------------
39 
41 
42  delete theProfile;
43  delete theGflashStep;
44 }
45 
46 G4bool GflashEMShowerModel::IsApplicable(const G4ParticleDefinition& particleType) {
47 
48  return ( &particleType == G4Electron::ElectronDefinition() ||
49  &particleType == G4Positron::PositronDefinition() );
50 
51 }
52 
53 // -----------------------------------------------------------------------------------
54 G4bool GflashEMShowerModel::ModelTrigger(const G4FastTrack & fastTrack ) {
55 
56  // Mininum energy cutoff to parameterize
57  if(fastTrack.GetPrimaryTrack()->GetKineticEnergy() < GeV) { return false; }
58  if(excludeDetectorRegion(fastTrack)) { return false; }
59 
60  // This will be changed accordingly when the way
61  // dealing with CaloRegion changes later.
62  G4VPhysicalVolume* pCurrentVolume = (fastTrack.GetPrimaryTrack()->GetTouchable())->GetVolume();
63  if( pCurrentVolume == nullptr) { return false; }
64 
65  G4LogicalVolume* lv = pCurrentVolume->GetLogicalVolume();
66  if(lv->GetRegion() != theRegion) { return false; }
67  return true;
68 }
69 
70 // -----------------------------------------------------------------------------------
71 void GflashEMShowerModel::DoIt(const G4FastTrack& fastTrack, G4FastStep& fastStep) {
72 
73  // Kill the parameterised particle:
74  fastStep.KillPrimaryTrack();
75  fastStep.ProposePrimaryTrackPathLength(0.0);
76 
77  //input variables for GflashEMShowerProfile with showerType = 1,5 (shower starts inside crystals)
78  G4double energy = fastTrack.GetPrimaryTrack()->GetKineticEnergy()/GeV;
79  G4double globalTime = fastTrack.GetPrimaryTrack()->GetStep()->GetPostStepPoint()->GetGlobalTime();
80  G4double charge = fastTrack.GetPrimaryTrack()->GetStep()->GetPreStepPoint()->GetCharge();
81  G4ThreeVector position = fastTrack.GetPrimaryTrack()->GetPosition() / cm;
82  G4ThreeVector momentum = fastTrack.GetPrimaryTrack()->GetMomentum()/GeV;
83  G4int showerType = Gflash::findShowerType(position);
84 
85  // Do actual parameterization. The result of parameterization is gflashHitList
86  theProfile->initialize(showerType,energy,globalTime,charge,position,momentum);
88 
89  //make hits
90  makeHits(fastTrack);
91 }
92 
93 void GflashEMShowerModel::makeHits(const G4FastTrack& fastTrack) {
94 
95  std::vector<GflashHit>& gflashHitList = theProfile->getGflashHitList();
96 
97  theGflashStep->SetTrack(const_cast<G4Track*>(fastTrack.GetPrimaryTrack()));
98 
99  theGflashStep->GetPostStepPoint()->SetProcessDefinedStep(const_cast<G4VProcess*>
100  (fastTrack.GetPrimaryTrack()->GetStep()->GetPostStepPoint()->GetProcessDefinedStep()));
101  theGflashNavigator->SetWorldVolume(G4TransportationManager::GetTransportationManager()->GetNavigatorForTracking()->GetWorldVolume());
102 
103  std::vector<GflashHit>::const_iterator spotIter = gflashHitList.begin();
104  std::vector<GflashHit>::const_iterator spotIterEnd = gflashHitList.end();
105 
106  for( ; spotIter != spotIterEnd; spotIter++){
107 
108  //put touchable for each hit so that touchable history keeps track of each step.
109  theGflashNavigator->LocateGlobalPointAndUpdateTouchableHandle(spotIter->getPosition(),G4ThreeVector(0,0,0),
110  theGflashTouchableHandle, false);
111  updateGflashStep(spotIter->getPosition(),spotIter->getTime());
112 
113  // Send G4Step information to Hit/Digi if the volume is sensitive
114  // Copied from G4SteppingManager.cc
115 
116  G4VPhysicalVolume* aCurrentVolume = theGflashStep->GetPreStepPoint()->GetPhysicalVolume();
117  if( aCurrentVolume == nullptr ) continue;
118 
119  G4LogicalVolume* lv = aCurrentVolume->GetLogicalVolume();
120  if(lv->GetRegion() != theRegion) continue;
121 
122  theGflashStep->GetPreStepPoint()->SetSensitiveDetector(aCurrentVolume->GetLogicalVolume()->GetSensitiveDetector());
123  G4VSensitiveDetector* aSensitive = theGflashStep->GetPreStepPoint()->GetSensitiveDetector();
124 
125  if( aSensitive == nullptr ) continue;
126 
127  theGflashStep->SetTotalEnergyDeposit(spotIter->getEnergy());
128  aSensitive->Hit(theGflashStep);
129  }
130 
131 }
132 
133 void GflashEMShowerModel::updateGflashStep(const G4ThreeVector& spotPosition,
134  G4double timeGlobal)
135 {
136  theGflashStep->GetPostStepPoint()->SetGlobalTime(timeGlobal);
137  theGflashStep->GetPreStepPoint()->SetPosition(spotPosition);
138  theGflashStep->GetPostStepPoint()->SetPosition(spotPosition);
139  theGflashStep->GetPreStepPoint()->SetTouchableHandle(theGflashTouchableHandle);
140 }
141 
142 // -----------------------------------------------------------------------------------
143 G4bool GflashEMShowerModel::excludeDetectorRegion(const G4FastTrack& fastTrack) {
144 
145  G4bool isExcluded=false;
146 
147  //exclude regions where geometry are complicated
148  //+- one supermodule around the EB/EE boundary: 1.479 +- 0.0174*5
149  G4double eta = fastTrack.GetPrimaryTrack()->GetPosition().pseudoRapidity() ;
150  if(std::fabs(eta) > 1.392 && std::fabs(eta) < 1.566) { return true; }
151 
152  return isExcluded;
153 }
154 
155 /*
156 G4int GflashEMShowerModel::findShowerType(const G4FastTrack& fastTrack)
157 {
158  // Initialization of longitudinal and lateral parameters for
159  // hadronic showers. Simulation of the intrinsic fluctuations
160 
161  // type of hadron showers subject to the shower starting point (ssp)
162  // showerType = -1 : default (invalid)
163  // showerType = 0 : ssp before EBRY (barrel crystal)
164  // showerType = 1 : ssp inside EBRY
165  // showerType = 2 : ssp after EBRY before HB
166  // showerType = 3 : ssp inside HB
167  // showerType = 4 : ssp before EFRY (endcap crystal)
168  // showerType = 5 : ssp inside EFRY
169  // showerType = 6 : ssp after EFRY before HE
170  // showerType = 7 : ssp inside HE
171 
172  G4TouchableHistory* touch = (G4TouchableHistory*)(fastTrack.GetPrimaryTrack()->GetTouchable());
173  G4LogicalVolume* lv = touch->GetVolume()->GetLogicalVolume();
174 
175  std::size_t pos1 = lv->GetName().find("EBRY");
176  std::size_t pos11 = lv->GetName().find("EWAL");
177  std::size_t pos12 = lv->GetName().find("EWRA");
178  std::size_t pos2 = lv->GetName().find("EFRY");
179 
180  G4ThreeVector position = fastTrack.GetPrimaryTrack()->GetPosition()/cm;
181  Gflash::CalorimeterNumber kCalor = Gflash::getCalorimeterNumber(position);
182 
183  G4int showerType = -1;
184 
185  //central
186  if (kCalor == Gflash::kESPM || kCalor == Gflash::kHB ) {
187 
188  G4double posRho = position.getRho();
189 
190  if(pos1 != std::string::npos || pos11 != std::string::npos || pos12 != std::string::npos ) {
191  showerType = 1;
192  }
193  else {
194  if(kCalor == Gflash::kESPM) {
195  showerType = 2;
196  if( posRho < Gflash::Rmin[Gflash::kESPM]+ Gflash::ROffCrystalEB ) showerType = 0;
197  }
198  else showerType = 3;
199  }
200 
201  }
202  //forward
203  else if (kCalor == Gflash::kENCA || kCalor == Gflash::kHE) {
204  if(pos2 != std::string::npos) {
205  showerType = 5;
206  }
207  else {
208  if(kCalor == Gflash::kENCA) {
209  showerType = 6;
210  if(fabs(position.getZ()) < Gflash::Zmin[Gflash::kENCA] + Gflash::ZOffCrystalEE) showerType = 4;
211  }
212  else showerType = 7;
213  }
214  //@@@need z-dependent correction on the mean energy reponse
215  }
216 
217  return showerType;
218 }
219 */
const double GeV
Definition: MathUtil.h:16
std::vector< GflashHit > & getGflashHitList()
G4bool excludeDetectorRegion(const G4FastTrack &fastTrack)
G4TouchableHandle theGflashTouchableHandle
const G4Region * theRegion
G4bool IsApplicable(const G4ParticleDefinition &) override
GflashEMShowerProfile * theProfile
GflashEMShowerModel(const G4String &name, G4Envelope *env, const edm::ParameterSet &parSet)
void initialize(int showerType, double energy, double globalTime, double charge, Gflash3Vector &position, Gflash3Vector &momentum)
int findShowerType(const Gflash3Vector &position)
void updateGflashStep(const G4ThreeVector &position, G4double time)
static const G4LogicalVolume * GetVolume(const std::string &name)
void DoIt(const G4FastTrack &, G4FastStep &) override
static int position[264][3]
Definition: ReadPGInfo.cc:509
G4Navigator * theGflashNavigator
void makeHits(const G4FastTrack &fastTrack)
G4bool ModelTrigger(const G4FastTrack &) override