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 (G4String name, G4Envelope *env, edm::ParameterSet parSet)
 
G4bool IsApplicable (const G4ParticleDefinition &)
 
G4bool ModelTrigger (const G4FastTrack &)
 
 ~GflashEMShowerModel ()
 

Private Member Functions

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

Private Attributes

G4Navigator * theGflashNavigator
 
G4Step * theGflashStep
 
G4TouchableHandle theGflashTouchableHandle
 
edm::ParameterSet theParSet
 
GflashEMShowerProfiletheProfile
 
bool theWatcherOn
 

Detailed Description

Definition at line 55 of file GflashEMShowerModel.h.

Constructor & Destructor Documentation

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

Definition at line 22 of file GflashEMShowerModel.cc.

References edm::ParameterSet::getParameter(), theGflashNavigator, theGflashStep, theGflashTouchableHandle, theProfile, and theWatcherOn.

23  : G4VFastSimulationModel(modelName, envelope), theParSet(parSet) {
24 
25  theWatcherOn = parSet.getParameter<bool>("watcherOn");
26 
27  theProfile = new GflashEMShowerProfile(parSet);
28 
29  theGflashStep = new G4Step();
30  theGflashTouchableHandle = new G4TouchableHistory();
31  theGflashNavigator = new G4Navigator();
32 
33 }
T getParameter(std::string const &) const
G4TouchableHandle theGflashTouchableHandle
edm::ParameterSet theParSet
GflashEMShowerProfile * theProfile
G4Navigator * theGflashNavigator
GflashEMShowerModel::~GflashEMShowerModel ( )

Definition at line 37 of file GflashEMShowerModel.cc.

References theGflashStep, and theProfile.

37  {
38 
39  if(theProfile) delete theProfile;
40  if(theGflashStep) delete theGflashStep;
41 }
GflashEMShowerProfile * theProfile

Member Function Documentation

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

Definition at line 88 of file GflashEMShowerModel.cc.

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

88  {
89 
90  // Kill the parameterised particle:
91  fastStep.KillPrimaryTrack();
92  fastStep.ProposePrimaryTrackPathLength(0.0);
93 
94  //input variables for GflashEMShowerProfile with showerType = 1,5 (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. The result of parameterization is gflashHitList
103  theProfile->initialize(showerType,energy,globalTime,charge,position,momentum);
105 
106  //make hits
107  makeHits(fastTrack);
108 }
double charge(const std::vector< uint8_t > &Ampls)
static int position[TOTALCHAMBERS][3]
Definition: ReadPGInfo.cc:509
int findShowerType(const Gflash3Vector position)
GflashEMShowerProfile * theProfile
void initialize(int showerType, double energy, double globalTime, double charge, Gflash3Vector &position, Gflash3Vector &momentum)
void makeHits(const G4FastTrack &fastTrack)
G4bool GflashEMShowerModel::excludeDetectorRegion ( const G4FastTrack &  fastTrack)
private

Definition at line 165 of file GflashEMShowerModel.cc.

References eta().

Referenced by ModelTrigger().

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

Definition at line 43 of file GflashEMShowerModel.cc.

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

Definition at line 110 of file GflashEMShowerModel.cc.

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

Referenced by DoIt().

110  {
111 
112  std::vector<GflashHit>& gflashHitList = theProfile->getGflashHitList();
113 
114  theGflashStep->SetTrack(const_cast<G4Track*>(fastTrack.GetPrimaryTrack()));
115 
116  theGflashStep->GetPostStepPoint()->SetProcessDefinedStep(const_cast<G4VProcess*>
117  (fastTrack.GetPrimaryTrack()->GetStep()->GetPostStepPoint()->GetProcessDefinedStep()));
118  theGflashNavigator->SetWorldVolume(G4TransportationManager::GetTransportationManager()->GetNavigatorForTracking()->GetWorldVolume());
119 
120  std::vector<GflashHit>::const_iterator spotIter = gflashHitList.begin();
121  std::vector<GflashHit>::const_iterator spotIterEnd = gflashHitList.end();
122 
123  for( ; spotIter != spotIterEnd; spotIter++){
124 
125  //put touchable for each hit so that touchable history keeps track of each step.
126  theGflashNavigator->LocateGlobalPointAndUpdateTouchableHandle(spotIter->getPosition(),G4ThreeVector(0,0,0),
127  theGflashTouchableHandle, false);
128  updateGflashStep(spotIter->getPosition(),spotIter->getTime());
129 
130  // if there is a watcher defined in a job and the flag is turned on
131  if(theWatcherOn) {
132  SteppingAction* userSteppingAction = (SteppingAction*) G4EventManager::GetEventManager()->GetUserSteppingAction();
133  userSteppingAction->m_g4StepSignal(theGflashStep);
134  }
135 
136  // Send G4Step information to Hit/Digi if the volume is sensitive
137  // Copied from G4SteppingManager.cc
138 
139  G4VPhysicalVolume* aCurrentVolume = theGflashStep->GetPreStepPoint()->GetPhysicalVolume();
140  if( aCurrentVolume == 0 ) continue;
141 
142  G4LogicalVolume* lv = aCurrentVolume->GetLogicalVolume();
143  if(lv->GetRegion()->GetName() != "CaloRegion") continue;
144 
145  theGflashStep->GetPreStepPoint()->SetSensitiveDetector(aCurrentVolume->GetLogicalVolume()->GetSensitiveDetector());
146  G4VSensitiveDetector* aSensitive = theGflashStep->GetPreStepPoint()->GetSensitiveDetector();
147 
148  if( aSensitive == 0 ) continue;
149 
150  theGflashStep->SetTotalEnergyDeposit(spotIter->getEnergy());
151  aSensitive->Hit(theGflashStep);
152  }
153 
154 }
std::vector< GflashHit > & getGflashHitList()
SimActivityRegistry::G4StepSignal m_g4StepSignal
void updateGflashStep(G4ThreeVector position, G4double time)
G4TouchableHandle theGflashTouchableHandle
GflashEMShowerProfile * theProfile
G4Navigator * theGflashNavigator
G4bool GflashEMShowerModel::ModelTrigger ( const G4FastTrack &  fastTrack)

Definition at line 51 of file GflashEMShowerModel.cc.

References excludeDetectorRegion().

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

Definition at line 156 of file GflashEMShowerModel.cc.

References theGflashStep, and theGflashTouchableHandle.

Referenced by makeHits().

157 {
158  theGflashStep->GetPostStepPoint()->SetGlobalTime(timeGlobal);
159  theGflashStep->GetPreStepPoint()->SetPosition(spotPosition);
160  theGflashStep->GetPostStepPoint()->SetPosition(spotPosition);
161  theGflashStep->GetPreStepPoint()->SetTouchableHandle(theGflashTouchableHandle);
162 }
G4TouchableHandle theGflashTouchableHandle

Member Data Documentation

G4Navigator* GflashEMShowerModel::theGflashNavigator
private

Definition at line 78 of file GflashEMShowerModel.h.

Referenced by GflashEMShowerModel(), and makeHits().

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

Definition at line 79 of file GflashEMShowerModel.h.

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

edm::ParameterSet GflashEMShowerModel::theParSet
private

Definition at line 72 of file GflashEMShowerModel.h.

GflashEMShowerProfile* GflashEMShowerModel::theProfile
private

Definition at line 75 of file GflashEMShowerModel.h.

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

bool GflashEMShowerModel::theWatcherOn
private

Definition at line 73 of file GflashEMShowerModel.h.

Referenced by GflashEMShowerModel(), and makeHits().