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 
21 
22 #include "G4String.hh"
23 #include <map>
24 #include <string>
25 #include <TH1F.h>
26 
27 class DDCompactView;
28 class DDFilteredView;
29 class G4LogicalVolume;
30 class G4Material;
31 class G4Step;
32 
33 class HCalSD : public CaloSD {
34 
35 public:
36 
37  HCalSD(G4String , const DDCompactView &, const SensitiveDetectorCatalog &,
38  edm::ParameterSet const &, const SimTrackManager*);
39  virtual ~HCalSD();
40  virtual bool ProcessHits(G4Step * , G4TouchableHistory * );
41  virtual double getEnergyDeposit(G4Step* );
42  virtual uint32_t setDetUnitId(G4Step* step);
44 
45 protected:
46 
47  virtual void initRun();
48  virtual bool filterHit(CaloG4Hit*, double);
49 
50 private:
51 
52  uint32_t setDetUnitId(int, const G4ThreeVector&, int, int);
53  std::vector<double> getDDDArray(const std::string&,
54  const DDsvalues_type&);
55  std::vector<G4String> getNames(DDFilteredView&);
56  bool isItHF(G4Step *);
57  bool isItHF(G4String);
58  bool isItFibre(G4LogicalVolume*);
59  bool isItFibre(G4String);
60  bool isItPMT(G4LogicalVolume*);
61  bool isItStraightBundle(G4LogicalVolume*);
62  bool isItConicalBundle(G4LogicalVolume*);
63  bool isItScintillator(G4Material*);
64  bool isItinFidVolume (G4ThreeVector&);
65  void getFromLibrary(G4Step * step, double weight);
66  void hitForFibre(G4Step * step, double weight);
67  void getFromParam(G4Step * step, double weight);
68  void getHitPMT(G4Step * step);
69  void getHitFibreBundle(G4Step * step, bool type);
70  int setTrackID(G4Step * step);
72  double layerWeight(int, const G4ThreeVector&, int, int);
73  void plotProfile(G4Step* step, const G4ThreeVector& pos,
74  double edep, double time, int id);
75  void plotHF(G4ThreeVector& pos, bool emType);
76 
87  double birk1, birk2, birk3, betaThr;
90  double deliveredLumi;
91  G4int mumPDG, mupPDG;
92  std::vector<double> layer0wt, gpar;
93  std::vector<int> hfLevels;
94  std::vector<G4String> hfNames, fibreNames, matNames;
95  std::vector<G4Material*> materials;
96  std::vector<G4LogicalVolume*> hfLV, fibreLV, pmtLV, fibre1LV, fibre2LV;
97  std::map<uint32_t,double> layerWeights;
98  TH1F *hit_[9], *time_[9], *dist_[9], *hzvem, *hzvhad;
99 
100 };
101 
102 #endif // HCalSD_h
HCalSD(G4String, const DDCompactView &, const SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
Definition: HCalSD.cc:38
type
Definition: HCALResponse.h:21
std::vector< double > layer0wt
Definition: HCalSD.h:92
bool useParam
Definition: HCalSD.h:88
double eminHitHE
Definition: HCalSD.h:89
HFShowerParam * showerParam
Definition: HCalSD.h:81
TH1F * time_[9]
Definition: HCalSD.h:98
void hitForFibre(G4Step *step, double weight)
Definition: HCalSD.cc:833
void getFromParam(G4Step *step, double weight)
Definition: HCalSD.cc:885
bool useLayerWt
Definition: HCalSD.h:86
Definition: CaloSD.h:42
bool isItStraightBundle(G4LogicalVolume *)
Definition: HCalSD.cc:735
virtual bool ProcessHits(G4Step *, G4TouchableHistory *)
Definition: HCalSD.cc:395
std::vector< double > gpar
Definition: HCalSD.h:92
bool useFibreBundle
Definition: HCalSD.h:86
bool isItinFidVolume(G4ThreeVector &)
Definition: HCalSD.cc:753
double betaThr
Definition: HCalSD.h:87
double deliveredLumi
Definition: HCalSD.h:90
HFShowerFibreBundle * showerBundle
Definition: HCalSD.h:83
std::vector< G4LogicalVolume * > fibre1LV
Definition: HCalSD.h:96
std::vector< G4LogicalVolume * > fibreLV
Definition: HCalSD.h:96
double eminHitHB
Definition: HCalSD.h:89
bool useShowerLibrary
Definition: HCalSD.h:88
std::vector< double > getDDDArray(const std::string &, const DDsvalues_type &)
Definition: HCalSD.cc:654
virtual uint32_t setDetUnitId(G4Step *step)
Definition: HCalSD.cc:589
double birk2
Definition: HCalSD.h:87
int setTrackID(G4Step *step)
Definition: HCalSD.cc:1054
void getFromLibrary(G4Step *step, double weight)
Definition: HCalSD.cc:770
bool usePMTHit
Definition: HCalSD.h:86
void setNumberingScheme(HcalNumberingScheme *)
Definition: HCalSD.cc:602
type of data representation of DDCompactView
Definition: DDCompactView.h:77
double birk1
Definition: HCalSD.h:87
bool isItHF(G4Step *)
Definition: HCalSD.cc:699
TH1F * hzvhad
Definition: HCalSD.h:98
double eminHitHF
Definition: HCalSD.h:89
bool isItConicalBundle(G4LogicalVolume *)
Definition: HCalSD.cc:741
Definition: HCalSD.h:33
TH1F * hit_[9]
Definition: HCalSD.h:98
virtual void initRun()
Definition: HCalSD.cc:610
bool testNumber
Definition: HCalSD.h:86
std::vector< G4LogicalVolume * > pmtLV
Definition: HCalSD.h:96
std::vector< G4LogicalVolume * > hfLV
Definition: HCalSD.h:96
bool isItScintillator(G4Material *)
Definition: HCalSD.cc:747
HFShower * hfshower
Definition: HCalSD.h:80
bool useBirk
Definition: HCalSD.h:86
double birk3
Definition: HCalSD.h:87
virtual double getEnergyDeposit(G4Step *)
Definition: HCalSD.cc:499
void plotHF(G4ThreeVector &pos, bool emType)
Definition: HCalSD.cc:1168
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:19
bool applyFidCut
Definition: HCalSD.h:88
TH1F * dist_[9]
Definition: HCalSD.h:98
virtual ~HCalSD()
Definition: HCalSD.cc:382
HcalNumberingScheme * numberingScheme
Definition: HCalSD.h:78
void getHitPMT(G4Step *step)
Definition: HCalSD.cc:924
HFDarkening * m_HFDarkening
Definition: HCalSD.h:85
std::vector< int > hfLevels
Definition: HCalSD.h:93
virtual std::vector< std::string > getNames()
void plotProfile(G4Step *step, const G4ThreeVector &pos, double edep, double time, int id)
Definition: HCalSD.cc:1121
HEDarkening * m_HEDarkening
Definition: HCalSD.h:84
bool isItPMT(G4LogicalVolume *)
Definition: HCalSD.cc:729
G4int mupPDG
Definition: HCalSD.h:91
double eminHitHO
Definition: HCalSD.h:89
TH1F * hzvem
Definition: HCalSD.h:98
void readWeightFromFile(std::string)
Definition: HCalSD.cc:1074
void getHitFibreBundle(G4Step *step, bool type)
Definition: HCalSD.cc:991
HcalNumberingFromDDD * numberingFromDDD
Definition: HCalSD.h:77
double layerWeight(int, const G4ThreeVector &, int, int)
Definition: HCalSD.cc:1100
std::vector< G4String > fibreNames
Definition: HCalSD.h:94
std::vector< G4Material * > materials
Definition: HCalSD.h:95
virtual bool filterHit(CaloG4Hit *, double)
Definition: HCalSD.cc:624
std::vector< G4String > hfNames
Definition: HCalSD.h:94
G4int mumPDG
Definition: HCalSD.h:91
std::map< uint32_t, double > layerWeights
Definition: HCalSD.h:97
HFShowerPMT * showerPMT
Definition: HCalSD.h:82
int weight
Definition: histoStyle.py:50
HFShowerLibrary * showerLibrary
Definition: HCalSD.h:79
std::vector< G4LogicalVolume * > fibre2LV
Definition: HCalSD.h:96
std::vector< G4String > matNames
Definition: HCalSD.h:94
bool isItFibre(G4LogicalVolume *)
Definition: HCalSD.cc:717
bool useHF
Definition: HCalSD.h:88