CMS 3D CMS Logo

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

#include <HGCScintSD.h>

Inheritance diagram for HGCScintSD:
CaloSD Observer< const BeginOfJob * > SensitiveCaloDetector Observer< const BeginOfRun * > Observer< const BeginOfEvent * > Observer< const BeginOfTrack * > Observer< const EndOfTrack * > Observer< const EndOfEvent * > SensitiveDetector

Public Member Functions

 HGCScintSD (const std::string &, const DDCompactView &, const SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
 
uint32_t setDetUnitId (const G4Step *step) override
 
 ~HGCScintSD () override=default
 
- Public Member Functions inherited from CaloSD
 CaloSD (const std::string &aSDname, const DDCompactView &cpv, const SensitiveDetectorCatalog &clg, edm::ParameterSet const &p, const SimTrackManager *, float timeSlice=1., bool ignoreTkID=false)
 
void clear () override
 
void clearHits () override
 
void DrawAll () override
 
void EndOfEvent (G4HCofThisEvent *eventHC) override
 
void fillHits (edm::PCaloHitContainer &, const std::string &) override
 
void Initialize (G4HCofThisEvent *HCE) override
 
void PrintAll () override
 
G4bool ProcessHits (G4Step *step, G4TouchableHistory *) override
 
bool ProcessHits (G4GFlashSpot *aSpot, G4TouchableHistory *) override
 
void reset () override
 
 ~CaloSD () override
 
- Public Member Functions inherited from SensitiveCaloDetector
 SensitiveCaloDetector (const std::string &iname, const DDCompactView &cpv, const SensitiveDetectorCatalog &clg, edm::ParameterSet const &p)
 
- Public Member Functions inherited from SensitiveDetector
void EndOfEvent (G4HCofThisEvent *eventHC) override
 
const std::vector< std::string > & getNames () const
 
void Initialize (G4HCofThisEvent *eventHC) override
 
bool isCaloSD () const
 
 SensitiveDetector (const std::string &iname, const DDCompactView &cpv, const SensitiveDetectorCatalog &, edm::ParameterSet const &p, bool calo)
 
 ~SensitiveDetector () override
 
- Public Member Functions inherited from Observer< const BeginOfRun * >
 Observer ()
 
void slotForUpdate (const BeginOfRun * iT)
 
virtual ~Observer ()
 
- Public Member Functions inherited from Observer< const BeginOfEvent * >
 Observer ()
 
void slotForUpdate (const BeginOfEvent * iT)
 
virtual ~Observer ()
 
- Public Member Functions inherited from Observer< const BeginOfTrack * >
 Observer ()
 
void slotForUpdate (const BeginOfTrack * iT)
 
virtual ~Observer ()
 
- Public Member Functions inherited from Observer< const EndOfTrack * >
 Observer ()
 
void slotForUpdate (const EndOfTrack * iT)
 
virtual ~Observer ()
 
- Public Member Functions inherited from Observer< const EndOfEvent * >
 Observer ()
 
void slotForUpdate (const EndOfEvent * iT)
 
virtual ~Observer ()
 
- Public Member Functions inherited from Observer< const BeginOfJob * >
 Observer ()
 
void slotForUpdate (const BeginOfJob * iT)
 
virtual ~Observer ()
 

Protected Member Functions

bool filterHit (CaloG4Hit *, double) override
 
double getEnergyDeposit (const G4Step *) override
 
void initRun () override
 
void update (const BeginOfJob *) override
 This routine will be called when the appropriate signal arrives. More...
 
- Protected Member Functions inherited from CaloSD
bool checkHit ()
 
CaloG4HitcreateNewHit (const G4Step *)
 
virtual void endEvent ()
 
double getAttenuation (const G4Step *aStep, double birk1, double birk2, double birk3) const
 
virtual uint16_t getDepth (const G4Step *)
 
virtual bool getFromLibrary (const G4Step *step)
 
int getNumberOfHits ()
 
double getResponseWt (const G4Track *)
 
virtual int getTrackID (const G4Track *)
 
bool hitExists (const G4Step *)
 
virtual void initEvent (const BeginOfEvent *)
 
void processHit (const G4Step *step)
 
void resetForNewPrimary (const G4Step *)
 
void setNumberCheckedHits (int val)
 
void setParameterized (bool val)
 
G4ThreeVector setToGlobal (const G4ThreeVector &, const G4VTouchable *) const
 
G4ThreeVector setToLocal (const G4ThreeVector &, const G4VTouchable *) const
 
virtual int setTrackID (const G4Step *)
 
void setUseMap (bool val)
 
void update (const BeginOfRun *) override
 This routine will be called when the appropriate signal arrives. More...
 
void update (const BeginOfEvent *) override
 This routine will be called when the appropriate signal arrives. More...
 
void update (const BeginOfTrack *trk) override
 This routine will be called when the appropriate signal arrives. More...
 
void update (const EndOfTrack *trk) override
 This routine will be called when the appropriate signal arrives. More...
 
void update (const ::EndOfEvent *) override
 
void updateHit (CaloG4Hit *)
 
- Protected Member Functions inherited from SensitiveDetector
TrackInformationcmsTrackInformation (const G4Track *aTrack)
 
Local3DPoint ConvertToLocal3DPoint (const G4ThreeVector &point) const
 
Local3DPoint FinalStepPosition (const G4Step *step, coordinates) const
 
Local3DPoint InitialStepPosition (const G4Step *step, coordinates) const
 
Local3DPoint LocalPostStepPosition (const G4Step *step) const
 
Local3DPoint LocalPreStepPosition (const G4Step *step) const
 
void NaNTrap (const G4Step *step) const
 
void setNames (const std::vector< std::string > &)
 
- Protected Member Functions inherited from Observer< const EndOfEvent * >
virtual void update (const EndOfEvent *)=0
 This routine will be called when the appropriate signal arrives. More...
 

Private Member Functions

bool isItinFidVolume (const G4ThreeVector &)
 
uint32_t setDetUnitId (int, int, int, int, G4ThreeVector &)
 

Private Attributes

double birk1_
 
double birk2_
 
double birk3_
 
double distanceFromEdge_
 
double eminHit_
 
bool fiducialCut_
 
HGCalGeometryMode::GeometryMode geom_mode_
 
const HGCalDDDConstantshgcons_
 
int levelT1_
 
int levelT2_
 
DetId::Detector mydet_
 
std::string nameX_
 
std::unique_ptr< HGCalNumberingSchemenumberingScheme_
 
double slopeMin_
 
bool storeAllG4Hits_
 
bool useBirk_
 
double weight_
 

Additional Inherited Members

- Protected Types inherited from SensitiveDetector
enum  coordinates { WorldCoordinates, LocalCoordinates }
 
- Protected Attributes inherited from CaloSD
CaloG4HitcurrentHit
 
CaloHitID currentID
 
float edepositEM
 
float edepositHAD
 
double eminHit
 
double energyCut
 
G4ThreeVector entranceLocal
 
G4ThreeVector entrancePoint
 
bool forceSave
 
float incidentEnergy
 
double kmaxIon
 
double kmaxNeutron
 
double kmaxProton
 
G4ThreeVector posGlobal
 
CaloHitID previousID
 
bool suppressHeavy
 
double tmaxHit
 

Detailed Description

Definition at line 20 of file HGCScintSD.h.

Constructor & Destructor Documentation

HGCScintSD::HGCScintSD ( const std::string &  name,
const DDCompactView cpv,
const SensitiveDetectorCatalog clg,
edm::ParameterSet const &  p,
const SimTrackManager manager 
)

Definition at line 37 of file HGCScintSD.cc.

References birk1_, birk2_, birk3_, distanceFromEdge_, eminHit_, fiducialCut_, DetId::Forward, g, edm::ParameterSet::getParameter(), DetId::HGCalHSc, MeV, mydet_, dataset::name, nameX_, numberingScheme_, CaloSD::setNumberCheckedHits(), CaloSD::setUseMap(), storeAllG4Hits_, and useBirk_.

42  : CaloSD(name,
43  cpv,
44  clg,
45  p,
46  manager,
47  (float)(p.getParameter<edm::ParameterSet>("HGCSD").getParameter<double>("TimeSliceUnit")),
48  p.getParameter<edm::ParameterSet>("HGCSD").getParameter<bool>("IgnoreTrackID")),
49  hgcons_(nullptr),
50  slopeMin_(0),
51  levelT1_(99),
52  levelT2_(99) {
53  numberingScheme_.reset(nullptr);
54 
55  edm::ParameterSet m_HGC = p.getParameter<edm::ParameterSet>("HGCScintSD");
56  eminHit_ = m_HGC.getParameter<double>("EminHit") * CLHEP::MeV;
57  fiducialCut_ = m_HGC.getParameter<bool>("FiducialCut");
58  distanceFromEdge_ = m_HGC.getParameter<double>("DistanceFromEdge");
59  useBirk_ = m_HGC.getParameter<bool>("UseBirkLaw");
60  birk1_ = m_HGC.getParameter<double>("BirkC1") * (CLHEP::g / (CLHEP::MeV * CLHEP::cm2));
61  birk2_ = m_HGC.getParameter<double>("BirkC2");
62  birk3_ = m_HGC.getParameter<double>("BirkC3");
63  storeAllG4Hits_ = m_HGC.getParameter<bool>("StoreAllG4Hits");
64 
65  if (storeAllG4Hits_) {
66  setUseMap(false);
68  }
69 
70  //this is defined in the hgcsens.xml
71  G4String myName = name;
73  nameX_ = "HGCal";
74  if (myName.find("HitsHEback") != std::string::npos) {
76  nameX_ = "HGCalHEScintillatorSensitive";
77  }
78 
79 #ifdef EDM_ML_DEBUG
80  edm::LogVerbatim("HGCSim") << "**************************************************"
81  << "\n"
82  << "* *"
83  << "\n"
84  << "* Constructing a HGCScintSD with name " << name << "\n"
85  << "* *"
86  << "\n"
87  << "**************************************************";
88 #endif
89  edm::LogVerbatim("HGCSim") << "HGCScintSD:: Threshold for storing hits: " << eminHit_ << " for " << nameX_
90  << " detector " << mydet_;
91  edm::LogVerbatim("HGCSim") << "Flag for storing individual Geant4 Hits " << storeAllG4Hits_;
92  edm::LogVerbatim("HGCSim") << "Fiducial volume cut with cut from eta/phi "
93  << "boundary " << fiducialCut_ << " at " << distanceFromEdge_;
94  edm::LogVerbatim("HGCSim") << "Use of Birks law is set to " << useBirk_
95  << " with three constants kB = " << birk1_ << ", C1 = " << birk2_ << ", C2 = " << birk3_;
96 }
double birk2_
Definition: HGCScintSD.h:51
T getParameter(std::string const &) const
int levelT2_
Definition: HGCScintSD.h:48
void setNumberCheckedHits(int val)
Definition: CaloSD.h:112
double birk1_
Definition: HGCScintSD.h:51
double distanceFromEdge_
Definition: HGCScintSD.h:47
int levelT1_
Definition: HGCScintSD.h:48
void setUseMap(bool val)
Definition: CaloSD.h:101
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e g
Definition: Activities.doc:4
bool fiducialCut_
Definition: HGCScintSD.h:49
bool useBirk_
Definition: HGCScintSD.h:50
const HGCalDDDConstants * hgcons_
Definition: HGCScintSD.h:42
const double MeV
CaloSD(const std::string &aSDname, const DDCompactView &cpv, const SensitiveDetectorCatalog &clg, edm::ParameterSet const &p, const SimTrackManager *, float timeSlice=1., bool ignoreTkID=false)
Definition: CaloSD.cc:27
double eminHit_
Definition: HGCScintSD.h:47
double slopeMin_
Definition: HGCScintSD.h:47
std::unique_ptr< HGCalNumberingScheme > numberingScheme_
Definition: HGCScintSD.h:43
std::string nameX_
Definition: HGCScintSD.h:45
DetId::Detector mydet_
Definition: HGCScintSD.h:44
double birk3_
Definition: HGCScintSD.h:51
bool storeAllG4Hits_
Definition: HGCScintSD.h:49
HGCScintSD::~HGCScintSD ( )
overridedefault

Member Function Documentation

bool HGCScintSD::filterHit ( CaloG4Hit aHit,
double  time 
)
overrideprotectedvirtual

Reimplemented from CaloSD.

Definition at line 209 of file HGCScintSD.cc.

References eminHit_, CaloG4Hit::getEnergyDeposit(), and CaloSD::tmaxHit.

209  {
210  return ((time <= tmaxHit) && (aHit->getEnergyDeposit() > eminHit_));
211 }
double tmaxHit
Definition: CaloSD.h:133
double eminHit_
Definition: HGCScintSD.h:47
double getEnergyDeposit() const
Definition: CaloG4Hit.h:77
double HGCScintSD::getEnergyDeposit ( const G4Step *  aStep)
overrideprotectedvirtual

Reimplemented from CaloSD.

Definition at line 98 of file HGCScintSD.cc.

References funct::abs(), birk1_, birk2_, birk3_, CaloSD::getAttenuation(), CaloSD::getResponseWt(), alignCSCRings::r, slopeMin_, useBirk_, weight_, and z.

98  {
99  double r = aStep->GetPreStepPoint()->GetPosition().perp();
100  double z = std::abs(aStep->GetPreStepPoint()->GetPosition().z());
101 #ifdef EDM_ML_DEBUG
102  G4int parCode = aStep->GetTrack()->GetDefinition()->GetPDGEncoding();
103  G4String parName = aStep->GetTrack()->GetDefinition()->GetParticleName();
104  G4LogicalVolume* lv = aStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume();
105  edm::LogVerbatim("HGCSim") << "HGCScintSD: Hit from standard path from " << lv->GetName() << " for Track "
106  << aStep->GetTrack()->GetTrackID() << " (" << parCode << ":" << parName << ") R = " << r
107  << " Z = " << z << " slope = " << r / z << ":" << slopeMin_;
108 #endif
109  // Apply fiducial cut
110  if (r < z * slopeMin_) {
111 #ifdef EDM_ML_DEBUG
112  edm::LogVerbatim("HGCSim") << "HGCScintSD: Fiducial Volume cut";
113 #endif
114  return 0.0;
115  }
116 
117  double wt1 = getResponseWt(aStep->GetTrack());
118  double wt2 = aStep->GetTrack()->GetWeight();
119  double wt3 = (useBirk_ ? getAttenuation(aStep, birk1_, birk2_, birk3_) : 1.0);
120  double destep = weight_ * wt1 * wt3 * (aStep->GetTotalEnergyDeposit());
121  if (wt2 > 0)
122  destep *= wt2;
123 #ifdef EDM_ML_DEBUG
124  edm::LogVerbatim("HGCalSim") << "HGCScintSD: weights= " << weight_ << ":" << wt1 << ":" << wt2 << ":" << wt3
125  << " Total weight " << weight_ * wt1 * wt2 * wt3
126  << " deStep: " << aStep->GetTotalEnergyDeposit() << ":" << destep;
127 #endif
128  return destep;
129 }
double birk2_
Definition: HGCScintSD.h:51
double birk1_
Definition: HGCScintSD.h:51
bool useBirk_
Definition: HGCScintSD.h:50
double weight_
Definition: HGCScintSD.h:51
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double slopeMin_
Definition: HGCScintSD.h:47
double getAttenuation(const G4Step *aStep, double birk1, double birk2, double birk3) const
Definition: CaloSD.cc:438
double getResponseWt(const G4Track *)
Definition: CaloSD.cc:623
double birk3_
Definition: HGCScintSD.h:51
void HGCScintSD::initRun ( )
overrideprotectedvirtual

Reimplemented from CaloSD.

Definition at line 207 of file HGCScintSD.cc.

207 {}
bool HGCScintSD::isItinFidVolume ( const G4ThreeVector &  pos)
private

Definition at line 218 of file HGCScintSD.cc.

References distanceFromEdge_, HGCalDDDConstants::distFromEdgeTrap(), fiducialCut_, and hgcons_.

Referenced by setDetUnitId().

218  {
219  if (fiducialCut_) {
220  return (hgcons_->distFromEdgeTrap(pos.x(), pos.y(), pos.z()) > distanceFromEdge_);
221  } else {
222  return true;
223  }
224 }
double distanceFromEdge_
Definition: HGCScintSD.h:47
bool fiducialCut_
Definition: HGCScintSD.h:49
const HGCalDDDConstants * hgcons_
Definition: HGCScintSD.h:42
double distFromEdgeTrap(double x, double y, double z) const
uint32_t HGCScintSD::setDetUnitId ( const G4Step *  step)
overridevirtual

Implements CaloSD.

Definition at line 131 of file HGCScintSD.cc.

References TauDecayModes::dec, triggerObjects_cff::id, isItinFidVolume(), levelT1_, and levelT2_.

131  {
132  const G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
133  const G4VTouchable* touch = preStepPoint->GetTouchable();
134 
135  //determine the exact position in global coordinates in the mass geometry
136  G4ThreeVector hitPoint = preStepPoint->GetPosition();
137  float globalZ = touch->GetTranslation(0).z();
138  int iz(globalZ > 0 ? 1 : -1);
139 
140  int layer, module, cell;
141  if ((touch->GetHistoryDepth() == levelT1_) || (touch->GetHistoryDepth() == levelT2_)) {
142  layer = touch->GetReplicaNumber(0);
143  module = -1;
144  cell = -1;
145 #ifdef EDM_ML_DEBUG
146  edm::LogVerbatim("HGCSim") << "DepthsTop: " << touch->GetHistoryDepth() << ":" << levelT1_ << ":" << levelT2_
147  << " name " << touch->GetVolume(0)->GetName() << " layer:module:cell " << layer << ":"
148  << module << ":" << cell;
149 #endif
150  } else {
151  layer = touch->GetReplicaNumber(3);
152  module = touch->GetReplicaNumber(2);
153  cell = touch->GetReplicaNumber(1);
154 #ifdef EDM_ML_DEBUG
155  edm::LogVerbatim("HGCSim") << "DepthsInside: " << touch->GetHistoryDepth() << " name "
156  << touch->GetVolume(0)->GetName() << " layer:module:cell " << layer << ":" << module
157  << ":" << cell;
158 #endif
159  }
160 #ifdef EDM_ML_DEBUG
161  G4Material* mat = aStep->GetPreStepPoint()->GetMaterial();
162  edm::LogVerbatim("HGCSim") << "Depths: " << touch->GetHistoryDepth() << " name " << touch->GetVolume(0)->GetName()
163  << ":" << touch->GetReplicaNumber(0) << " " << touch->GetVolume(1)->GetName() << ":"
164  << touch->GetReplicaNumber(1) << " " << touch->GetVolume(2)->GetName() << ":"
165  << touch->GetReplicaNumber(2) << " " << touch->GetVolume(3)->GetName() << ":"
166  << touch->GetReplicaNumber(3) << " " << touch->GetVolume(4)->GetName() << ":"
167  << touch->GetReplicaNumber(4) << " "
168  << " layer:module:cell " << layer << ":" << module << ":" << cell << " Material "
169  << mat->GetName() << ":" << mat->GetRadlen();
170 #endif
171  // The following statement should be examined later before elimination
172  if (aStep->GetPreStepPoint()->GetMaterial()->GetRadlen() > 100000.)
173  return 0;
174 
175  uint32_t id = setDetUnitId(layer, module, cell, iz, hitPoint);
176  if (!isItinFidVolume(hitPoint)) {
177 #ifdef EDM_ML_DEBUG
178  edm::LogVerbatim("HGCSim") << "ID " << std::hex << id << std::dec << " " << HGCScintillatorDetId(id)
179  << " is rejected by fiducilal volume cut";
180 #endif
181  id = 0;
182  }
183  return id;
184 }
int levelT2_
Definition: HGCScintSD.h:48
uint32_t setDetUnitId(const G4Step *step) override
Definition: HGCScintSD.cc:131
int levelT1_
Definition: HGCScintSD.h:48
bool isItinFidVolume(const G4ThreeVector &)
Definition: HGCScintSD.cc:218
Definition: vlib.h:208
uint32_t HGCScintSD::setDetUnitId ( int  layer,
int  module,
int  cell,
int  iz,
G4ThreeVector &  pos 
)
private

Definition at line 213 of file HGCScintSD.cc.

References triggerObjects_cff::id, numberingScheme_, and weight_.

213  {
214  uint32_t id = numberingScheme_ ? numberingScheme_->getUnitID(layer, module, cell, iz, pos, weight_) : 0;
215  return id;
216 }
double weight_
Definition: HGCScintSD.h:51
std::unique_ptr< HGCalNumberingScheme > numberingScheme_
Definition: HGCScintSD.h:43
Definition: vlib.h:208
void HGCScintSD::update ( const BeginOfJob )
overrideprotectedvirtual

This routine will be called when the appropriate signal arrives.

Implements Observer< const BeginOfJob * >.

Definition at line 186 of file HGCScintSD.cc.

References Exception, geom_mode_, HGCalDDDConstants::geomMode(), edm::EventSetup::get(), hgcons_, levelT1_, levelT2_, HGCalDDDConstants::levelTop(), HGCalDDDConstants::minSlope(), mydet_, nameX_, numberingScheme_, and slopeMin_.

Referenced by progressbar.ProgressBar::__next__(), MatrixUtil.Matrix::__setitem__(), MatrixUtil.Steps::__setitem__(), Vispa.Gui.VispaWidget.VispaWidget::autosize(), Vispa.Views.LineDecayView.LineDecayContainer::createObject(), Vispa.Views.LineDecayView.LineDecayContainer::deselectAllObjects(), Vispa.Gui.VispaWidgetOwner.VispaWidgetOwner::deselectAllWidgets(), Vispa.Gui.VispaWidget.VispaWidget::enableAutosizing(), progressbar.ProgressBar::finish(), Vispa.Gui.MenuWidget.MenuWidget::leaveEvent(), Vispa.Gui.VispaWidgetOwner.VispaWidgetOwner::mouseMoveEvent(), Vispa.Gui.MenuWidget.MenuWidget::mouseMoveEvent(), Vispa.Views.LineDecayView.LineDecayContainer::mouseMoveEvent(), Vispa.Gui.VispaWidgetOwner.VispaWidgetOwner::mouseReleaseEvent(), Vispa.Views.LineDecayView.LineDecayContainer::objectMoved(), MatrixUtil.Steps::overwrite(), Vispa.Views.LineDecayView.LineDecayContainer::removeObject(), Vispa.Gui.ConnectableWidget.ConnectableWidget::removePorts(), Vispa.Gui.FindDialog.FindDialog::reset(), Vispa.Gui.PortConnection.PointToPointConnection::select(), Vispa.Gui.VispaWidget.VispaWidget::select(), Vispa.Views.LineDecayView.LineDecayContainer::select(), Vispa.Gui.VispaWidget.VispaWidget::setText(), Vispa.Gui.VispaWidget.VispaWidget::setTitle(), Vispa.Gui.ZoomableWidget.ZoomableWidget::setZoom(), Vispa.Views.LineDecayView.LineDecayContainer::setZoom(), and Vispa.Gui.PortConnection.PointToPointConnection::updateConnection().

186  {
187  const edm::EventSetup* es = (*job)();
189  es->get<IdealGeometryRecord>().get(nameX_, hdc);
190  if (hdc.isValid()) {
191  hgcons_ = hdc.product();
194  levelT1_ = hgcons_->levelTop(0);
195  levelT2_ = hgcons_->levelTop(1);
196 #ifdef EDM_ML_DEBUG
197  edm::LogVerbatim("HGCSim") << "HGCScintSD::Initialized with mode " << geom_mode_ << " Slope cut " << slopeMin_
198  << " top Level " << levelT1_ << ":" << levelT2_;
199 #endif
200 
202  } else {
203  throw cms::Exception("Unknown", "HGCScintSD") << "Cannot find HGCalDDDConstants for " << nameX_ << "\n";
204  }
205 }
HGCalGeometryMode::GeometryMode geom_mode_
Definition: HGCScintSD.h:46
int levelT2_
Definition: HGCScintSD.h:48
int levelT1_
Definition: HGCScintSD.h:48
const HGCalDDDConstants * hgcons_
Definition: HGCScintSD.h:42
double minSlope() const
HGCalGeometryMode::GeometryMode geomMode() const
double slopeMin_
Definition: HGCScintSD.h:47
std::unique_ptr< HGCalNumberingScheme > numberingScheme_
Definition: HGCScintSD.h:43
std::string nameX_
Definition: HGCScintSD.h:45
T get() const
Definition: EventSetup.h:71
DetId::Detector mydet_
Definition: HGCScintSD.h:44
int levelTop(int ind=0) const

Member Data Documentation

double HGCScintSD::birk1_
private

Definition at line 51 of file HGCScintSD.h.

Referenced by getEnergyDeposit(), and HGCScintSD().

double HGCScintSD::birk2_
private

Definition at line 51 of file HGCScintSD.h.

Referenced by getEnergyDeposit(), and HGCScintSD().

double HGCScintSD::birk3_
private

Definition at line 51 of file HGCScintSD.h.

Referenced by getEnergyDeposit(), and HGCScintSD().

double HGCScintSD::distanceFromEdge_
private

Definition at line 47 of file HGCScintSD.h.

Referenced by HGCScintSD(), and isItinFidVolume().

double HGCScintSD::eminHit_
private

Definition at line 47 of file HGCScintSD.h.

Referenced by filterHit(), and HGCScintSD().

bool HGCScintSD::fiducialCut_
private

Definition at line 49 of file HGCScintSD.h.

Referenced by HGCScintSD(), and isItinFidVolume().

HGCalGeometryMode::GeometryMode HGCScintSD::geom_mode_
private

Definition at line 46 of file HGCScintSD.h.

Referenced by update().

const HGCalDDDConstants* HGCScintSD::hgcons_
private

Definition at line 42 of file HGCScintSD.h.

Referenced by isItinFidVolume(), and update().

int HGCScintSD::levelT1_
private

Definition at line 48 of file HGCScintSD.h.

Referenced by setDetUnitId(), and update().

int HGCScintSD::levelT2_
private

Definition at line 48 of file HGCScintSD.h.

Referenced by setDetUnitId(), and update().

DetId::Detector HGCScintSD::mydet_
private

Definition at line 44 of file HGCScintSD.h.

Referenced by HGCScintSD(), and update().

std::string HGCScintSD::nameX_
private

Definition at line 45 of file HGCScintSD.h.

Referenced by HGCScintSD(), and update().

std::unique_ptr<HGCalNumberingScheme> HGCScintSD::numberingScheme_
private

Definition at line 43 of file HGCScintSD.h.

Referenced by HGCScintSD(), setDetUnitId(), and update().

double HGCScintSD::slopeMin_
private

Definition at line 47 of file HGCScintSD.h.

Referenced by getEnergyDeposit(), and update().

bool HGCScintSD::storeAllG4Hits_
private

Definition at line 49 of file HGCScintSD.h.

Referenced by HGCScintSD().

bool HGCScintSD::useBirk_
private

Definition at line 50 of file HGCScintSD.h.

Referenced by getEnergyDeposit(), and HGCScintSD().

double HGCScintSD::weight_
private

Definition at line 51 of file HGCScintSD.h.

Referenced by getEnergyDeposit(), and setDetUnitId().