CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
GflashEMShowerModel Class Reference

#include <GflashEMShowerModel.h>

Inheritance diagram for GflashEMShowerModel:

Public Member Functions

void DoIt (const G4FastTrack &, G4FastStep &) override
 
 GflashEMShowerModel (const G4String &name, G4Envelope *env, const edm::ParameterSet &parSet)
 
G4bool IsApplicable (const G4ParticleDefinition &) override
 
G4bool ModelTrigger (const G4FastTrack &) override
 
 ~GflashEMShowerModel () override
 

Private Member Functions

G4bool excludeDetectorRegion (const G4FastTrack &fastTrack)
 
void makeHits (const G4FastTrack &fastTrack)
 
void updateGflashStep (const G4ThreeVector &position, G4double time)
 

Private Attributes

G4Navigator * theGflashNavigator
 
G4Step * theGflashStep
 
G4TouchableHandle theGflashTouchableHandle
 
edm::ParameterSet theParSet
 
GflashEMShowerProfiletheProfile
 
const G4Region * theRegion
 

Detailed Description

Definition at line 28 of file GflashEMShowerModel.h.

Constructor & Destructor Documentation

GflashEMShowerModel::GflashEMShowerModel ( const G4String &  name,
G4Envelope *  env,
const edm::ParameterSet parSet 
)

Definition at line 24 of file GflashEMShowerModel.cc.

References theGflashNavigator, theGflashStep, theGflashTouchableHandle, theProfile, and theRegion.

27  : G4VFastSimulationModel(modelName, envelope), theParSet(parSet) {
28  theProfile = new GflashEMShowerProfile(parSet);
29  theRegion = const_cast<const G4Region *>(envelope);
30 
31  theGflashStep = new G4Step();
32  theGflashTouchableHandle = new G4TouchableHistory();
33  theGflashNavigator = new G4Navigator();
34 }
G4TouchableHandle theGflashTouchableHandle
edm::ParameterSet theParSet
const G4Region * theRegion
GflashEMShowerProfile * theProfile
G4Navigator * theGflashNavigator
GflashEMShowerModel::~GflashEMShowerModel ( )
override

Definition at line 38 of file GflashEMShowerModel.cc.

References theGflashStep, and theProfile.

38  {
39  delete theProfile;
40  delete theGflashStep;
41 }
GflashEMShowerProfile * theProfile

Member Function Documentation

void GflashEMShowerModel::DoIt ( const G4FastTrack &  fastTrack,
G4FastStep &  fastStep 
)
override

Definition at line 72 of file GflashEMShowerModel.cc.

References ALCARECOTkAlJpsiMuMu_cff::charge, Gflash::findShowerType(), GeV, GflashEMShowerProfile::initialize(), makeHits(), GflashEMShowerProfile::parameterization(), position, and theProfile.

72  {
73  // Kill the parameterised particle:
74  fastStep.KillPrimaryTrack();
75  fastStep.ProposePrimaryTrackPathLength(0.0);
76 
77  // input variables for GflashEMShowerProfile with showerType = 1,5 (shower
78  // starts inside crystals)
79  G4double energy = fastTrack.GetPrimaryTrack()->GetKineticEnergy() / GeV;
80  G4double globalTime = fastTrack.GetPrimaryTrack()->GetStep()->GetPostStepPoint()->GetGlobalTime();
81  G4double charge = fastTrack.GetPrimaryTrack()->GetStep()->GetPreStepPoint()->GetCharge();
82  G4ThreeVector position = fastTrack.GetPrimaryTrack()->GetPosition() / cm;
83  G4ThreeVector momentum = fastTrack.GetPrimaryTrack()->GetMomentum() / GeV;
84  G4int showerType = Gflash::findShowerType(position);
85 
86  // Do actual parameterization. The result of parameterization is gflashHitList
87  theProfile->initialize(showerType, energy, globalTime, charge, position, momentum);
89 
90  // make hits
91  makeHits(fastTrack);
92 }
const double GeV
Definition: MathUtil.h:16
GflashEMShowerProfile * theProfile
void initialize(int showerType, double energy, double globalTime, double charge, Gflash3Vector &position, Gflash3Vector &momentum)
int findShowerType(const Gflash3Vector &position)
static int position[264][3]
Definition: ReadPGInfo.cc:509
void makeHits(const G4FastTrack &fastTrack)
G4bool GflashEMShowerModel::excludeDetectorRegion ( const G4FastTrack &  fastTrack)
private

Definition at line 144 of file GflashEMShowerModel.cc.

References PVValHelper::eta.

Referenced by ModelTrigger().

144  {
145  G4bool isExcluded = false;
146 
147  // exclude regions where geometry are complicated
148  //+- one supermodule around the EB/EE boundary: 1.479 +- 0.0174*5
149  G4double eta = fastTrack.GetPrimaryTrack()->GetPosition().pseudoRapidity();
150  if (std::fabs(eta) > 1.392 && std::fabs(eta) < 1.566) {
151  return true;
152  }
153 
154  return isExcluded;
155 }
G4bool GflashEMShowerModel::IsApplicable ( const G4ParticleDefinition &  particleType)
override

Definition at line 43 of file GflashEMShowerModel.cc.

43  {
44  return (&particleType == G4Electron::ElectronDefinition() || &particleType == G4Positron::PositronDefinition());
45 }
void GflashEMShowerModel::makeHits ( const G4FastTrack &  fastTrack)
private

Definition at line 94 of file GflashEMShowerModel.cc.

References GflashEMShowerProfile::getGflashHitList(), theGflashNavigator, theGflashStep, theGflashTouchableHandle, theProfile, theRegion, and updateGflashStep().

