CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 &)
 
 GflashEMShowerModel (const G4String &name, G4Envelope *env, const edm::ParameterSet &parSet)
 
G4bool IsApplicable (const G4ParticleDefinition &)
 
G4bool ModelTrigger (const G4FastTrack &)
 
virtual ~GflashEMShowerModel ()
 

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
 
bool theWatcherOn
 

Detailed Description

Definition at line 55 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  //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 }
G4TouchableHandle theGflashTouchableHandle
edm::ParameterSet theParSet
const G4Region * theRegion
GflashEMShowerProfile * theProfile
G4Navigator * theGflashNavigator
GflashEMShowerModel::~GflashEMShowerModel ( )
virtual

Definition at line 42 of file GflashEMShowerModel.cc.

References theGflashStep, and theProfile.

42  {
43 
44  delete theProfile;
45  delete theGflashStep;
46 }
GflashEMShowerProfile * theProfile

Member Function Documentation

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

Definition at line 90 of file GflashEMShowerModel.cc.

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

90  {
91 
92  // Kill the parameterised particle:
93  fastStep.KillPrimaryTrack();
94  fastStep.ProposePrimaryTrackPathLength(0.0);
95 
96  //input variables for GflashEMShowerProfile with showerType = 1,5 (shower starts inside crystals)
97  G4double energy = fastTrack.GetPrimaryTrack()->GetKineticEnergy()/GeV;
98  G4double globalTime = fastTrack.GetPrimaryTrack()->GetStep()->GetPostStepPoint()->GetGlobalTime();
99  G4double charge = fastTrack.GetPrimaryTrack()->GetStep()->GetPreStepPoint()->GetCharge();
100  G4ThreeVector position = fastTrack.GetPrimaryTrack()->GetPosition() / cm;
101  G4ThreeVector momentum = fastTrack.GetPrimaryTrack()->GetMomentum()/GeV;
102  G4int showerType = Gflash::findShowerType(position);
103 
104  // Do actual parameterization. The result of parameterization is gflashHitList
105  theProfile->initialize(showerType,energy,globalTime,charge,position,momentum);
107 
108  //make hits
109  makeHits(fastTrack);
110 }
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 168 of file GflashEMShowerModel.cc.

References eta().

Referenced by ModelTrigger().

168  {
169 
170  G4bool isExcluded=false;
171 
172  //exclude regions where geometry are complicated
173  //+- one supermodule around the EB/EE boundary: 1.479 +- 0.0174*5
174  G4double eta = fastTrack.GetPrimaryTrack()->GetPosition().pseudoRapidity() ;
175  if(std::fabs(eta) > 1.392 && std::fabs(eta) < 1.566) { return true; }
176 
177  return isExcluded;
178 }
T eta() const
G4bool GflashEMShowerModel::IsApplicable ( const G4ParticleDefinition &  particleType)

Definition at line 48 of file GflashEMShowerModel.cc.

48  {
49 
50  return ( &particleType == G4Electron::ElectronDefinition() ||
51  &particleType == G4Positron::PositronDefinition() );
52 
53 }
tuple particleType
Definition: autophobj.py:28
void GflashEMShowerModel::makeHits ( const G4FastTrack &  fastTrack)
private

Definition at line 112 of file GflashEMShowerModel.cc.

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

Referenced by DoIt().

112  {
113 
114  std::vector<GflashHit>& gflashHitList = theProfile->getGflashHitList();
115 
116  theGflashStep->SetTrack(const_cast<G4Track*>(fastTrack.GetPrimaryTrack()));
117 
118  theGflashStep->GetPostStepPoint()->SetProcessDefinedStep(const_cast<G4VProcess*>
119  (fastTrack.GetPrimaryTrack()->GetStep()->GetPostStepPoint()->GetProcessDefinedStep()));
120  theGflashNavigator->SetWorldVolume(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 
127  //put touchable for each hit so that touchable history keeps track of each step.
128  theGflashNavigator->LocateGlobalPointAndUpdateTouchableHandle(spotIter->getPosition(),G4ThreeVector(0,0,0),
129  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 == 0 ) continue;
143 
144  G4LogicalVolume* lv = aCurrentVolume->GetLogicalVolume();
145  if(lv->GetRegion() != theRegion) continue;
146 
147  theGflashStep->GetPreStepPoint()->SetSensitiveDetector(aCurrentVolume->GetLogicalVolume()->GetSensitiveDetector());
148  G4VSensitiveDetector* aSensitive = theGflashStep->GetPreStepPoint()->GetSensitiveDetector();
149 
150  if( aSensitive == 0 ) continue;
151 
152  theGflashStep->SetTotalEnergyDeposit(spotIter->getEnergy());
153  aSensitive->Hit(theGflashStep);
154  }
155 
156 }
std::vector< GflashHit > & getGflashHitList()
SimActivityRegistry::G4StepSignal m_g4StepSignal
G4TouchableHandle theGflashTouchableHandle
const G4Region * theRegion
GflashEMShowerProfile * theProfile
void updateGflashStep(const G4ThreeVector &position, G4double time)
G4Navigator * theGflashNavigator
G4bool GflashEMShowerModel::ModelTrigger ( const G4FastTrack &  fastTrack)

