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