Referenced by DoIt().

94  {
95  std::vector<GflashHit> &gflashHitList = theProfile->getGflashHitList();
96 
97  theGflashStep->SetTrack(const_cast<G4Track *>(fastTrack.GetPrimaryTrack()));
98 
99  theGflashStep->GetPostStepPoint()->SetProcessDefinedStep(
100  const_cast<G4VProcess *>(fastTrack.GetPrimaryTrack()->GetStep()->GetPostStepPoint()->GetProcessDefinedStep()));
101  theGflashNavigator->SetWorldVolume(
102  G4TransportationManager::GetTransportationManager()->GetNavigatorForTracking()->GetWorldVolume());
103 
104  std::vector<GflashHit>::const_iterator spotIter = gflashHitList.begin();
105  std::vector<GflashHit>::const_iterator spotIterEnd = gflashHitList.end();
106 
107  for (; spotIter != spotIterEnd; spotIter++) {
108  // put touchable for each hit so that touchable history keeps track of each
109  // step.
110  theGflashNavigator->LocateGlobalPointAndUpdateTouchableHandle(
111  spotIter->getPosition(), G4ThreeVector(0, 0, 0), theGflashTouchableHandle, false);
112  updateGflashStep(spotIter->getPosition(), spotIter->getTime());
113 
114  // Send G4Step information to Hit/Digi if the volume is sensitive
115  // Copied from G4SteppingManager.cc
116 
117  G4VPhysicalVolume *aCurrentVolume = theGflashStep->GetPreStepPoint()->GetPhysicalVolume();
118  if (aCurrentVolume == nullptr)
119  continue;
120 
121  G4LogicalVolume *lv = aCurrentVolume->GetLogicalVolume();
122  if (lv->GetRegion() != theRegion)
123  continue;
124 
125  theGflashStep->GetPreStepPoint()->SetSensitiveDetector(aCurrentVolume->GetLogicalVolume()->GetSensitiveDetector());
126  G4VSensitiveDetector *aSensitive = theGflashStep->GetPreStepPoint()->GetSensitiveDetector();
127 
128  if (aSensitive == nullptr)
129  continue;
130 
131  theGflashStep->SetTotalEnergyDeposit(spotIter->getEnergy());
132  aSensitive->Hit(theGflashStep);
133  }
134 }
std::vector< GflashHit > & getGflashHitList()
G4TouchableHandle theGflashTouchableHandle
const G4Region * theRegion
GflashEMShowerProfile * theProfile
void updateGflashStep(const G4ThreeVector &position, G4double time)
G4Navigator * theGflashNavigator
G4bool GflashEMShowerModel::ModelTrigger ( const G4FastTrack &  fastTrack)
override

Definition at line 48 of file GflashEMShowerModel.cc.

References excludeDetectorRegion(), GetVolume(), GeV, and theRegion.

48  {
49  // Mininum energy cutoff to parameterize
50  if (fastTrack.GetPrimaryTrack()->GetKineticEnergy() < GeV) {
51  return false;
52  }
53  if (excludeDetectorRegion(fastTrack)) {
54  return false;
55  }
56 
57  // This will be changed accordingly when the way
58  // dealing with CaloRegion changes later.
59  G4VPhysicalVolume *pCurrentVolume = (fastTrack.GetPrimaryTrack()->GetTouchable())->GetVolume();
60  if (pCurrentVolume == nullptr) {
61  return false;
62  }
63 
64  G4LogicalVolume *lv = pCurrentVolume->GetLogicalVolume();
65  if (lv->GetRegion() != theRegion) {
66  return false;
67  }
68  return true;
69 }
const double GeV
Definition: MathUtil.h:16
G4bool excludeDetectorRegion(const G4FastTrack &fastTrack)
const G4Region * theRegion
static const G4LogicalVolume * GetVolume(const std::string &name)
void GflashEMShowerModel::updateGflashStep ( const G4ThreeVector &  position,
G4double  time 
)
private

Definition at line 136 of file GflashEMShowerModel.cc.

References theGflashStep, and theGflashTouchableHandle.

Referenced by makeHits().

136  {
137  theGflashStep->GetPostStepPoint()->SetGlobalTime(timeGlobal);
138  theGflashStep->GetPreStepPoint()->SetPosition(spotPosition);
139  theGflashStep->GetPostStepPoint()->SetPosition(spotPosition);
140  theGflashStep->GetPreStepPoint()->SetTouchableHandle(theGflashTouchableHandle);
141 }
G4TouchableHandle theGflashTouchableHandle

Member Data Documentation

G4Navigator* GflashEMShowerModel::theGflashNavigator
private

Definition at line 50 of file GflashEMShowerModel.h.

Referenced by GflashEMShowerModel(), and makeHits().

G4Step* GflashEMShowerModel::theGflashStep
private
G4TouchableHandle GflashEMShowerModel::theGflashTouchableHandle
private

Definition at line 51 of file GflashEMShowerModel.h.

Referenced by GflashEMShowerModel(), makeHits(), and updateGflashStep().

edm::ParameterSet GflashEMShowerModel::theParSet
private

Definition at line 43 of file GflashEMShowerModel.h.

GflashEMShowerProfile* GflashEMShowerModel::theProfile
private

Definition at line 45 of file GflashEMShowerModel.h.

Referenced by DoIt(), GflashEMShowerModel(), makeHits(), and ~GflashEMShowerModel().

const G4Region* GflashEMShowerModel::theRegion
private

Definition at line 47 of file GflashEMShowerModel.h.

Referenced by GflashEMShowerModel(), makeHits(), and ModelTrigger().