CMS 3D CMS Logo

HFShowerLibrary.h
Go to the documentation of this file.
1 #ifndef SimG4CMS_HFShowerLibrary_h
2 #define SimG4CMS_HFShowerLibrary_h 1
3 // File: HFShowerLibrary.h
5 // Description: Gets information from a shower library
7 
14 
15 #include "G4ThreeVector.hh"
16 
17 //ROOT
18 #include "TROOT.h"
19 #include "TFile.h"
20 #include "TTree.h"
21 
22 #include <string>
23 #include <memory>
24 
25 class G4Step;
26 class G4ParticleTable;
27 
29 public:
30  //Constructor and Destructor
32  const HcalDDDSimConstants *hcons,
33  const HcalSimulationParameters *hps,
34  edm::ParameterSet const &p);
36 
37 public:
38  struct Hit {
39  Hit() {}
40  G4ThreeVector position;
41  int depth;
42  double time;
43  };
44 
45  std::vector<Hit> getHits(const G4Step *aStep, bool &ok, double weight, bool onlyLong = false);
46  std::vector<Hit> fillHits(const G4ThreeVector &p,
47  const G4ThreeVector &v,
48  int parCode,
49  double parEnergy,
50  bool &ok,
51  double weight,
52  double time,
53  bool onlyLong = false);
54 
55  struct Params {
56  double probMax_;
57  double backProb_;
58  double dphi_;
60  bool verbose_;
62  };
63  struct FileParams {
69  bool cacheBranches_ = false;
70  };
72 
73 private:
75  const HcalSimulationParameters *hps,
76  edm::ParameterSet const &hfShower,
77  edm::ParameterSet const &hfShowerLibrary);
78  bool rInside(double r) const;
79  HFShowerPhotonCollection getRecord(int, int) const;
80 
81  struct VersionInfo {
82  float libVers_;
83  float listVersion_;
84  };
85  VersionInfo loadEventInfo(TBranch *, int fileVersion);
88  void storePhoton(HFShowerPhoton const &iPhoton, HFShowerPhotonCollection &iPhotons) const;
89 
90  enum class FileFormat { kOld, kNew, kNewV3 };
91 
92  struct BranchCache;
93  struct BranchReader {
94  BranchReader() : branch_(nullptr), offset_(0), format_(FileFormat::kOld) {}
95  BranchReader(TBranch *iBranch, FileFormat iFormat, size_t iReadOffset, size_t maxRecordsToCache)
96  : branch_(iBranch), offset_(iReadOffset), format_(iFormat) {
97  if (0 < maxRecordsToCache) {
98  doCaching(maxRecordsToCache);
99  branch_ = nullptr;
100  }
101  }
102 
104 
105  std::size_t numberOfRecords() const;
106 
107  private:
108  static HFShowerPhotonCollection getRecordOldForm(TBranch *, int iEntry);
109  static HFShowerPhotonCollection getRecordNewForm(TBranch *, int iEntry);
110  static HFShowerPhotonCollection getRecordNewFormV3(TBranch *, int iEntry);
111  void doCaching(size_t maxRecords);
112 
113  static std::shared_ptr<BranchCache const> makeCache(BranchReader &,
114  size_t maxRecordsToCache,
115  std::string const &iFileName,
116  std::string const &iBranchName);
117 
118  TBranch *branch_;
119  std::shared_ptr<BranchCache const> cache_;
120  size_t offset_;
122  };
123 
124  struct BranchCache {
125  explicit BranchCache(BranchReader &, size_t maxRecordsToCache);
126  [[nodiscard]] HFShowerPhotonCollection getRecord(int) const;
127 
128  private:
130  std::vector<std::size_t> offsets_;
131  };
132 
134  std::unique_ptr<TFile> hf_;
135 
138 
141  std::vector<double> pmom_;
142 
145  double dphi_, rMin_, rMax_;
146  std::vector<double> gpar_;
147 };
148 #endif
static HFShowerPhotonCollection getRecordOldForm(TBranch *, int iEntry)
BranchCache(BranchReader &, size_t maxRecordsToCache)
std::vector< Hit > getHits(const G4Step *aStep, bool &ok, double weight, bool onlyLong=false)
void storePhoton(HFShowerPhoton const &iPhoton, HFShowerPhotonCollection &iPhotons) const
HFShowerPhotonCollection getRecord(int, int) const
std::size_t numberOfRecords() const
Definition: weight.py:1
HFShowerPhotonCollection extrapolate(int, double)
std::vector< double > gpar_
void doCaching(size_t maxRecords)
BranchReader hadBranch_
std::unique_ptr< TFile > hf_
VersionInfo loadEventInfo(TBranch *, int fileVersion)
std::vector< HFShowerPhoton > HFShowerPhotonCollection
HFShowerPhotonCollection interpolate(int, double)
std::shared_ptr< BranchCache const > cache_
static std::shared_ptr< BranchCache const > makeCache(BranchReader &, size_t maxRecordsToCache, std::string const &iFileName, std::string const &iBranchName)
std::vector< double > pmom_
G4ThreeVector position
HFShowerPhotonCollection getRecord(int) const
HFShowerLibrary(const std::string &name, const HcalDDDSimConstants *hcons, const HcalSimulationParameters *hps, edm::ParameterSet const &p)
BranchReader emBranch_
static HFShowerPhotonCollection getRecordNewForm(TBranch *, int iEntry)
std::vector< std::size_t > offsets_
static HFShowerPhotonCollection getRecordNewFormV3(TBranch *, int iEntry)
HFShowerPhotonCollection photons_
BranchReader(TBranch *iBranch, FileFormat iFormat, size_t iReadOffset, size_t maxRecordsToCache)
std::vector< Hit > fillHits(const G4ThreeVector &p, const G4ThreeVector &v, int parCode, double parEnergy, bool &ok, double weight, double time, bool onlyLong=false)
bool rInside(double r) const
HFShowerPhotonCollection getRecord(int) const