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  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  return true;
70 }
71 
72 // -----------------------------------------------------------------------------------
73 void GflashEMShowerModel::DoIt(const G4FastTrack& fastTrack, G4FastStep& fastStep) {
74 
75  // Kill the parameterised particle:
76  fastStep.KillPrimaryTrack();
77  fastStep.ProposePrimaryTrackPathLength(0.0);
78 
79  //input variables for GflashEMShowerProfile with showerType = 1,5 (shower starts inside crystals)
80  G4double energy = fastTrack.GetPrimaryTrack()->GetKineticEnergy()/GeV;
81  G4double globalTime = fastTrack.GetPrimaryTrack()->GetStep()->GetPostStepPoint()->GetGlobalTime();
82  G4double charge = fastTrack.GetPrimaryTrack()->GetStep()->GetPreStepPoint()->GetCharge();
83  G4ThreeVector position = fastTrack.GetPrimaryTrack()->GetPosition() / cm;
84  G4ThreeVector momentum = fastTrack.GetPrimaryTrack()->GetMomentum()/GeV;
85  G4int showerType = Gflash::findShowerType(position);
86 
87  // Do actual parameterization. The result of parameterization is gflashHitList
88  theProfile->initialize(showerType,energy,globalTime,charge,position,momentum);
90 
91  //make hits
92  makeHits(fastTrack);
93 }
94 
95 void GflashEMShowerModel::makeHits(const G4FastTrack& fastTrack) {
96 
97  std::vector<GflashHit>& gflashHitList = theProfile->getGflashHitList();
98 
99  theGflashStep->SetTrack(const_cast<G4Track*>(fastTrack.GetPrimaryTrack()));
100 
101  theGflashStep->GetPostStepPoint()->SetProcessDefinedStep(const_cast<G4VProcess*>
102  (fastTrack.GetPrimaryTrack()->GetStep()->GetPostStepPoint()->GetProcessDefinedStep()));
103  theGflashNavigator->SetWorldVolume(G4TransportationManager::GetTransportationManager()->GetNavigatorForTracking()->GetWorldVolume());
104 
105  std::vector<GflashHit>::const_iterator spotIter = gflashHitList.begin();
106  std::vector<GflashHit>::const_iterator spotIterEnd = gflashHitList.end();
107 
108  for( ; spotIter != spotIterEnd; spotIter++){
109 
110  //put touchable for each hit so that touchable history keeps track of each step.
111  theGflashNavigator->LocateGlobalPointAndUpdateTouchableHandle(spotIter->getPosition(),G4ThreeVector(0,0,0),
112  theGflashTouchableHandle, false);
113  updateGflashStep(spotIter->getPosition(),spotIter->getTime());
114 
115  // Send G4Step information to Hit/Digi if the volume is sensitive
116  // Copied from G4SteppingManager.cc
117 
118  G4VPhysicalVolume* aCurrentVolume = theGflashStep->GetPreStepPoint()->GetPhysicalVolume();
119  if( aCurrentVolume == 0 ) continue;
120 
121  G4LogicalVolume* lv = aCurrentVolume->GetLogicalVolume();
122  if(lv->GetRegion() != theRegion) continue;
123 
124  theGflashStep->GetPreStepPoint()->SetSensitiveDetector(aCurrentVolume->GetLogicalVolume()->GetSensitiveDetector());
125  G4VSensitiveDetector* aSensitive = theGflashStep->GetPreStepPoint()->GetSensitiveDetector();
126 
127  if( aSensitive == 0 ) continue;
128 
129  theGflashStep->SetTotalEnergyDeposit(spotIter->getEnergy());
130  aSensitive->Hit(theGflashStep);
131  }
132 
133 }
134 
135 void GflashEMShowerModel::updateGflashStep(const G4ThreeVector& spotPosition,
136  G4double timeGlobal)
137 {
138  theGflashStep->GetPostStepPoint()->SetGlobalTime(timeGlobal);
139  theGflashStep->GetPreStepPoint()->SetPosition(spotPosition);
140  theGflashStep->GetPostStepPoint()->SetPosition(spotPosition);
141  theGflashStep->GetPreStepPoint()->SetTouchableHandle(theGflashTouchableHandle);
142 }
143 
144 // -----------------------------------------------------------------------------------
145 G4bool GflashEMShowerModel::excludeDetectorRegion(const G4FastTrack& fastTrack) {
146 
147  G4bool isExcluded=false;
148 
149  //exclude regions where geometry are complicated
150  //+- one supermodule around the EB/EE boundary: 1.479 +- 0.0174*5
151  G4double eta = fastTrack.GetPrimaryTrack()->GetPosition().pseudoRapidity() ;
152  if(std::fabs(eta) > 1.392 && std::fabs(eta) < 1.566) { return true; }
153 
154  return isExcluded;
155 }
156 
157 /*
158 G4int GflashEMShowerModel::findShowerType(const G4FastTrack& fastTrack)
159 {
160  // Initialization of longitudinal and lateral parameters for
161  // hadronic showers. Simulation of the intrinsic fluctuations
162 
163  // type of hadron showers subject to the shower starting point (ssp)
164  // showerType = -1 : default (invalid)
165  // showerType = 0 : ssp before EBRY (barrel crystal)
166  // showerType = 1 : ssp inside EBRY
167  // showerType = 2 : ssp after EBRY before HB
168  // showerType = 3 : ssp inside HB
169  // showerType = 4 : ssp before EFRY (endcap crystal)
170  // showerType = 5 : ssp inside EFRY
171  // showerType = 6 : ssp after EFRY before HE
172  // showerType = 7 : ssp inside HE
173 
174  G4TouchableHistory* touch = (G4TouchableHistory*)(fastTrack.GetPrimaryTrack()->GetTouchable());
175  G4LogicalVolume* lv = touch->GetVolume()->GetLogicalVolume();
176 
177  std::size_t pos1 = lv->GetName().find("EBRY");
178  std::size_t pos11 = lv->GetName().find("EWAL");
179  std::size_t pos12 = lv->GetName().find("EWRA");
180  std::size_t pos2 = lv->GetName().find("EFRY");
181 
182  G4ThreeVector position = fastTrack.GetPrimaryTrack()->GetPosition()/cm;
183  Gflash::CalorimeterNumber kCalor = Gflash::getCalorimeterNumber(position);
184 
185  G4int showerType = -1;
186 
187  //central
188  if (kCalor == Gflash::kESPM || kCalor == Gflash::kHB ) {
189 
190  G4double posRho = position.getRho();
191 
192  if(pos1 != std::string::npos || pos11 != std::string::npos || pos12 != std::string::npos ) {
193  showerType = 1;
194  }
195  else {
196  if(kCalor == Gflash::kESPM) {
197  showerType = 2;
198  if( posRho < Gflash::Rmin[Gflash::kESPM]+ Gflash::ROffCrystalEB ) showerType = 0;
199  }
200  else showerType = 3;
201  }
202 
203  }
204  //forward
205  else if (kCalor == Gflash::kENCA || kCalor == Gflash::kHE) {
206  if(pos2 != std::string::npos) {
207  showerType = 5;
208  }
209  else {
210  if(kCalor == Gflash::kENCA) {
211  showerType = 6;
212  if(fabs(position.getZ()) < Gflash::Zmin[Gflash::kENCA] + Gflash::ZOffCrystalEE) showerType = 4;
213  }
214  else showerType = 7;
215  }
216  //@@@need z-dependent correction on the mean energy reponse
217  }
218 
219  return showerType;
220 }
221 */
const double GeV
Definition: MathUtil.h:16
std::vector< GflashHit > & getGflashHitList()
void DoIt(const G4FastTrack &, G4FastStep &)
G4bool excludeDetectorRegion(const G4FastTrack &fastTrack)
G4TouchableHandle theGflashTouchableHandle
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 &)