CMS 3D CMS Logo

All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 edm::EventSetup&,
43  edm::ParameterSet const&,
44  const SimTrackManager*);
45  ~HCalSD() override = default;
46  uint32_t setDetUnitId(const G4Step* step) override;
48 
49 protected:
50  double getEnergyDeposit(const G4Step*) override;
51  bool getFromLibrary(const G4Step*) override;
52  using CaloSD::update;
53  void update(const BeginOfJob*) override;
54  void initRun() override;
55  bool filterHit(CaloG4Hit*, double) override;
56 
57 private:
58  void fillLogVolumeVector(const std::string&, const std::vector<std::string>&, std::vector<const G4LogicalVolume*>&);
59  uint32_t setDetUnitId(int, const G4ThreeVector&, int, int);
61  bool isItHF(const G4Step*);
62  bool isItHF(const G4String&);
63  bool isItFibre(const G4LogicalVolume*);
64  bool isItFibre(const G4String&);
65  bool isItPMT(const G4LogicalVolume*);
66  bool isItStraightBundle(const G4LogicalVolume*);
67  bool isItConicalBundle(const G4LogicalVolume*);
68  bool isItScintillator(const G4Material*);
69  bool isItinFidVolume(const G4ThreeVector&);
70  void getFromHFLibrary(const G4Step* step, bool& isKilled);
71  void hitForFibre(const G4Step* step);
72  void getFromParam(const G4Step* step, bool& isKilled);
73  void getHitPMT(const G4Step* step);
74  void getHitFibreBundle(const G4Step* step, bool type);
75  void readWeightFromFile(const std::string&);
76  double layerWeight(int, const G4ThreeVector&, int, int);
77  void plotProfile(const G4Step* step, const G4ThreeVector& pos, double edep, double time, int id);
78  void plotHF(const G4ThreeVector& pos, bool emType);
80 
81  std::unique_ptr<HcalNumberingFromDDD> numberingFromDDD;
82  std::unique_ptr<HcalNumberingScheme> numberingScheme;
83  std::unique_ptr<HFShowerLibrary> showerLibrary;
84  std::unique_ptr<HFShower> hfshower;
85  std::unique_ptr<HFShowerParam> showerParam;
86  std::unique_ptr<HFShowerPMT> showerPMT;
87  std::unique_ptr<HFShowerFibreBundle> showerBundle;
88 
93  std::unique_ptr<HFDarkening> m_HFDarkening;
94  std::unique_ptr<HcalTestNS> m_HcalTestNS;
95 
96  bool isHF;
104  double weight_;
105  int depth_;
106  std::vector<double> gpar;
107  std::vector<int> hfLevels;
108  std::vector<std::string> hfNames;
109  std::vector<std::string> fibreNames;
110  std::vector<std::string> matNames;
111  std::vector<const G4Material*> materials;
112  std::vector<const G4LogicalVolume*> hfLV, fibreLV, pmtLV, fibre1LV, fibre2LV;
113  std::map<uint32_t, double> layerWeights;
114  TH1F *hit_[9], *time_[9], *dist_[9], *hzvem, *hzvhad;
115 };
116 
117 #endif // HCalSD_h
type
Definition: HCALResponse.h:21
std::vector< std::string > matNames
Definition: HCalSD.h:110
void readWeightFromFile(const std::string &)
Definition: HCalSD.cc:928
std::vector< const G4LogicalVolume * > hfLV
Definition: HCalSD.h:112
bool useParam
Definition: HCalSD.h:101
double eminHitHE
Definition: HCalSD.h:102
TH1F * time_[9]
Definition: HCalSD.h:114
HCalSD(const std::string &, const edm::EventSetup &, const SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
Definition: HCalSD.cc:44
void hitForFibre(const G4Step *step)
Definition: HCalSD.cc:735
void fillLogVolumeVector(const std::string &, const std::vector< std::string > &, std::vector< const G4LogicalVolume * > &)
Definition: HCalSD.cc:324
bool useLayerWt
Definition: HCalSD.h:98
Definition: CaloSD.h:38
std::vector< std::string > fibreNames
Definition: HCalSD.h:109
std::unique_ptr< HFShowerParam > showerParam
Definition: HCalSD.h:85
double weight_
Definition: HCalSD.h:104
std::vector< double > gpar
Definition: HCalSD.h:106
bool useFibreBundle
Definition: HCalSD.h:98
const HcalDDDSimConstants * hcalConstants_
Definition: HCalSD.h:89
double betaThr
Definition: HCalSD.h:100
double deliveredLumi
Definition: HCalSD.h:103
double eminHitHB
Definition: HCalSD.h:102
bool useShowerLibrary
Definition: HCalSD.h:101
std::vector< std::string > hfNames
Definition: HCalSD.h:108
void plotProfile(const G4Step *step, const G4ThreeVector &pos, double edep, double time, int id)
Definition: HCalSD.cc:972
double birk2
Definition: HCalSD.h:100
bool usePMTHit
Definition: HCalSD.h:98
std::vector< const G4LogicalVolume * > fibre2LV
Definition: HCalSD.h:112
void modifyDepth(HcalNumberingFromDDD::HcalID &id)
Definition: HCalSD.cc:1043
void update(const BeginOfJob *) override
This routine will be called when the appropriate signal arrives.
Definition: HCalSD.cc:557
void setNumberingScheme(HcalNumberingScheme *)
Definition: HCalSD.cc:550
bool agingFlagHE
Definition: HCalSD.h:97
double birk1
Definition: HCalSD.h:100
std::unique_ptr< HFShowerPMT > showerPMT
Definition: HCalSD.h:86
TH1F * hzvhad
Definition: HCalSD.h:114
bool isItConicalBundle(const G4LogicalVolume *)
Definition: HCalSD.cc:660
double eminHitHF
Definition: HCalSD.h:102
Definition: HCalSD.h:38
std::unique_ptr< HcalNumberingFromDDD > numberingFromDDD
Definition: HCalSD.h:81
TH1F * hit_[9]
Definition: HCalSD.h:114
bool agingFlagHB
Definition: HCalSD.h:97
bool testNumber
Definition: HCalSD.h:99
~HCalSD() override=default
bool useBirk
Definition: HCalSD.h:98
std::unique_ptr< HFDarkening > m_HFDarkening
Definition: HCalSD.h:93
double birk3
Definition: HCalSD.h:100
const HBHEDarkening * m_HEDarkening
Definition: HCalSD.h:92
bool applyFidCut
Definition: HCalSD.h:101
TH1F * dist_[9]
Definition: HCalSD.h:114
std::unique_ptr< HcalTestNS > m_HcalTestNS
Definition: HCalSD.h:94
void plotHF(const G4ThreeVector &pos, bool emType)
Definition: HCalSD.cc:1032
std::vector< int > hfLevels
Definition: HCalSD.h:107
void getFromParam(const G4Step *step, bool &isKilled)
Definition: HCalSD.cc:777
void getFromHFLibrary(const G4Step *step, bool &isKilled)
Definition: HCalSD.cc:692
uint32_t setDetUnitId(const G4Step *step) override
Definition: HCalSD.cc:538
bool isItStraightBundle(const G4LogicalVolume *)
Definition: HCalSD.cc:652
bool testNS_
Definition: HCalSD.h:99
bool getFromLibrary(const G4Step *) override
Definition: HCalSD.cc:349
double eminHitHO
Definition: HCalSD.h:102
TH1F * hzvem
Definition: HCalSD.h:114
std::vector< const G4LogicalVolume * > pmtLV
Definition: HCalSD.h:112
bool filterHit(CaloG4Hit *, double) override
Definition: HCalSD.cc:561
bool neutralDensity
Definition: HCalSD.h:99
bool isHF
Definition: HCalSD.h:96
void getHitFibreBundle(const G4Step *step, bool type)
Definition: HCalSD.cc:866
double getEnergyDeposit(const G4Step *) override
Definition: HCalSD.cc:406
std::unique_ptr< HFShowerFibreBundle > showerBundle
Definition: HCalSD.h:87
bool isItPMT(const G4LogicalVolume *)
Definition: HCalSD.cc:644
bool isItinFidVolume(const G4ThreeVector &)
Definition: HCalSD.cc:676
void getHitPMT(const G4Step *step)
Definition: HCalSD.cc:807
std::unique_ptr< HFShowerLibrary > showerLibrary
Definition: HCalSD.h:83
double layerWeight(int, const G4ThreeVector &, int, int)
Definition: HCalSD.cc:952
int depth_
Definition: HCalSD.h:105
std::map< uint32_t, double > layerWeights
Definition: HCalSD.h:113
bool isItHF(const G4Step *)
Definition: HCalSD.cc:607
bool isItScintillator(const G4Material *)
Definition: HCalSD.cc:668
std::unique_ptr< HcalNumberingScheme > numberingScheme
Definition: HCalSD.h:82
std::unique_ptr< HFShower > hfshower
Definition: HCalSD.h:84
std::vector< const G4LogicalVolume * > fibreLV
Definition: HCalSD.h:112
step
Definition: StallMonitor.cc:94
std::vector< const G4LogicalVolume * > fibre1LV
Definition: HCalSD.h:112
tmp
align.sh
Definition: createJobs.py:716
std::vector< const G4Material * > materials
Definition: HCalSD.h:111
void initRun() override
Definition: HCalSD.cc:559
const HBHEDarkening * m_HBDarkening
Definition: HCalSD.h:91
void update(const BeginOfRun *) override
This routine will be called when the appropriate signal arrives.
Definition: CaloSD.cc:456
bool useHF
Definition: HCalSD.h:101
bool isItFibre(const G4LogicalVolume *)
Definition: HCalSD.cc:628
const HcalDDDSimulationConstants * hcalSimConstants_
Definition: HCalSD.h:90