CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 
13 
14 #include "G4ParticleTable.hh"
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 DDCompactView;
26 class G4Step;
27 
29 
30 public:
31 
32  //Constructor and Destructor
34  edm::ParameterSet const & p);
36 
37 public:
38 
39  struct Hit {
40  Hit() {}
41  G4ThreeVector position;
42  int depth;
43  double time;
44  };
45 
46  void initRun(G4ParticleTable * theParticleTable);
47  std::vector<Hit> getHits(G4Step * aStep, bool &ok, double weight,
48  bool onlyLong=false);
49 
50 protected:
51 
52  bool rInside(double r);
53  void getRecord(int, int);
54  void loadEventInfo(TBranch *);
55  void interpolate(int, double);
56  void extrapolate(int, double);
57  void storePhoton(int j);
58  std::vector<double> getDDDArray(const std::string&, const DDsvalues_type&,
59  int&);
60 
61 private:
62 
64  TFile * hf;
65  TBranch *emBranch, *hadBranch;
66 
70  std::vector<double> pmom;
71 
72  double probMax, backProb;
73  double dphi, rMin, rMax;
74  std::vector<double> gpar;
75 
79 
80  int npe;
84 
85 };
86 #endif
void getRecord(int, int)
std::vector< Hit > getHits(G4Step *aStep, bool &ok, double weight, bool onlyLong=false)
HFShowerPhotonCollection pe
type of data representation of DDCompactView
Definition: DDCompactView.h:77
HFShowerLibrary(std::string &name, const DDCompactView &cpv, edm::ParameterSet const &p)
std::vector< HFShowerPhoton > HFShowerPhotonCollection
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
void loadEventInfo(TBranch *)
int j
Definition: DBlmapReader.cc:9
std::vector< double > gpar
std::vector< double > getDDDArray(const std::string &, const DDsvalues_type &, int &)
HFShowerPhotonCollection * photo
TBranch * emBranch
std::vector< double > pmom
G4ThreeVector position
void initRun(G4ParticleTable *theParticleTable)
void extrapolate(int, double)
void storePhoton(int j)
int weight
Definition: histoStyle.py:50
HFShowerPhotonCollection photon
bool rInside(double r)
void interpolate(int, double)
TBranch * hadBranch