CMS 3D CMS Logo

HCalSD.h
Go to the documentation of this file.
1 #ifndef SimG4CMS_HCalSD_h
2 #define SimG4CMS_HCalSD_h
3 // File: HCalSD.h
5 // Description: Stores hits of Hadron calorimeter in appropriate container
6 // Use in your sensitive detector builder:
7 // HCalSD* hcalSD = new HCalSD(SDname, new CaloNumberingScheme());
9 
26 
27 #include "G4String.hh"
28 #include <map>
29 #include <string>
30 
31 class DDFilteredView;
32 class G4LogicalVolume;
33 class G4Material;
34 class G4Step;
35 class HcalTestNS;
36 class TH1F;
37 
38 class HCalSD : public CaloSD, public Observer<const BeginOfJob*> {
39 public:
40  HCalSD(const std::string&,
41  const HcalDDDSimConstants*,
42  const HcalDDDRecConstants*,
44  const HBHEDarkening*,
45  const HBHEDarkening*,
47  edm::ParameterSet const&,
48  const SimTrackManager*);
49  ~HCalSD() override = default;
50  uint32_t setDetUnitId(const G4Step* step) override;
52 
53 protected:
54  double getEnergyDeposit(const G4Step*) override;
55  bool getFromLibrary(const G4Step*) override;
56  using CaloSD::update;
57  void update(const BeginOfJob*) override;
58  void initRun() override;
59  bool filterHit(CaloG4Hit*, double) override;
60  void initEvent(const BeginOfEvent*) override;
61  void endEvent() override;
62 
63 private:
64  void fillLogVolumeVector(const std::string&, const std::vector<std::string>&, std::vector<const G4LogicalVolume*>&);
65  uint32_t setDetUnitId(int, const G4ThreeVector&, int, int);
67  bool isItHF(const G4Step*);
68  bool isItHF(const G4String&);
69  bool isItFibre(const G4LogicalVolume*);
70  bool isItFibre(const G4String&);
71  bool isItPMT(const G4LogicalVolume*);
72  bool isItStraightBundle(const G4LogicalVolume*);
73  bool isItConicalBundle(const G4LogicalVolume*);
74  bool isItScintillator(const G4Material*);
75  bool isItinFidVolume(const G4ThreeVector&);
76  void getFromHFLibrary(const G4Step* step, bool& isKilled);
77  void hitForFibre(const G4Step* step);
78  void getFromParam(const G4Step* step, bool& isKilled);
79  void getHitPMT(const G4Step* step);
80  void getHitFibreBundle(const G4Step* step, bool type);
81  void readWeightFromFile(const std::string&);
82  double layerWeight(int, const G4ThreeVector&, int, int);
83  void plotProfile(const G4Step* step, const G4ThreeVector& pos, double edep, double time, int id);
84  void plotHF(const G4ThreeVector& pos, bool emType);
86  void printVolume(const G4VTouchable* touch) const;
87 
88  std::unique_ptr<HcalNumberingFromDDD> numberingFromDDD;
89  std::unique_ptr<HcalNumberingScheme> numberingScheme;
90  std::unique_ptr<HFShowerLibrary> showerLibrary;
91  std::unique_ptr<HFShower> hfshower;
92  std::unique_ptr<HFShowerParam> showerParam;
93  std::unique_ptr<HFShowerPMT> showerPMT;
94  std::unique_ptr<HFShowerFibreBundle> showerBundle;
95 
100  std::unique_ptr<HFDarkening> m_HFDarkening;
101  std::unique_ptr<HcalTestNS> m_HcalTestNS;
102 
103  static constexpr double maxZ_ = 10000.0;
104  static constexpr double minRoff_ = -1500.0;
105  static constexpr double maxRoff_ = 450.0;
106  static constexpr double slopeHE_ = 0.4;
107  bool isHF;
115  double weight_;
116  int depth_;
117  std::vector<double> gpar;
118  std::vector<int> hfLevels;
119  std::vector<std::string> hfNames;
120  std::vector<std::string> fibreNames;
121  std::vector<std::string> matNames;
122  std::vector<const G4Material*> materials;
123  std::vector<const G4LogicalVolume*> hfLV, fibreLV, pmtLV, fibre1LV, fibre2LV;
124  std::map<uint32_t, double> layerWeights;
125  TH1F *hit_[9], *time_[9], *dist_[9], *hzvem, *hzvhad;
126  std::vector<int> detNull_;
127 };
128 
129 #endif // HCalSD_h
std::vector< std::string > matNames
Definition: HCalSD.h:121
void readWeightFromFile(const std::string &)
Definition: HCalSD.cc:978
std::vector< const G4LogicalVolume * > hfLV
Definition: HCalSD.h:123
bool useParam
Definition: HCalSD.h:112
double eminHitHE
Definition: HCalSD.h:113
TH1F * time_[9]
Definition: HCalSD.h:125
void hitForFibre(const G4Step *step)
Definition: HCalSD.cc:785
void initEvent(const BeginOfEvent *) override
Definition: HCalSD.cc:1109
bool useLayerWt
Definition: HCalSD.h:109
Definition: CaloSD.h:40
std::vector< std::string > fibreNames
Definition: HCalSD.h:120
static constexpr double slopeHE_
Definition: HCalSD.h:106
std::unique_ptr< HFShowerParam > showerParam
Definition: HCalSD.h:92
double weight_
Definition: HCalSD.h:115
std::vector< double > gpar
Definition: HCalSD.h:117
bool useFibreBundle
Definition: HCalSD.h:109
const HcalDDDSimConstants * hcalConstants_
Definition: HCalSD.h:96
double betaThr
Definition: HCalSD.h:111
double deliveredLumi
Definition: HCalSD.h:114
void endEvent() override
Definition: HCalSD.cc:1115
double eminHitHB
Definition: HCalSD.h:113
bool useShowerLibrary
Definition: HCalSD.h:112
std::vector< std::string > hfNames
Definition: HCalSD.h:119
void plotProfile(const G4Step *step, const G4ThreeVector &pos, double edep, double time, int id)
Definition: HCalSD.cc:1022
double birk2
Definition: HCalSD.h:111
bool usePMTHit
Definition: HCalSD.h:109
std::vector< const G4LogicalVolume * > fibre2LV
Definition: HCalSD.h:123
void modifyDepth(HcalNumberingFromDDD::HcalID &id)
Definition: HCalSD.cc:1095
void update(const BeginOfJob *) override
This routine will be called when the appropriate signal arrives.
Definition: HCalSD.cc:576
void setNumberingScheme(HcalNumberingScheme *)
Definition: HCalSD.cc:569
bool agingFlagHE
Definition: HCalSD.h:108
double birk1
Definition: HCalSD.h:111
std::unique_ptr< HFShowerPMT > showerPMT
Definition: HCalSD.h:93
TH1F * hzvhad
Definition: HCalSD.h:125
bool isItConicalBundle(const G4LogicalVolume *)
Definition: HCalSD.cc:710
HCalSD(const std::string &, const HcalDDDSimConstants *, const HcalDDDRecConstants *, const HcalSimulationConstants *, const HBHEDarkening *, const HBHEDarkening *, const SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
Definition: HCalSD.cc:43
double eminHitHF
Definition: HCalSD.h:113
Definition: HCalSD.h:38
static constexpr double minRoff_
Definition: HCalSD.h:104
std::unique_ptr< HcalNumberingFromDDD > numberingFromDDD
Definition: HCalSD.h:88
TH1F * hit_[9]
Definition: HCalSD.h:125
bool agingFlagHB
Definition: HCalSD.h:108
bool testNumber
Definition: HCalSD.h:110
~HCalSD() override=default
bool useBirk
Definition: HCalSD.h:109
std::unique_ptr< HFDarkening > m_HFDarkening
Definition: HCalSD.h:100
double birk3
Definition: HCalSD.h:111
const HBHEDarkening * m_HEDarkening
Definition: HCalSD.h:99
bool applyFidCut
Definition: HCalSD.h:112
TH1F * dist_[9]
Definition: HCalSD.h:125
std::unique_ptr< HcalTestNS > m_HcalTestNS
Definition: HCalSD.h:101
void plotHF(const G4ThreeVector &pos, bool emType)
Definition: HCalSD.cc:1084
std::vector< int > hfLevels
Definition: HCalSD.h:118
void getFromParam(const G4Step *step, bool &isKilled)
Definition: HCalSD.cc:827
void getFromHFLibrary(const G4Step *step, bool &isKilled)
Definition: HCalSD.cc:742
uint32_t setDetUnitId(const G4Step *step) override
Definition: HCalSD.cc:536
bool isItStraightBundle(const G4LogicalVolume *)
Definition: HCalSD.cc:702
bool testNS_
Definition: HCalSD.h:110
bool getFromLibrary(const G4Step *) override
Definition: HCalSD.cc:335
double eminHitHO
Definition: HCalSD.h:113
TH1F * hzvem
Definition: HCalSD.h:125
std::vector< const G4LogicalVolume * > pmtLV
Definition: HCalSD.h:123
bool filterHit(CaloG4Hit *, double) override
Definition: HCalSD.cc:580
bool neutralDensity
Definition: HCalSD.h:110
bool isHF
Definition: HCalSD.h:107
void fillLogVolumeVector(const std::string &, const std::vector< std::string > &, std::vector< const G4LogicalVolume *> &)
Definition: HCalSD.cc:310
static constexpr double maxZ_
Definition: HCalSD.h:103
void getHitFibreBundle(const G4Step *step, bool type)
Definition: HCalSD.cc:916
double getEnergyDeposit(const G4Step *) override
Definition: HCalSD.cc:402
std::unique_ptr< HFShowerFibreBundle > showerBundle
Definition: HCalSD.h:94
bool isItPMT(const G4LogicalVolume *)
Definition: HCalSD.cc:694
bool isItinFidVolume(const G4ThreeVector &)
Definition: HCalSD.cc:726
void printVolume(const G4VTouchable *touch) const
Definition: HCalSD.cc:1124
void getHitPMT(const G4Step *step)
Definition: HCalSD.cc:857
std::vector< int > detNull_
Definition: HCalSD.h:126
std::unique_ptr< HFShowerLibrary > showerLibrary
Definition: HCalSD.h:90
double layerWeight(int, const G4ThreeVector &, int, int)
Definition: HCalSD.cc:1002
const HcalSimulationConstants * hcalSimConstants_
Definition: HCalSD.h:97
int depth_
Definition: HCalSD.h:116
std::map< uint32_t, double > layerWeights
Definition: HCalSD.h:124
bool isItHF(const G4Step *)
Definition: HCalSD.cc:657
bool isItScintillator(const G4Material *)
Definition: HCalSD.cc:718
std::unique_ptr< HcalNumberingScheme > numberingScheme
Definition: HCalSD.h:89
std::unique_ptr< HFShower > hfshower
Definition: HCalSD.h:91
std::vector< const G4LogicalVolume * > fibreLV
Definition: HCalSD.h:123
step
Definition: StallMonitor.cc:83
std::vector< const G4LogicalVolume * > fibre1LV
Definition: HCalSD.h:123
tmp
align.sh
Definition: createJobs.py:716
std::vector< const G4Material * > materials
Definition: HCalSD.h:122
void initRun() override
Definition: HCalSD.cc:578
const HBHEDarkening * m_HBDarkening
Definition: HCalSD.h:98
static constexpr double maxRoff_
Definition: HCalSD.h:105
void update(const BeginOfRun *) override
This routine will be called when the appropriate signal arrives.
Definition: CaloSD.cc:748
bool useHF
Definition: HCalSD.h:112
bool isItFibre(const G4LogicalVolume *)
Definition: HCalSD.cc:678