CMS 3D CMS Logo

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