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 // V.Ivanchenko rename the class, cleanup, and move
5 // to SimG4Core/Application - 2012/08/14
6 //
9 
12 
13 #include "G4Electron.hh"
14 #include "G4Positron.hh"
15 #include "G4VProcess.hh"
16 #include "G4VPhysicalVolume.hh"
17 #include "G4LogicalVolume.hh"
18 #include "G4TransportationManager.hh"
19 #include "G4EventManager.hh"
20 #include "G4FastSimulationManager.hh"
21 #include "G4TouchableHandle.hh"
22 #include "G4VSensitiveDetector.hh"
23 #include <CLHEP/Units/SystemOfUnits.h>
24 using CLHEP::cm;
25 using CLHEP::GeV;
26 
28  G4Envelope* envelope,
29  const edm::ParameterSet& parSet)
30  : G4VFastSimulationModel(modelName, envelope), theParSet(parSet) {
31  theWatcherOn = parSet.getParameter<bool>("watcherOn");
32 
33  theProfile = new GflashEMShowerProfile(parSet);
34  theRegion = reinterpret_cast<G4Region*>(envelope);
35 
36  theGflashStep = new G4Step();
37  theGflashTouchableHandle = new G4TouchableHistory();
38  theGflashNavigator = new G4Navigator();
39 }
40 
41 // --------------------------------------------------------------------------
42 
44  delete theProfile;
45  delete theGflashStep;
46 }
47 
48 G4bool GFlashEMShowerModel::IsApplicable(const G4ParticleDefinition& particleType) {
49  return (&particleType == G4Electron::Electron() || &particleType == G4Positron::Positron());
50 }
51 
52 // --------------------------------------------------------------------------
53 G4bool GFlashEMShowerModel::ModelTrigger(const G4FastTrack& fastTrack) {
54  // Mininum energy cutoff to parameterize
55  if (fastTrack.GetPrimaryTrack()->GetKineticEnergy() < Gflash::energyCutOff) {
56  return false;
57  }
58  if (excludeDetectorRegion(fastTrack)) {
59  return false;
60  }
61 
62  // This will be changed accordingly when the way
63  // dealing with CaloRegion changes later.
64  const G4TouchableHistory* touch = static_cast<const G4TouchableHistory*>(fastTrack.GetPrimaryTrack()->GetTouchable());
65  const G4VPhysicalVolume* pCurrentVolume = touch->GetVolume();
66  if (pCurrentVolume == nullptr) {
67  return false;
68  }
69 
70  const G4LogicalVolume* lv = pCurrentVolume->GetLogicalVolume();
71  if (lv->GetRegion() != theRegion) {
72  return false;
73  }
74  return true;
75 }
76 
77 // ---------------------------------------------------------------------------
78 void GFlashEMShowerModel::DoIt(const G4FastTrack& fastTrack, G4FastStep& fastStep) {
79  // Kill the parameterised particle:
80  fastStep.KillPrimaryTrack();
81  fastStep.ProposePrimaryTrackPathLength(0.0);
82 
83  // Input variables for GFlashEMShowerProfile with showerType = 1,5
84  // Shower starts inside crystals
85  G4double energy = fastTrack.GetPrimaryTrack()->GetKineticEnergy() / GeV;
86  G4double globalTime = fastTrack.GetPrimaryTrack()->GetStep()->GetPostStepPoint()->GetGlobalTime();
87  G4double charge = fastTrack.GetPrimaryTrack()->GetStep()->GetPreStepPoint()->GetCharge();
88  G4ThreeVector position = fastTrack.GetPrimaryTrack()->GetPosition() / cm;
89  G4ThreeVector momentum = fastTrack.GetPrimaryTrack()->GetMomentum() / GeV;
90  G4int showerType = Gflash::findShowerType(position);
91 
92  // Do actual parameterization
93  // The result of parameterization is gflashHitList
94  theProfile->initialize(showerType, energy, globalTime, charge, position, momentum);
96 
97  // Make hits
98  makeHits(fastTrack);
99 }
100 
101 // ---------------------------------------------------------------------------
102 void GFlashEMShowerModel::makeHits(const G4FastTrack& fastTrack) {
103  std::vector<GflashHit>& gflashHitList = theProfile->getGflashHitList();
104 
105  theGflashStep->SetTrack(const_cast<G4Track*>(fastTrack.GetPrimaryTrack()));
106 
107  theGflashStep->GetPostStepPoint()->SetProcessDefinedStep(
108  const_cast<G4VProcess*>(fastTrack.GetPrimaryTrack()->GetStep()->GetPostStepPoint()->GetProcessDefinedStep()));
109  theGflashNavigator->SetWorldVolume(
110  G4TransportationManager::GetTransportationManager()->GetNavigatorForTracking()->GetWorldVolume());
111 
112  std::vector<GflashHit>::const_iterator spotIter = gflashHitList.begin();
113  std::vector<GflashHit>::const_iterator spotIterEnd = gflashHitList.end();
114 
115  for (; spotIter != spotIterEnd; ++spotIter) {
116  // Put touchable for each hit so that touchable history
117  // keeps track of each step.
118  theGflashNavigator->LocateGlobalPointAndUpdateTouchableHandle(
119  spotIter->getPosition(), G4ThreeVector(0, 0, 0), theGflashTouchableHandle, false);
120  updateGflashStep(spotIter->getPosition(), spotIter->getTime());
121 
122  // If there is a watcher defined in a job and the flag is turned on
123  if (theWatcherOn) {
124  SteppingAction* userSteppingAction = (SteppingAction*)G4EventManager::GetEventManager()->GetUserSteppingAction();
125  userSteppingAction->m_g4StepSignal(theGflashStep);
126  }
127 
128  // Send G4Step information to Hit/Digi if the volume is sensitive
129  // Copied from G4SteppingManager.cc
130 
131  G4VPhysicalVolume* aCurrentVolume = theGflashStep->GetPreStepPoint()->GetPhysicalVolume();
132  if (aCurrentVolume == nullptr) {
133  continue;
134  }
135 
136  G4LogicalVolume* lv = aCurrentVolume->GetLogicalVolume();
137  if (lv->GetRegion() != theRegion) {
138  continue;
139  }
140 
141  theGflashStep->GetPreStepPoint()->SetSensitiveDetector(aCurrentVolume->GetLogicalVolume()->GetSensitiveDetector());
142  G4VSensitiveDetector* aSensitive = theGflashStep->GetPreStepPoint()->GetSensitiveDetector();
143 
144  if (aSensitive == nullptr) {
145  continue;
146  }
147 
148  theGflashStep->SetTotalEnergyDeposit(spotIter->getEnergy());
149  aSensitive->Hit(theGflashStep);
150  }
151 }
152 
153 // ---------------------------------------------------------------------------
154 void GFlashEMShowerModel::updateGflashStep(const G4ThreeVector& spotPosition, G4double timeGlobal) {
155  theGflashStep->GetPostStepPoint()->SetGlobalTime(timeGlobal);
156  theGflashStep->GetPreStepPoint()->SetPosition(spotPosition);
157  theGflashStep->GetPostStepPoint()->SetPosition(spotPosition);
158  theGflashStep->GetPreStepPoint()->SetTouchableHandle(theGflashTouchableHandle);
159 }
160 
161 // ---------------------------------------------------------------------------
162 G4bool GFlashEMShowerModel::excludeDetectorRegion(const G4FastTrack& fastTrack) {
163  //exclude regions where geometry are complicated
164  //+- one supermodule around the EB/EE boundary: 1.479 +- 0.0174*5
165  G4double eta = fastTrack.GetPrimaryTrack()->GetPosition().pseudoRapidity();
166  return (std::fabs(eta) > 1.392 && std::fabs(eta) < 1.566);
167 }
168 
169 // ---------------------------------------------------------------------------
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
void makeHits(const G4FastTrack &fastTrack)
std::vector< GflashHit > & getGflashHitList()
G4TouchableHandle theGflashTouchableHandle
SimActivityRegistry::G4StepSignal m_g4StepSignal
void DoIt(const G4FastTrack &, G4FastStep &) override
GflashEMShowerProfile * theProfile
void initialize(int showerType, double energy, double globalTime, double charge, Gflash3Vector &position, Gflash3Vector &momentum)
const double energyCutOff
G4bool IsApplicable(const G4ParticleDefinition &) override
int findShowerType(const Gflash3Vector &position)
GFlashEMShowerModel(const G4String &name, G4Envelope *env, const edm::ParameterSet &parSet)
G4bool ModelTrigger(const G4FastTrack &) override
static int position[264][3]
Definition: ReadPGInfo.cc:289
void updateGflashStep(const G4ThreeVector &position, G4double time)
G4bool excludeDetectorRegion(const G4FastTrack &fastTrack)
G4Navigator * theGflashNavigator