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 29 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 edm::ParameterSet::getParameter(), theGflashNavigator, theGflashStep, theGflashTouchableHandle, theProfile, theRegion, and theWatcherOn.

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 }
T getParameter(std::string const &) const
G4TouchableHandle theGflashTouchableHandle
GflashEMShowerProfile * theProfile
const G4Region * theRegion
edm::ParameterSet theParSet
G4Navigator * theGflashNavigator
GFlashEMShowerModel::~GFlashEMShowerModel ( )
virtual

Definition at line 42 of file GFlashEMShowerModel.cc.

References theGflashStep, and theProfile.

43 {
44  delete theProfile;
45  delete theGflashStep;
46 }
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.

90 {
91  // Kill the parameterised particle:
92  fastStep.KillPrimaryTrack();
93  fastStep.ProposePrimaryTrackPathLength(0.0);
94 
95  // Input variables for GFlashEMShowerProfile with showerType = 1,5
96  // Shower starts inside crystals
97  G4double energy = fastTrack.GetPrimaryTrack()->GetKineticEnergy()/GeV;
98  G4double globalTime =
99  fastTrack.GetPrimaryTrack()->GetStep()->GetPostStepPoint()->GetGlobalTime();
100  G4double charge =
101  fastTrack.GetPrimaryTrack()->GetStep()->GetPreStepPoint()->GetCharge();
102  G4ThreeVector position = fastTrack.GetPrimaryTrack()->GetPosition()/cm;
103  G4ThreeVector momentum = fastTrack.GetPrimaryTrack()->GetMomentum()/GeV;
104  G4int showerType = Gflash::findShowerType(position);
105 
106  // Do actual parameterization
107  // The result of parameterization is gflashHitList
108  theProfile->initialize(showerType,energy,globalTime,charge,
109  position,momentum);
111 
112  // Make hits
113  makeHits(fastTrack);
114 }
void makeHits(const G4FastTrack &fastTrack)
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)
G4bool GFlashEMShowerModel::excludeDetectorRegion ( const G4FastTrack &  fastTrack)
private

Definition at line 175 of file GFlashEMShowerModel.cc.

References eta().

Referenced by ModelTrigger().

176 {
177  G4bool isExcluded=false;
178 
179  //exclude regions where geometry are complicated
180  //+- one supermodule around the EB/EE boundary: 1.479 +- 0.0174*5
181  G4double eta = fastTrack.GetPrimaryTrack()->GetPosition().pseudoRapidity();
182  if(std::fabs(eta) > 1.392 && std::fabs(eta) < 1.566) { return true; }
183 
184  return isExcluded;
185 }
T eta() const
G4bool GFlashEMShowerModel::IsApplicable ( const G4ParticleDefinition &  particleType)

Definition at line 49 of file GFlashEMShowerModel.cc.

50 {
51  return ( &particleType == G4Electron::Electron() ||
52  &particleType == G4Positron::Positron() );
53 }
void GFlashEMShowerModel::makeHits ( const G4FastTrack &  fastTrack)
private

Definition at line 117 of file GFlashEMShowerModel.cc.

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

Referenced by DoIt().

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

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 
83  return true;
84 
85 }
const G4Region * theRegion
G4bool excludeDetectorRegion(const G4FastTrack &fastTrack)
void GFlashEMShowerModel::updateGflashStep ( const G4ThreeVector &  position,
G4double  time 
)
private

Definition at line 165 of file GFlashEMShowerModel.cc.

References theGflashStep, and theGflashTouchableHandle.

Referenced by makeHits().

167 {
168  theGflashStep->GetPostStepPoint()->SetGlobalTime(timeGlobal);
169  theGflashStep->GetPreStepPoint()->SetPosition(spotPosition);
170  theGflashStep->GetPostStepPoint()->SetPosition(spotPosition);
171  theGflashStep->GetPreStepPoint()->SetTouchableHandle(theGflashTouchableHandle);
172 }
G4TouchableHandle theGflashTouchableHandle

Member Data Documentation

G4Navigator* GFlashEMShowerModel::theGflashNavigator
private

Definition at line 56 of file GFlashEMShowerModel.h.

Referenced by GFlashEMShowerModel(), and makeHits().

G4Step* GFlashEMShowerModel::theGflashStep
private
G4TouchableHandle GFlashEMShowerModel::theGflashTouchableHandle
private

Definition at line 57 of file GFlashEMShowerModel.h.

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

edm::ParameterSet GFlashEMShowerModel::theParSet
private

Definition at line 48 of file GFlashEMShowerModel.h.

GflashEMShowerProfile* GFlashEMShowerModel::theProfile
private

Definition at line 51 of file GFlashEMShowerModel.h.

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

const G4Region* GFlashEMShowerModel::theRegion
private

Definition at line 53 of file GFlashEMShowerModel.h.

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

bool GFlashEMShowerModel::theWatcherOn
private

Definition at line 49 of file GFlashEMShowerModel.h.

Referenced by GFlashEMShowerModel(), and makeHits().