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 
29  theProfile = new GflashEMShowerProfile(parSet);
30  theRegion = const_cast<const G4Region*>(envelope);
31 
32  theGflashStep = new G4Step();
33  theGflashTouchableHandle = new G4TouchableHistory();
34  theGflashNavigator = new G4Navigator();
35 
36 }
G4TouchableHandle theGflashTouchableHandle
edm::ParameterSet theParSet
const G4Region * theRegion
GflashEMShowerProfile * theProfile
G4Navigator * theGflashNavigator
GflashEMShowerModel::~GflashEMShowerModel ( )
override

Definition at line 40 of file GflashEMShowerModel.cc.

References theGflashStep, and theProfile.

40  {
41 
42  delete theProfile;
43  delete theGflashStep;
44 }
GflashEMShowerProfile * theProfile

Member Function Documentation

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

Definition at line 71 of file GflashEMShowerModel.cc.

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

71  {
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 starts inside crystals)
78  G4double energy = fastTrack.GetPrimaryTrack()->GetKineticEnergy()/GeV;
79  G4double globalTime = fastTrack.GetPrimaryTrack()->GetStep()->GetPostStepPoint()->GetGlobalTime();
80  G4double charge = fastTrack.GetPrimaryTrack()->GetStep()->GetPreStepPoint()->GetCharge();
81  G4ThreeVector position = fastTrack.GetPrimaryTrack()->GetPosition() / cm;
82  G4ThreeVector momentum = fastTrack.GetPrimaryTrack()->GetMomentum()/GeV;
83  G4int showerType = Gflash::findShowerType(position);
84 
85  // Do actual parameterization. The result of parameterization is gflashHitList
86  theProfile->initialize(showerType,energy,globalTime,charge,position,momentum);
88 
89  //make hits
90  makeHits(fastTrack);
91 }
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 143 of file GflashEMShowerModel.cc.

References PVValHelper::eta.

Referenced by ModelTrigger().

143  {
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) { return true; }
151 
152  return isExcluded;
153 }
G4bool GflashEMShowerModel::IsApplicable ( const G4ParticleDefinition &  particleType)
override

Definition at line 46 of file GflashEMShowerModel.cc.

46  {
47 
48  return ( &particleType == G4Electron::ElectronDefinition() ||
49  &particleType == G4Positron::PositronDefinition() );
50 
51 }
void GflashEMShowerModel::makeHits ( const G4FastTrack &  fastTrack)
private

Definition at line 93 of file GflashEMShowerModel.cc.

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

Referenced by DoIt().

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

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

54  {
55 
56  // Mininum energy cutoff to parameterize
57  if(fastTrack.GetPrimaryTrack()->GetKineticEnergy() < GeV) { return false; }
58  if(excludeDetectorRegion(fastTrack)) { return false; }
59 
60  // This will be changed accordingly when the way
61  // dealing with CaloRegion changes later.
62  G4VPhysicalVolume* pCurrentVolume = (fastTrack.GetPrimaryTrack()->GetTouchable())->GetVolume();
63  if( pCurrentVolume == nullptr) { return false; }
64 
65  G4LogicalVolume* lv = pCurrentVolume->GetLogicalVolume();
66  if(lv->GetRegion() != theRegion) { return false; }
67  return true;
68 }
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 133 of file GflashEMShowerModel.cc.

References theGflashStep, and theGflashTouchableHandle.

Referenced by makeHits().

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

Member Data Documentation

G4Navigator* GflashEMShowerModel::theGflashNavigator
private

Definition at line 54 of file GflashEMShowerModel.h.

Referenced by GflashEMShowerModel(), and makeHits().

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

Definition at line 55 of file GflashEMShowerModel.h.

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

edm::ParameterSet GflashEMShowerModel::theParSet
private

Definition at line 47 of file GflashEMShowerModel.h.

GflashEMShowerProfile* GflashEMShowerModel::theProfile
private

Definition at line 49 of file GflashEMShowerModel.h.

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

const G4Region* GflashEMShowerModel::theRegion
private

Definition at line 51 of file GflashEMShowerModel.h.

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