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 "G4EventManager.hh"
13 #include "G4FastSimulationManager.hh"
14 #include "G4LogicalVolume.hh"
15 #include "G4Positron.hh"
16 #include "G4TouchableHandle.hh"
17 #include "G4TransportationManager.hh"
18 #include "G4VPhysicalVolume.hh"
19 #include "G4VProcess.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  theProfile = new GflashEMShowerProfile(parSet);
29  theRegion = const_cast<const G4Region *>(envelope);
30 
31  theGflashStep = new G4Step();
32  theGflashTouchableHandle = new G4TouchableHistory();
33  theGflashNavigator = new G4Navigator();
34 }
35 
36 // -----------------------------------------------------------------------------------
37 
39  delete theProfile;
40  delete theGflashStep;
41 }
42 
43 G4bool GflashEMShowerModel::IsApplicable(const G4ParticleDefinition &particleType) {
44  return (&particleType == G4Electron::ElectronDefinition() || &particleType == G4Positron::PositronDefinition());
45 }
46 
47 // -----------------------------------------------------------------------------------
48 G4bool GflashEMShowerModel::ModelTrigger(const G4FastTrack &fastTrack) {
49  // Mininum energy cutoff to parameterize
50  if (fastTrack.GetPrimaryTrack()->GetKineticEnergy() < GeV) {
51  return false;
52  }
53  if (excludeDetectorRegion(fastTrack)) {
54  return false;
55  }
56 
57  // This will be changed accordingly when the way
58  // dealing with CaloRegion changes later.
59  G4VPhysicalVolume *pCurrentVolume = (fastTrack.GetPrimaryTrack()->GetTouchable())->GetVolume();
60  if (pCurrentVolume == nullptr) {
61  return false;
62  }
63 
64  G4LogicalVolume *lv = pCurrentVolume->GetLogicalVolume();
65  if (lv->GetRegion() != theRegion) {
66  return false;
67  }
68  return true;
69 }
70 
71 // -----------------------------------------------------------------------------------
72 void GflashEMShowerModel::DoIt(const G4FastTrack &fastTrack, G4FastStep &fastStep) {
73  // Kill the parameterised particle:
74  fastStep.KillPrimaryTrack();
75  fastStep.ProposePrimaryTrackPathLength(0.0);
76 
77  // input variables for GflashEMShowerProfile with showerType = 1,5 (shower
78  // starts inside crystals)
79  G4double energy = fastTrack.GetPrimaryTrack()->GetKineticEnergy() / GeV;
80  G4double globalTime = fastTrack.GetPrimaryTrack()->GetStep()->GetPostStepPoint()->GetGlobalTime();
81  G4double charge = fastTrack.GetPrimaryTrack()->GetStep()->GetPreStepPoint()->GetCharge();
82  G4ThreeVector position = fastTrack.GetPrimaryTrack()->GetPosition() / cm;
83  G4ThreeVector momentum = fastTrack.GetPrimaryTrack()->GetMomentum() / GeV;
84  G4int showerType = Gflash::findShowerType(position);
85 
86  // Do actual parameterization. The result of parameterization is gflashHitList
87  theProfile->initialize(showerType, energy, globalTime, charge, position, momentum);
89 
90  // make hits
91  makeHits(fastTrack);
92 }
93 
94 void GflashEMShowerModel::makeHits(const G4FastTrack &fastTrack) {
95  std::vector<GflashHit> &gflashHitList = theProfile->getGflashHitList();
96 
97  theGflashStep->SetTrack(const_cast<G4Track *>(fastTrack.GetPrimaryTrack()));
98 
99  theGflashStep->GetPostStepPoint()->SetProcessDefinedStep(
100  const_cast<G4VProcess *>(fastTrack.GetPrimaryTrack()->GetStep()->GetPostStepPoint()->GetProcessDefinedStep()));
101  theGflashNavigator->SetWorldVolume(
102  G4TransportationManager::GetTransportationManager()->GetNavigatorForTracking()->GetWorldVolume());
103 
104  std::vector<GflashHit>::const_iterator spotIter = gflashHitList.begin();
105  std::vector<GflashHit>::const_iterator spotIterEnd = gflashHitList.end();
106 
107  for (; spotIter != spotIterEnd; spotIter++) {
108  // put touchable for each hit so that touchable history keeps track of each
109  // step.
110  theGflashNavigator->LocateGlobalPointAndUpdateTouchableHandle(
111  spotIter->getPosition(), G4ThreeVector(0, 0, 0), theGflashTouchableHandle, false);
112  updateGflashStep(spotIter->getPosition(), spotIter->getTime());
113 
114  // Send G4Step information to Hit/Digi if the volume is sensitive
115  // Copied from G4SteppingManager.cc
116 
117  G4VPhysicalVolume *aCurrentVolume = theGflashStep->GetPreStepPoint()->GetPhysicalVolume();
118  if (aCurrentVolume == nullptr)
119  continue;
120 
121  G4LogicalVolume *lv = aCurrentVolume->GetLogicalVolume();
122  if (lv->GetRegion() != theRegion)
123  continue;
124 
125  theGflashStep->GetPreStepPoint()->SetSensitiveDetector(aCurrentVolume->GetLogicalVolume()->GetSensitiveDetector());
126  G4VSensitiveDetector *aSensitive = theGflashStep->GetPreStepPoint()->GetSensitiveDetector();
127 
128  if (aSensitive == nullptr)
129  continue;
130 
131  theGflashStep->SetTotalEnergyDeposit(spotIter->getEnergy());
132  aSensitive->Hit(theGflashStep);
133  }
134 }
135 
136 void GflashEMShowerModel::updateGflashStep(const G4ThreeVector &spotPosition, G4double timeGlobal) {
137  theGflashStep->GetPostStepPoint()->SetGlobalTime(timeGlobal);
138  theGflashStep->GetPreStepPoint()->SetPosition(spotPosition);
139  theGflashStep->GetPostStepPoint()->SetPosition(spotPosition);
140  theGflashStep->GetPreStepPoint()->SetTouchableHandle(theGflashTouchableHandle);
141 }
142 
143 // -----------------------------------------------------------------------------------
144 G4bool GflashEMShowerModel::excludeDetectorRegion(const G4FastTrack &fastTrack) {
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) {
151  return true;
152  }
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 =
175 (G4TouchableHistory*)(fastTrack.GetPrimaryTrack()->GetTouchable());
176  G4LogicalVolume* lv = touch->GetVolume()->GetLogicalVolume();
177 
178  std::size_t pos1 = lv->GetName().find("EBRY");
179  std::size_t pos11 = lv->GetName().find("EWAL");
180  std::size_t pos12 = lv->GetName().find("EWRA");
181  std::size_t pos2 = lv->GetName().find("EFRY");
182 
183  G4ThreeVector position = fastTrack.GetPrimaryTrack()->GetPosition()/cm;
184  Gflash::CalorimeterNumber kCalor = Gflash::getCalorimeterNumber(position);
185 
186  G4int showerType = -1;
187 
188  //central
189  if (kCalor == Gflash::kESPM || kCalor == Gflash::kHB ) {
190 
191  G4double posRho = position.getRho();
192 
193  if(pos1 != std::string::npos || pos11 != std::string::npos || pos12 !=
194 std::string::npos ) { showerType = 1;
195  }
196  else {
197  if(kCalor == Gflash::kESPM) {
198  showerType = 2;
199  if( posRho < Gflash::Rmin[Gflash::kESPM]+ Gflash::ROffCrystalEB )
200 showerType = 0;
201  }
202  else showerType = 3;
203  }
204 
205  }
206  //forward
207  else if (kCalor == Gflash::kENCA || kCalor == Gflash::kHE) {
208  if(pos2 != std::string::npos) {
209  showerType = 5;
210  }
211  else {
212  if(kCalor == Gflash::kENCA) {
213  showerType = 6;
214  if(fabs(position.getZ()) < Gflash::Zmin[Gflash::kENCA] +
215 Gflash::ZOffCrystalEE) showerType = 4;
216  }
217  else showerType = 7;
218  }
219  //@@@need z-dependent correction on the mean energy reponse
220  }
221 
222  return showerType;
223 }
224 */
GflashEMShowerProfile.h
GflashHit.h
GflashEMShowerModel::DoIt
void DoIt(const G4FastTrack &, G4FastStep &) override
Definition: GflashEMShowerModel.cc:72
GflashEMShowerProfile::initialize
void initialize(int showerType, double energy, double globalTime, double charge, Gflash3Vector &position, Gflash3Vector &momentum)
Definition: GflashEMShowerProfile.cc:37
GflashEMShowerProfile::parameterization
void parameterization()
Definition: GflashEMShowerProfile.cc:42
GflashEMShowerModel::~GflashEMShowerModel
~GflashEMShowerModel() override
Definition: GflashEMShowerModel.cc:38
SteppingAction.h
GflashEMShowerModel.h
Gflash::findShowerType
int findShowerType(const Gflash3Vector &position)
Definition: GflashNameSpace.cc:41
HLTEgPhaseIITestSequence_cff.modelName
modelName
Definition: HLTEgPhaseIITestSequence_cff.py:16
GflashEMShowerProfile
Definition: GflashEMShowerProfile.h:15
PVValHelper::eta
Definition: PVValidationHelpers.h:69
GflashEMShowerModel::theRegion
const G4Region * theRegion
Definition: GflashEMShowerModel.h:47
HCALHighEnergyHPDFilter_cfi.energy
energy
Definition: HCALHighEnergyHPDFilter_cfi.py:5
CLHEP
Definition: CocoaGlobals.h:27
ALCARECOTkAlJpsiMuMu_cff.charge
charge
Definition: ALCARECOTkAlJpsiMuMu_cff.py:47
edm::ParameterSet
Definition: ParameterSet.h:36
GeV
const double GeV
Definition: MathUtil.h:16
GetVolume
static const G4LogicalVolume * GetVolume(const std::string &name)
Definition: TrackingMaterialProducer.cc:41
GflashEMShowerModel::IsApplicable
G4bool IsApplicable(const G4ParticleDefinition &) override
Definition: GflashEMShowerModel.cc:43
position
static int position[264][3]
Definition: ReadPGInfo.cc:289
GflashEMShowerModel::excludeDetectorRegion
G4bool excludeDetectorRegion(const G4FastTrack &fastTrack)
Definition: GflashEMShowerModel.cc:144
GflashEMShowerModel::theGflashStep
G4Step * theGflashStep
Definition: GflashEMShowerModel.h:49
GflashEMShowerModel::ModelTrigger
G4bool ModelTrigger(const G4FastTrack &) override
Definition: GflashEMShowerModel.cc:48
GflashEMShowerModel::makeHits
void makeHits(const G4FastTrack &fastTrack)
Definition: GflashEMShowerModel.cc:94
GflashEMShowerProfile::getGflashHitList
std::vector< GflashHit > & getGflashHitList()
Definition: GflashEMShowerProfile.h:32
GflashEMShowerModel::theGflashNavigator
G4Navigator * theGflashNavigator
Definition: GflashEMShowerModel.h:50
GflashEMShowerModel::theGflashTouchableHandle
G4TouchableHandle theGflashTouchableHandle
Definition: GflashEMShowerModel.h:51
PbPb_ZMuSkimMuonDPG_cff.particleType
particleType
Definition: PbPb_ZMuSkimMuonDPG_cff.py:27
GflashEMShowerModel::updateGflashStep
void updateGflashStep(const G4ThreeVector &position, G4double time)
Definition: GflashEMShowerModel.cc:136
GflashEMShowerModel::GflashEMShowerModel
GflashEMShowerModel(const G4String &name, G4Envelope *env, const edm::ParameterSet &parSet)
Definition: GflashEMShowerModel.cc:24
GflashEMShowerModel::theProfile
GflashEMShowerProfile * theProfile
Definition: GflashEMShowerModel.h:45