CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 GflashEMShowerModel::GflashEMShowerModel(const G4String& modelName,
23  G4Envelope* envelope,
24  const edm::ParameterSet& parSet)
25  : G4VFastSimulationModel(modelName, envelope), theParSet(parSet) {
26 
27  //theWatcherOn = parSet.getParameter<bool>("watcherOn");
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  G4TouchableHistory* touch =
63  (G4TouchableHistory*)(fastTrack.GetPrimaryTrack()->GetTouchable());
64  G4VPhysicalVolume* pCurrentVolume = touch->GetVolume();
65  if( pCurrentVolume == 0) { return false; }
66 
67  G4LogicalVolume* lv = pCurrentVolume->GetLogicalVolume();
68  if(lv->GetRegion() != theRegion) { return false; }
69  //std::cout << "GflashEMShowerModel::ModelTrigger: LV "
70  // << lv->GetRegion()->GetName() << std::endl;
71 
72  // The parameterization starts inside crystals
73  //std::size_t pos1 = lv->GetName().find("EBRY");
74  //std::size_t pos2 = lv->GetName().find("EFRY");
75 
76  //std::size_t pos3 = lv->GetName().find("HVQ");
77  //std::size_t pos4 = lv->GetName().find("HF");
78  //if(pos1 == std::string::npos && pos2 == std::string::npos &&
79  // pos3 == std::string::npos && pos4 == std::string::npos) return false;
80  //@@@for now, HF is not a part of Gflash Envelopes
81  //if(pos1 == std::string::npos && pos2 == std::string::npos ) return false;
82 
83  return true;
84 
85 }
86 
87 // -----------------------------------------------------------------------------------
88 void GflashEMShowerModel::DoIt(const G4FastTrack& fastTrack, G4FastStep& fastStep) {
89 
90  // Kill the parameterised particle:
91  fastStep.KillPrimaryTrack();
92  fastStep.ProposePrimaryTrackPathLength(0.0);
93 
94  //input variables for GflashEMShowerProfile with showerType = 1,5 (shower starts inside crystals)
95  G4double energy = fastTrack.GetPrimaryTrack()->GetKineticEnergy()/GeV;
96  G4double globalTime = fastTrack.GetPrimaryTrack()->GetStep()->GetPostStepPoint()->GetGlobalTime();
97  G4double charge = fastTrack.GetPrimaryTrack()->GetStep()->GetPreStepPoint()->GetCharge();
98  G4ThreeVector position = fastTrack.GetPrimaryTrack()->GetPosition() / cm;
99  G4ThreeVector momentum = fastTrack.GetPrimaryTrack()->GetMomentum()/GeV;
100  G4int showerType = Gflash::findShowerType(position);
101 
102  // Do actual parameterization. The result of parameterization is gflashHitList
103  theProfile->initialize(showerType,energy,globalTime,charge,position,momentum);
105 
106  //make hits
107  makeHits(fastTrack);
108 }
109 
110 void GflashEMShowerModel::makeHits(const G4FastTrack& fastTrack) {
111 
112  std::vector<GflashHit>& gflashHitList = theProfile->getGflashHitList();
113 
114  theGflashStep->SetTrack(const_cast<G4Track*>(fastTrack.GetPrimaryTrack()));
115 
116  theGflashStep->GetPostStepPoint()->SetProcessDefinedStep(const_cast<G4VProcess*>
117  (fastTrack.GetPrimaryTrack()->GetStep()->GetPostStepPoint()->GetProcessDefinedStep()));
118  theGflashNavigator->SetWorldVolume(G4TransportationManager::GetTransportationManager()->GetNavigatorForTracking()->GetWorldVolume());
119 
120  std::vector<GflashHit>::const_iterator spotIter = gflashHitList.begin();
121  std::vector<GflashHit>::const_iterator spotIterEnd = gflashHitList.end();
122 
123  for( ; spotIter != spotIterEnd; spotIter++){
124 
125  //put touchable for each hit so that touchable history keeps track of each step.
126  theGflashNavigator->LocateGlobalPointAndUpdateTouchableHandle(spotIter->getPosition(),G4ThreeVector(0,0,0),
127  theGflashTouchableHandle, false);
128  updateGflashStep(spotIter->getPosition(),spotIter->getTime());
129 
130  // if there is a watcher defined in a job and the flag is turned on
131  if(theWatcherOn) {
132  SteppingAction* userSteppingAction = (SteppingAction*) G4EventManager::GetEventManager()->GetUserSteppingAction();
133  userSteppingAction->m_g4StepSignal(theGflashStep);
134  }
135 
136  // Send G4Step information to Hit/Digi if the volume is sensitive
137  // Copied from G4SteppingManager.cc
138 
139  G4VPhysicalVolume* aCurrentVolume = theGflashStep->GetPreStepPoint()->GetPhysicalVolume();
140  if( aCurrentVolume == 0 ) continue;
141 
142  G4LogicalVolume* lv = aCurrentVolume->GetLogicalVolume();
143  if(lv->GetRegion() != theRegion) continue;
144 
145  theGflashStep->GetPreStepPoint()->SetSensitiveDetector(aCurrentVolume->GetLogicalVolume()->GetSensitiveDetector());
146  G4VSensitiveDetector* aSensitive = theGflashStep->GetPreStepPoint()->GetSensitiveDetector();
147 
148  if( aSensitive == 0 ) continue;
149 
150  theGflashStep->SetTotalEnergyDeposit(spotIter->getEnergy());
151  aSensitive->Hit(theGflashStep);
152  }
153 
154 }
155 
156 void GflashEMShowerModel::updateGflashStep(const G4ThreeVector& spotPosition,
157  G4double timeGlobal)
158 {
159  theGflashStep->GetPostStepPoint()->SetGlobalTime(timeGlobal);
160  theGflashStep->GetPreStepPoint()->SetPosition(spotPosition);
161  theGflashStep->GetPostStepPoint()->SetPosition(spotPosition);
162  theGflashStep->GetPreStepPoint()->SetTouchableHandle(theGflashTouchableHandle);
163 }
164 
165 // -----------------------------------------------------------------------------------
166 G4bool GflashEMShowerModel::excludeDetectorRegion(const G4FastTrack& fastTrack) {
167 
168  G4bool isExcluded=false;
169 
170  //exclude regions where geometry are complicated
171  //+- one supermodule around the EB/EE boundary: 1.479 +- 0.0174*5
172  G4double eta = fastTrack.GetPrimaryTrack()->GetPosition().pseudoRapidity() ;
173  if(std::fabs(eta) > 1.392 && std::fabs(eta) < 1.566) { return true; }
174 
175  return isExcluded;
176 }
177 
178 /*
179 G4int GflashEMShowerModel::findShowerType(const G4FastTrack& fastTrack)
180 {
181  // Initialization of longitudinal and lateral parameters for
182  // hadronic showers. Simulation of the intrinsic fluctuations
183 
184  // type of hadron showers subject to the shower starting point (ssp)
185  // showerType = -1 : default (invalid)
186  // showerType = 0 : ssp before EBRY (barrel crystal)
187  // showerType = 1 : ssp inside EBRY
188  // showerType = 2 : ssp after EBRY before HB
189  // showerType = 3 : ssp inside HB
190  // showerType = 4 : ssp before EFRY (endcap crystal)
191  // showerType = 5 : ssp inside EFRY
192  // showerType = 6 : ssp after EFRY before HE
193  // showerType = 7 : ssp inside HE
194 
195  G4TouchableHistory* touch = (G4TouchableHistory*)(fastTrack.GetPrimaryTrack()->GetTouchable());
196  G4LogicalVolume* lv = touch->GetVolume()->GetLogicalVolume();
197 
198  std::size_t pos1 = lv->GetName().find("EBRY");
199  std::size_t pos11 = lv->GetName().find("EWAL");
200  std::size_t pos12 = lv->GetName().find("EWRA");
201  std::size_t pos2 = lv->GetName().find("EFRY");
202 
203  G4ThreeVector position = fastTrack.GetPrimaryTrack()->GetPosition()/cm;
204  Gflash::CalorimeterNumber kCalor = Gflash::getCalorimeterNumber(position);
205 
206  G4int showerType = -1;
207 
208  //central
209  if (kCalor == Gflash::kESPM || kCalor == Gflash::kHB ) {
210 
211  G4double posRho = position.getRho();
212 
213  if(pos1 != std::string::npos || pos11 != std::string::npos || pos12 != std::string::npos ) {
214  showerType = 1;
215  }
216  else {
217  if(kCalor == Gflash::kESPM) {
218  showerType = 2;
219  if( posRho < Gflash::Rmin[Gflash::kESPM]+ Gflash::ROffCrystalEB ) showerType = 0;
220  }
221  else showerType = 3;
222  }
223 
224  }
225  //forward
226  else if (kCalor == Gflash::kENCA || kCalor == Gflash::kHE) {
227  if(pos2 != std::string::npos) {
228  showerType = 5;
229  }
230  else {
231  if(kCalor == Gflash::kENCA) {
232  showerType = 6;
233  if(fabs(position.getZ()) < Gflash::Zmin[Gflash::kENCA] + Gflash::ZOffCrystalEE) showerType = 4;
234  }
235  else showerType = 7;
236  }
237  //@@@need z-dependent correction on the mean energy reponse
238  }
239 
240  return showerType;
241 }
242 */
std::vector< GflashHit > & getGflashHitList()
SimActivityRegistry::G4StepSignal m_g4StepSignal
void DoIt(const G4FastTrack &, G4FastStep &)
G4bool excludeDetectorRegion(const G4FastTrack &fastTrack)
G4TouchableHandle theGflashTouchableHandle
T eta() const
double charge(const std::vector< uint8_t > &Ampls)
static int position[TOTALCHAMBERS][3]
Definition: ReadPGInfo.cc:509
const G4Region * theRegion
int findShowerType(const Gflash3Vector position)
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)
G4bool IsApplicable(const G4ParticleDefinition &)
void updateGflashStep(const G4ThreeVector &position, G4double time)
G4Navigator * theGflashNavigator
void makeHits(const G4FastTrack &fastTrack)
G4bool ModelTrigger(const G4FastTrack &)