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 
24 
25 #include "G4String.hh"
26 #include <map>
27 #include <string>
28 
29 class DDCompactView;
30 class DDFilteredView;
31 class G4LogicalVolume;
32 class G4Material;
33 class G4Step;
34 class HcalTestNS;
35 class TH1F;
36 
37 class HCalSD : public CaloSD, public Observer<const BeginOfJob *> {
38 
39 public:
40 
41  HCalSD(const std::string& , const DDCompactView &, const SensitiveDetectorCatalog &,
42  edm::ParameterSet const &, const SimTrackManager*);
43  ~HCalSD() override = default;
44  uint32_t setDetUnitId(const G4Step* step) override;
46 
47 protected:
48 
49  double getEnergyDeposit(const G4Step*) override;
50  bool getFromLibrary(const G4Step*) override;
51  void update(const BeginOfJob *) override;
52  void initRun() override;
53  bool filterHit(CaloG4Hit*, double) override;
54 
55 private:
56 
57  void fillLogVolumeVector(const std::string&, const std::string&,
58  const DDCompactView&,
59  std::vector<const G4LogicalVolume*>&,
60  std::vector<G4String>&);
61 
62  uint32_t setDetUnitId(int, const G4ThreeVector&, int, int);
64  std::vector<double> getDDDArray(const std::string&,
65  const DDsvalues_type&);
66  std::vector<G4String> getNames(DDFilteredView&);
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,
84  double edep, double time, int id);
85  void plotHF(const G4ThreeVector& pos, bool emType);
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 
99  std::unique_ptr<HFDarkening> m_HFDarkening;
100  std::unique_ptr<HcalTestNS> m_HcalTestNS;
101 
102  bool isHF;
110  double weight_;
111  int depth_;
112  std::vector<double> gpar;
113  std::vector<int> hfLevels;
114  std::vector<G4String> hfNames, fibreNames, matNames;
115  std::vector<const G4Material*> materials;
116  std::vector<const G4LogicalVolume*> hfLV, fibreLV, pmtLV, fibre1LV, fibre2LV;
117  std::map<uint32_t,double> layerWeights;
118  TH1F *hit_[9], *time_[9], *dist_[9], *hzvem, *hzvhad;
119 
120 };
121 
122 #endif // HCalSD_h
type
Definition: HCALResponse.h:21
void readWeightFromFile(const std::string &)
Definition: HCalSD.cc:960
std::vector< const G4LogicalVolume * > hfLV
Definition: HCalSD.h:116
bool useParam
Definition: HCalSD.h:107
double eminHitHE
Definition: HCalSD.h:108
TH1F * time_[9]
Definition: HCalSD.h:118
void hitForFibre(const G4Step *step)
Definition: HCalSD.cc:767
bool useLayerWt
Definition: HCalSD.h:104
Definition: CaloSD.h:37
std::unique_ptr< HFShowerParam > showerParam
Definition: HCalSD.h:92
double weight_
Definition: HCalSD.h:110
std::vector< double > gpar
Definition: HCalSD.h:112
bool useFibreBundle
Definition: HCalSD.h:104
double betaThr
Definition: HCalSD.h:106
double deliveredLumi
Definition: HCalSD.h:109
double eminHitHB
Definition: HCalSD.h:108
bool useShowerLibrary
Definition: HCalSD.h:107
std::vector< double > getDDDArray(const std::string &, const DDsvalues_type &)
Definition: HCalSD.cc:615
void plotProfile(const G4Step *step, const G4ThreeVector &pos, double edep, double time, int id)
Definition: HCalSD.cc:1008
double birk2
Definition: HCalSD.h:106
bool usePMTHit
Definition: HCalSD.h:104
std::vector< const G4LogicalVolume * > fibre2LV
Definition: HCalSD.h:116
void modifyDepth(HcalNumberingFromDDD::HcalID &id)
Definition: HCalSD.cc:1064
void update(const BeginOfJob *) override
This routine will be called when the appropriate signal arrives.
Definition: HCalSD.cc:526
void setNumberingScheme(HcalNumberingScheme *)
Definition: HCalSD.cc:519
Compact representation of the geometrical detector hierarchy.
Definition: DDCompactView.h:83
bool agingFlagHE
Definition: HCalSD.h:103
double birk1
Definition: HCalSD.h:106
std::unique_ptr< HFShowerPMT > showerPMT
Definition: HCalSD.h:93
TH1F * hzvhad
Definition: HCalSD.h:118
bool isItConicalBundle(const G4LogicalVolume *)
Definition: HCalSD.cc:697
double eminHitHF
Definition: HCalSD.h:108
Definition: HCalSD.h:37
std::unique_ptr< HcalNumberingFromDDD > numberingFromDDD
Definition: HCalSD.h:88
TH1F * hit_[9]
Definition: HCalSD.h:118
HCalSD(const std::string &, const DDCompactView &, const SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
Definition: HCalSD.cc:48
bool agingFlagHB
Definition: HCalSD.h:103
bool testNumber
Definition: HCalSD.h:105
const std::vector< std::string > & getNames() const
~HCalSD() override=default
bool useBirk
Definition: HCalSD.h:104
std::unique_ptr< HFDarkening > m_HFDarkening
Definition: HCalSD.h:99
double birk3
Definition: HCalSD.h:106
const HBHEDarkening * m_HEDarkening
Definition: HCalSD.h:98
std::vector< std::pair< unsigned int, DDValue > > DDsvalues_type
std::maps an index to a DDValue. The index corresponds to the index assigned to the name of the std::...
Definition: DDsvalues.h:20
bool applyFidCut
Definition: HCalSD.h:107
TH1F * dist_[9]
Definition: HCalSD.h:118
std::unique_ptr< HcalTestNS > m_HcalTestNS
Definition: HCalSD.h:100
void plotHF(const G4ThreeVector &pos, bool emType)
Definition: HCalSD.cc:1055
std::vector< int > hfLevels
Definition: HCalSD.h:113
void getFromParam(const G4Step *step, bool &isKilled)
Definition: HCalSD.cc:809
void getFromHFLibrary(const G4Step *step, bool &isKilled)
Definition: HCalSD.cc:724
uint32_t setDetUnitId(const G4Step *step) override
Definition: HCalSD.cc:506
bool isItStraightBundle(const G4LogicalVolume *)
Definition: HCalSD.cc:692
bool testNS_
Definition: HCalSD.h:105
bool getFromLibrary(const G4Step *) override
Definition: HCalSD.cc:316
double eminHitHO
Definition: HCalSD.h:108
TH1F * hzvem
Definition: HCalSD.h:118
std::vector< const G4LogicalVolume * > pmtLV
Definition: HCalSD.h:116
bool filterHit(CaloG4Hit *, double) override
Definition: HCalSD.cc:571
bool neutralDensity
Definition: HCalSD.h:105
bool isHF
Definition: HCalSD.h:102
void getHitFibreBundle(const G4Step *step, bool type)
Definition: HCalSD.cc:900
double getEnergyDeposit(const G4Step *) override
Definition: HCalSD.cc:382
std::unique_ptr< HFShowerFibreBundle > showerBundle
Definition: HCalSD.h:94
bool isItPMT(const G4LogicalVolume *)
Definition: HCalSD.cc:687
bool isItinFidVolume(const G4ThreeVector &)
Definition: HCalSD.cc:707
void getHitPMT(const G4Step *step)
Definition: HCalSD.cc:839
std::unique_ptr< HFShowerLibrary > showerLibrary
Definition: HCalSD.h:90
double layerWeight(int, const G4ThreeVector &, int, int)
Definition: HCalSD.cc:986
std::vector< std::vector< double > > tmp
Definition: MVATrainer.cc:100
int depth_
Definition: HCalSD.h:111
bool isItHF(const G4Step *)
Definition: HCalSD.cc:660
std::vector< G4String > fibreNames
Definition: HCalSD.h:114
bool isItScintillator(const G4Material *)
Definition: HCalSD.cc:702
std::unique_ptr< HcalNumberingScheme > numberingScheme
Definition: HCalSD.h:89
std::vector< G4String > hfNames
Definition: HCalSD.h:114
std::unique_ptr< HFShower > hfshower
Definition: HCalSD.h:91
std::vector< const G4LogicalVolume * > fibreLV
Definition: HCalSD.h:116
std::map< uint32_t, double > layerWeights
Definition: HCalSD.h:117
step
std::vector< const G4LogicalVolume * > fibre1LV
Definition: HCalSD.h:116
std::vector< const G4Material * > materials
Definition: HCalSD.h:115
void initRun() override
Definition: HCalSD.cc:563
const HBHEDarkening * m_HBDarkening
Definition: HCalSD.h:97
std::vector< G4String > matNames
Definition: HCalSD.h:114
bool useHF
Definition: HCalSD.h:107
bool isItFibre(const G4LogicalVolume *)
Definition: HCalSD.cc:677
HcalDDDSimConstants * hcalConstants
Definition: HCalSD.h:96
void fillLogVolumeVector(const std::string &, const std::string &, const DDCompactView &, std::vector< const G4LogicalVolume * > &, std::vector< G4String > &)
Definition: HCalSD.cc:285