CMS 3D CMS Logo

HGCalSD.h
Go to the documentation of this file.
1 #ifndef SimG4CMS_HGCalSD_h
2 #define SimG4CMS_HGCalSD_h
3 // File: HGCalSD.h
5 // Description: Stores hits of the High Granularity Calorimeter (HGC) in the
6 // appropriate container (post TDR version)
8 
13 
14 #include <string>
15 
16 class DDCompactView;
17 class HGCalDDDConstants;
18 class G4LogicalVolume;
19 class G4Step;
20 
21 class HGCalSD : public CaloSD, public Observer<const BeginOfJob *> {
22 public:
23  HGCalSD(const std::string &,
24  const DDCompactView &,
26  edm::ParameterSet const &,
27  const SimTrackManager *);
28  ~HGCalSD() override = default;
29 
30  uint32_t setDetUnitId(const G4Step *step) override;
31 
32 protected:
33  double getEnergyDeposit(const G4Step *) override;
34  using CaloSD::update;
35  void update(const BeginOfJob *) override;
36  void initRun() override;
37  bool filterHit(CaloG4Hit *, double) override;
38 
39 private:
40  uint32_t setDetUnitId(int, int, int, int, G4ThreeVector &);
41  bool isItinFidVolume(const G4ThreeVector &);
42 
44  std::unique_ptr<HGCalNumberingScheme> numberingScheme_;
45  std::unique_ptr<HGCMouseBite> mouseBite_;
54  const double tan30deg_;
55  std::vector<double> angles_;
56 };
57 
58 #endif // HGCalSD_h
double distanceFromEdge_
Definition: HGCalSD.h:49
int levelT2_
Definition: HGCalSD.h:51
double slopeMin_
Definition: HGCalSD.h:49
Definition: CaloSD.h:38
std::unique_ptr< HGCalNumberingScheme > numberingScheme_
Definition: HGCalSD.h:44
std::string nameX_
Definition: HGCalSD.h:47
double mouseBiteCut_
Definition: HGCalSD.h:50
Compact representation of the geometrical detector hierarchy.
Definition: DDCompactView.h:80
HGCalSD(const std::string &, const DDCompactView &, const SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
Definition: HGCalSD.cc:30
void initRun() override
Definition: HGCalSD.cc:218
double getEnergyDeposit(const G4Step *) override
Definition: HGCalSD.cc:96
bool rejectMB_
Definition: HGCalSD.h:53
DetId::Detector mydet_
Definition: HGCalSD.h:46
uint32_t setDetUnitId(const G4Step *step) override
Definition: HGCalSD.cc:127
bool storeAllG4Hits_
Definition: HGCalSD.h:52
std::vector< double > angles_
Definition: HGCalSD.h:55
const HGCalDDDConstants * hgcons_
Definition: HGCalSD.h:43
const double tan30deg_
Definition: HGCalSD.h:54
bool waferRot_
Definition: HGCalSD.h:53
bool isItinFidVolume(const G4ThreeVector &)
Definition: HGCalSD.cc:233
~HGCalSD() override=default
int cornerMinMask_
Definition: HGCalSD.h:51
Detector
Definition: DetId.h:26
HGCalGeometryMode::GeometryMode geom_mode_
Definition: HGCalSD.h:48
double weight_
Definition: HGCalSD.h:50
bool fiducialCut_
Definition: HGCalSD.h:53
bool filterHit(CaloG4Hit *, double) override
Definition: HGCalSD.cc:220
step
Definition: StallMonitor.cc:94
int levelT1_
Definition: HGCalSD.h:51
std::unique_ptr< HGCMouseBite > mouseBite_
Definition: HGCalSD.h:45
double eminHit_
Definition: HGCalSD.h:49
void update(const BeginOfRun *) override
This routine will be called when the appropriate signal arrives.
Definition: CaloSD.cc:460
void update(const BeginOfJob *) override
This routine will be called when the appropriate signal arrives.
Definition: HGCalSD.cc:191