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