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 = const_cast<const 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  //std::cout << "GFlashEMShowerModel::ModelTrigger: LV "
73  // << lv->GetRegion()->GetName() << std::endl;
74 
75  // The parameterization starts inside crystals
76  //std::size_t pos1 = lv->GetName().find("EBRY");
77  //std::size_t pos2 = lv->GetName().find("EFRY");
78 
79  //std::size_t pos3 = lv->GetName().find("HVQ");
80  //std::size_t pos4 = lv->GetName().find("HF");
81  //if(pos1 == std::string::npos && pos2 == std::string::npos &&
82  // pos3 == std::string::npos && pos4 == std::string::npos) return false;
83 
84  return true;
85 }
86 
87 // ---------------------------------------------------------------------------
88 void GFlashEMShowerModel::DoIt(const G4FastTrack& fastTrack, G4FastStep& fastStep) {
89  // Kill the parameterised particle:
90  fastStep.KillPrimaryTrack();
91  fastStep.ProposePrimaryTrackPathLength(0.0);
92 
93  // Input variables for GFlashEMShowerProfile with showerType = 1,5
94  // Shower starts inside crystals
95  G4double energy = fastTrack.GetPrimaryTrack()->GetKineticEnergy() / GeV;
96  G4double globalTime = fastTrack.GetPrimaryTrack()->GetStep()->GetPostStepPoint()->GetGlobalTime();
97  G4double charge = fastTrack.GetPrimaryTrack()->GetStep()->GetPreStepPoint()->GetCharge();
98  G4ThreeVector position = fastTrack.GetPrimaryTrack()->GetPosition() / cm;
99  G4ThreeVector momentum = fastTrack.GetPrimaryTrack()->GetMomentum() / GeV;
100  G4int showerType = Gflash::findShowerType(position);
101 
102  // Do actual parameterization
103  // The result of parameterization is gflashHitList
104  theProfile->initialize(showerType, energy, globalTime, charge, position, momentum);
106 
107  // Make hits
108  makeHits(fastTrack);
109 }
110 
111 // ---------------------------------------------------------------------------
112 void GFlashEMShowerModel::makeHits(const G4FastTrack& fastTrack) {
113  std::vector<GflashHit>& gflashHitList = theProfile->getGflashHitList();
114 
115  theGflashStep->SetTrack(const_cast<G4Track*>(fastTrack.GetPrimaryTrack()));
116 
117  theGflashStep->GetPostStepPoint()->SetProcessDefinedStep(
118  const_cast<G4VProcess*>(fastTrack.GetPrimaryTrack()->GetStep()->GetPostStepPoint()->GetProcessDefinedStep()));
119  theGflashNavigator->SetWorldVolume(
120  G4TransportationManager::GetTransportationManager()->GetNavigatorForTracking()->GetWorldVolume());
121 
122  std::vector<GflashHit>::const_iterator spotIter = gflashHitList.begin();
123  std::vector<GflashHit>::const_iterator spotIterEnd = gflashHitList.end();
124 
125  for (; spotIter != spotIterEnd; spotIter++) {
126  // Put touchable for each hit so that touchable history
127  // keeps track of each step.
128  theGflashNavigator->LocateGlobalPointAndUpdateTouchableHandle(
129  spotIter->getPosition(), G4ThreeVector(0, 0, 0), theGflashTouchableHandle, false);
130  updateGflashStep(spotIter->getPosition(), spotIter->getTime());
131 
132  // If there is a watcher defined in a job and the flag is turned on
133  if (theWatcherOn) {
134  SteppingAction* userSteppingAction = (SteppingAction*)G4EventManager::GetEventManager()->GetUserSteppingAction();
135  userSteppingAction->m_g4StepSignal(theGflashStep);
136  }
137 
138  // Send G4Step information to Hit/Digi if the volume is sensitive
139  // Copied from G4SteppingManager.cc
140 
141  G4VPhysicalVolume* aCurrentVolume = theGflashStep->GetPreStepPoint()->GetPhysicalVolume();
142  if (aCurrentVolume == nullptr) {
143  continue;
144  }
145 
146  G4LogicalVolume* lv = aCurrentVolume->GetLogicalVolume();
147  if (lv->GetRegion() != theRegion) {
148  continue;
149  }
150 
151  theGflashStep->GetPreStepPoint()->SetSensitiveDetector(aCurrentVolume->GetLogicalVolume()->GetSensitiveDetector());
152  G4VSensitiveDetector* aSensitive = theGflashStep->GetPreStepPoint()->GetSensitiveDetector();
153 
154  if (aSensitive == nullptr) {
155  continue;
156  }
157 
158  theGflashStep->SetTotalEnergyDeposit(spotIter->getEnergy());
159  aSensitive->Hit(theGflashStep);
160  }
161 }
162 
163 // ---------------------------------------------------------------------------
164 void GFlashEMShowerModel::updateGflashStep(const G4ThreeVector& spotPosition, G4double timeGlobal) {
165  theGflashStep->GetPostStepPoint()->SetGlobalTime(timeGlobal);
166  theGflashStep->GetPreStepPoint()->SetPosition(spotPosition);
167  theGflashStep->GetPostStepPoint()->SetPosition(spotPosition);
168  theGflashStep->GetPreStepPoint()->SetTouchableHandle(theGflashTouchableHandle);
169 }
170 
171 // ---------------------------------------------------------------------------
172 G4bool GFlashEMShowerModel::excludeDetectorRegion(const G4FastTrack& fastTrack) {
173  G4bool isExcluded = false;
174 
175  //exclude regions where geometry are complicated
176  //+- one supermodule around the EB/EE boundary: 1.479 +- 0.0174*5
177  G4double eta = fastTrack.GetPrimaryTrack()->GetPosition().pseudoRapidity();
178  if (std::fabs(eta) > 1.392 && std::fabs(eta) < 1.566) {
179  return true;
180  }
181 
182  return isExcluded;
183 }
184 
185 // ---------------------------------------------------------------------------
GflashEMShowerProfile.h
GflashHit.h
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
SteppingAction.h
GFlashEMShowerModel::DoIt
void DoIt(const G4FastTrack &, G4FastStep &) override
Definition: GFlashEMShowerModel.cc:88
SteppingAction
Definition: SteppingAction.h:31
GFlashEMShowerModel.h
GFlashEMShowerModel::theProfile
GflashEMShowerProfile * theProfile
Definition: GFlashEMShowerModel.h:47
Gflash::findShowerType
int findShowerType(const Gflash3Vector &position)
Definition: GflashNameSpace.cc:41
GFlashEMShowerModel::updateGflashStep
void updateGflashStep(const G4ThreeVector &position, G4double time)
Definition: GFlashEMShowerModel.cc:164
HLTEgPhaseIITestSequence_cff.modelName
modelName
Definition: HLTEgPhaseIITestSequence_cff.py:16
GflashEMShowerProfile
Definition: GflashEMShowerProfile.h:15
PVValHelper::eta
Definition: PVValidationHelpers.h:69
Gflash::energyCutOff
const double energyCutOff
Definition: GflashNameSpace.h:44
GFlashEMShowerModel::IsApplicable
G4bool IsApplicable(const G4ParticleDefinition &) override
Definition: GFlashEMShowerModel.cc:46
HCALHighEnergyHPDFilter_cfi.energy
energy
Definition: HCALHighEnergyHPDFilter_cfi.py:5
ALCARECOTkAlJpsiMuMu_cff.charge
charge
Definition: ALCARECOTkAlJpsiMuMu_cff.py:47
edm::ParameterSet
Definition: ParameterSet.h:47
GeV
const double GeV
Definition: MathUtil.h:16
GFlashEMShowerModel::theRegion
const G4Region * theRegion
Definition: GFlashEMShowerModel.h:49
position
static int position[264][3]
Definition: ReadPGInfo.cc:289
GFlashEMShowerModel::theGflashStep
G4Step * theGflashStep
Definition: GFlashEMShowerModel.h:51
nanoDQM_cff.Electron
Electron
Definition: nanoDQM_cff.py:62
SteppingAction::m_g4StepSignal
SimActivityRegistry::G4StepSignal m_g4StepSignal
Definition: SteppingAction.h:38
GFlashEMShowerModel::makeHits
void makeHits(const G4FastTrack &fastTrack)
Definition: GFlashEMShowerModel.cc:112
GFlashEMShowerModel::theGflashNavigator
G4Navigator * theGflashNavigator
Definition: GFlashEMShowerModel.h:52
GFlashEMShowerModel::excludeDetectorRegion
G4bool excludeDetectorRegion(const G4FastTrack &fastTrack)
Definition: GFlashEMShowerModel.cc:172
GFlashEMShowerModel::ModelTrigger
G4bool ModelTrigger(const G4FastTrack &) override
Definition: GFlashEMShowerModel.cc:51
GflashEMShowerProfile::getGflashHitList
std::vector< GflashHit > & getGflashHitList()
Definition: GflashEMShowerProfile.h:32
GFlashEMShowerModel::theGflashTouchableHandle
G4TouchableHandle theGflashTouchableHandle
Definition: GFlashEMShowerModel.h:53
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
PbPb_ZMuSkimMuonDPG_cff.particleType
particleType
Definition: PbPb_ZMuSkimMuonDPG_cff.py:27
GFlashEMShowerModel::GFlashEMShowerModel
GFlashEMShowerModel(const G4String &name, G4Envelope *env, const edm::ParameterSet &parSet)
Definition: GFlashEMShowerModel.cc:25
GFlashEMShowerModel::~GFlashEMShowerModel
~GFlashEMShowerModel() override
Definition: GFlashEMShowerModel.cc:41
GFlashEMShowerModel::theWatcherOn
bool theWatcherOn
Definition: GFlashEMShowerModel.h:45