Definition at line 56 of file GflashEMShowerModel.cc.

References excludeDetectorRegion(), GeV, and theRegion.

56  {
57 
58  // Mininum energy cutoff to parameterize
59  if(fastTrack.GetPrimaryTrack()->GetKineticEnergy() < GeV) { return false; }
60  if(excludeDetectorRegion(fastTrack)) { return false; }
61 
62  // This will be changed accordingly when the way
63  // dealing with CaloRegion changes later.
64  G4TouchableHistory* touch =
65  (G4TouchableHistory*)(fastTrack.GetPrimaryTrack()->GetTouchable());
66  G4VPhysicalVolume* pCurrentVolume = touch->GetVolume();
67  if( pCurrentVolume == 0) { return false; }
68 
69  G4LogicalVolume* lv = pCurrentVolume->GetLogicalVolume();
70  if(lv->GetRegion() != theRegion) { return false; }
71  //std::cout << "GflashEMShowerModel::ModelTrigger: LV "
72  // << lv->GetRegion()->GetName() << std::endl;
73 
74  // The parameterization starts inside crystals
75  //std::size_t pos1 = lv->GetName().find("EBRY");
76  //std::size_t pos2 = lv->GetName().find("EFRY");
77 
78  //std::size_t pos3 = lv->GetName().find("HVQ");
79  //std::size_t pos4 = lv->GetName().find("HF");
80  //if(pos1 == std::string::npos && pos2 == std::string::npos &&
81  // pos3 == std::string::npos && pos4 == std::string::npos) return false;
82  //@@@for now, HF is not a part of Gflash Envelopes
83  //if(pos1 == std::string::npos && pos2 == std::string::npos ) return false;
84 
85  return true;
86 
87 }
const double GeV
Definition: MathUtil.h:16
G4bool excludeDetectorRegion(const G4FastTrack &fastTrack)
const G4Region * theRegion
void GflashEMShowerModel::updateGflashStep ( const G4ThreeVector &  position,
G4double  time 
)
private

Definition at line 158 of file GflashEMShowerModel.cc.

References theGflashStep, and theGflashTouchableHandle.

Referenced by makeHits().

160 {
161  theGflashStep->GetPostStepPoint()->SetGlobalTime(timeGlobal);
162  theGflashStep->GetPreStepPoint()->SetPosition(spotPosition);
163  theGflashStep->GetPostStepPoint()->SetPosition(spotPosition);
164  theGflashStep->GetPreStepPoint()->SetTouchableHandle(theGflashTouchableHandle);
165 }
G4TouchableHandle theGflashTouchableHandle

Member Data Documentation

G4Navigator* GflashEMShowerModel::theGflashNavigator
private

Definition at line 81 of file GflashEMShowerModel.h.

Referenced by GflashEMShowerModel(), and makeHits().

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

Definition at line 82 of file GflashEMShowerModel.h.

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

edm::ParameterSet GflashEMShowerModel::theParSet
private

Definition at line 73 of file GflashEMShowerModel.h.

GflashEMShowerProfile* GflashEMShowerModel::theProfile
private

Definition at line 76 of file GflashEMShowerModel.h.

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

const G4Region* GflashEMShowerModel::theRegion
private

Definition at line 78 of file GflashEMShowerModel.h.

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

bool GflashEMShowerModel::theWatcherOn
private

Definition at line 74 of file GflashEMShowerModel.h.

Referenced by makeHits().