CMS 3D CMS Logo

HGCalTB16SD01.h
Go to the documentation of this file.
1 #ifndef SimG4CMS_HGCalTB16SD01_h
2 #define SimG4CMS_HGCalTB16SD01_h
3 // File: HGCalTB16SD01.h
5 // Description: Stores hits of Beam counters for Fermilab TB16 in appropriate
6 // containers
8 
10 
11 #include <string>
12 
13 class DDCompactView;
14 class G4Step;
15 class G4Material;
16 
17 class HGCalTB16SD01 : public CaloSD {
18 
19 public:
20 
21  HGCalTB16SD01(const std::string& , const DDCompactView &,
23  const SimTrackManager*);
24  ~HGCalTB16SD01() override = default;
25  uint32_t setDetUnitId(const G4Step* step) override;
26  static uint32_t packIndex(int det, int lay, int x, int y);
27  static void unpackIndex(const uint32_t & idx, int& det, int& lay,
28  int& x, int& y);
29 
30 protected:
31 
32  double getEnergyDeposit(const G4Step*) override;
33 
34 private:
35  void initialize(const G4StepPoint* point);
36 
38  bool useBirk_;
39  double birk1_, birk2_, birk3_;
41  G4Material* matScin_;
42 };
43 
44 #endif // HGCalTB16SD01_h
double getEnergyDeposit(const G4Step *) override
Definition: CaloSD.h:37
Compact representation of the geometrical detector hierarchy.
Definition: DDCompactView.h:83
static void unpackIndex(const uint32_t &idx, int &det, int &lay, int &x, int &y)
std::string matName_
Definition: HGCalTB16SD01.h:37
uint32_t setDetUnitId(const G4Step *step) override
~HGCalTB16SD01() override=default
HGCalTB16SD01(const std::string &, const DDCompactView &, const SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
G4Material * matScin_
Definition: HGCalTB16SD01.h:41
void initialize(const G4StepPoint *point)
step
static uint32_t packIndex(int det, int lay, int x, int y)
*vegas h *****************************************************used in the default bin number in original ***version of VEGAS is ***a higher bin number might help to derive a more precise ***grade subtle point
Definition: invegas.h:5