CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups 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 //
6 
9 
10 #include "G4Electron.hh"
11 #include "G4EventManager.hh"
12 #include "G4FastSimulationManager.hh"
13 #include "G4LogicalVolume.hh"
14 #include "G4Positron.hh"
15 #include "G4TouchableHandle.hh"
16 #include "G4TransportationManager.hh"
17 #include "G4VPhysicalVolume.hh"
18 #include "G4VProcess.hh"
19 #include "G4VSensitiveDetector.hh"
20 
21 using namespace CLHEP;
22 
23 GflashEMShowerModel::GflashEMShowerModel(const G4String &modelName,
24  G4Envelope *envelope,
25  const edm::ParameterSet &parSet)
26  : G4VFastSimulationModel(modelName, envelope), theParSet(parSet) {
27  theProfile = new GflashEMShowerProfile(parSet);
28  theRegion = reinterpret_cast<G4Region *>(envelope);
29 
30  theGflashStep = new G4Step();
31  theGflashTouchableHandle = new G4TouchableHistory();
32  theGflashNavigator = new G4Navigator();
33 }
34 
35 // -----------------------------------------------------------------------------------
36 
38  delete theProfile;
39  delete theGflashStep;
40 }
41 
42 G4bool GflashEMShowerModel::IsApplicable(const G4ParticleDefinition &particleType) {
43  return (&particleType == G4Electron::ElectronDefinition() || &particleType == G4Positron::PositronDefinition());
44 }
45 
46 // -----------------------------------------------------------------------------------
47 G4bool GflashEMShowerModel::ModelTrigger(const G4FastTrack &fastTrack) {
48  // Mininum energy cutoff to parameterize
49  if (fastTrack.GetPrimaryTrack()->GetKineticEnergy() < GeV) {
50  return false;
51  }
52  if (excludeDetectorRegion(fastTrack)) {
53  return false;
54  }
55 
56  // This will be changed accordingly when the way
57  // dealing with CaloRegion changes later.
58  const G4VPhysicalVolume *pCurrentVolume = (fastTrack.GetPrimaryTrack()->GetTouchable())->GetVolume();
59  if (pCurrentVolume == nullptr) {
60  return false;
61  }
62 
63  const G4LogicalVolume *lv = pCurrentVolume->GetLogicalVolume();
64  if (lv->GetRegion() != theRegion) {
65  return false;
66  }
67  return true;
68 }
69 
70 // -----------------------------------------------------------------------------------
71 void GflashEMShowerModel::DoIt(const G4FastTrack &fastTrack, G4FastStep &fastStep) {
72  // Kill the parameterised particle:
73  fastStep.KillPrimaryTrack();
74  fastStep.ProposePrimaryTrackPathLength(0.0);
75 
76  // input variables for GflashEMShowerProfile with showerType = 1,5 (shower
77  // 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  std::vector<GflashHit> &gflashHitList = theProfile->getGflashHitList();
95 
96  theGflashStep->SetTrack(const_cast<G4Track *>(fastTrack.GetPrimaryTrack()));
97 
98  theGflashStep->GetPostStepPoint()->SetProcessDefinedStep(
99  const_cast<G4VProcess *>(fastTrack.GetPrimaryTrack()->GetStep()->GetPostStepPoint()->GetProcessDefinedStep()));
100  theGflashNavigator->SetWorldVolume(
101  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  // put touchable for each hit so that touchable history keeps track of each
108  // step.
109  theGflashNavigator->LocateGlobalPointAndUpdateTouchableHandle(
110  spotIter->getPosition(), G4ThreeVector(0, 0, 0), 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)
118  continue;
119 
120  G4LogicalVolume *lv = aCurrentVolume->GetLogicalVolume();
121  if (lv->GetRegion() != theRegion)
122  continue;
123 
124  theGflashStep->GetPreStepPoint()->SetSensitiveDetector(aCurrentVolume->GetLogicalVolume()->GetSensitiveDetector());
125  G4VSensitiveDetector *aSensitive = theGflashStep->GetPreStepPoint()->GetSensitiveDetector();
126 
127  if (aSensitive == nullptr)
128  continue;
129 
130  theGflashStep->SetTotalEnergyDeposit(spotIter->getEnergy());
131  aSensitive->Hit(theGflashStep);
132  }
133 }
134 
135 void GflashEMShowerModel::updateGflashStep(const G4ThreeVector &spotPosition, G4double timeGlobal) {
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  // exclude regions where geometry are complicated
145  //+- one supermodule around the EB/EE boundary: 1.479 +- 0.0174*5
146  G4double eta = fastTrack.GetPrimaryTrack()->GetPosition().pseudoRapidity();
147  return (std::fabs(eta) > 1.392 && std::fabs(eta) < 1.566);
148 }
std::vector< GflashHit > & getGflashHitList()
G4bool excludeDetectorRegion(const G4FastTrack &fastTrack)
G4TouchableHandle theGflashTouchableHandle
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:289
G4Navigator * theGflashNavigator
void makeHits(const G4FastTrack &fastTrack)
G4bool ModelTrigger(const G4FastTrack &) override