CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HFShower.h
Go to the documentation of this file.
1 #ifndef SimG4CMS_HFShower_h
2 #define SimG4CMS_HFShower_h
3 // File: HFShower.h
5 // Description: Generates hits for HF with Cerenkov photon code
7 
13 
14 #include "G4ThreeVector.hh"
15 #include "G4String.hh"
16 
17 class DDCompactView;
18 class G4Step;
19 
20 #include <vector>
21 
22 class HFShower {
23 
24 public:
25 
26  HFShower(std::string & name, const DDCompactView & cpv,
27  edm::ParameterSet const & p, int chk=0);
28  virtual ~HFShower();
29 
30 public:
31 
32  struct Hit {
33  Hit() {}
34  int depth;
35  double time;
36  double wavelength;
37  double momentum;
38  G4ThreeVector position;
39  };
40 
41  void initRun(G4ParticleTable *, HcalDDDSimConstants*);
42  std::vector<Hit> getHits(G4Step * aStep, double weight);
43  std::vector<Hit> getHits(G4Step * aStep, bool forLibrary);
44  std::vector<Hit> getHits(G4Step * aStep, bool forLibraryProducer, double zoffset);
45 
46 
47 private:
48 
49  std::vector<double> getDDDArray(const std::string &, const DDsvalues_type &, int &);
51 
52 private:
53 
56 
57  int chkFibre;
58  double probMax;
59  std::vector<double> gpar;
60 
61 };
62 
63 #endif // HFShower_h
double wavelength
Definition: HFShower.h:36
HFFibre * fibre
Definition: HFShower.h:55
double momentum
Definition: HFShower.h:37
void initRun(G4ParticleTable *, HcalDDDSimConstants *)
Definition: HFShower.cc:454
bool applyFidCut
Definition: HFShower.h:50
G4ThreeVector position
Definition: HFShower.h:38
type of data representation of DDCompactView
Definition: DDCompactView.h:77
HFCherenkov * cherenkov
Definition: HFShower.h:54
int chkFibre
Definition: HFShower.h:57
virtual ~HFShower()
Definition: HFShower.h:41
double probMax
Definition: HFShower.h:58
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
std::vector< Hit > getHits(G4Step *aStep, double weight)
Definition: HFShower.cc:47
std::vector< double > gpar
Definition: HFShower.h:59
double time
Definition: HFShower.h:35
HFShower(const RandomEngineAndDistribution *engine, HDShowerParametrization *myParam, EcalHitMaker *myGrid, HcalHitMaker *myHcalHitMaker, int onECAL, double epart)
Definition: HFShower.cc:30
std::vector< double > getDDDArray(const std::string &, const DDsvalues_type &, int &)
Definition: HFShower.cc:416