CMS 3D CMS Logo

DreamSD.h
Go to the documentation of this file.
1 #ifndef SimG4CMS_DreamSD_h
2 #define SimG4CMS_DreamSD_h
3 
5 
7 
8 #include "G4String.hh"
9 #include "G4PhysicsOrderedFreeVector.hh"
10 
11 #include <map>
12 
13 const int MAXPHOTONS = 500; // Maximum number of photons we can store
14 
15 class G4LogicalVolume;
16 
17 class TTree;
18 
19 class DreamSD : public CaloSD {
20 
21 public:
22 
24  edm::ParameterSet const &, const SimTrackManager*);
25  ~DreamSD() override {}
26  // bool ProcessHits(G4Step * step,G4TouchableHistory * tHistory) override;
27  uint32_t setDetUnitId(const G4Step*) override;
28 
29 protected:
30 
31  // G4bool getStepInfo(G4Step* aStep) override;
32  double getEnergyDeposit(const G4Step* ) override;
33  void initRun() override;
34 
35 private:
36 
37  typedef std::pair<double,double> Doubles;
38  typedef std::map<G4LogicalVolume*,Doubles> DimensionMap;
39 
40  void initMap(const std::string&, const DDCompactView &);
41  double curve_LY(const G4Step*, int);
42  const double crystalLength(G4LogicalVolume*) const;
43  const double crystalWidth(G4LogicalVolume*) const;
44 
46  double cherenkovDeposit_(const G4Step* aStep );
48  double getAverageNumberOfPhotons_(const double charge,
49  const double beta,
50  const G4Material* aMaterial,
51  const G4MaterialPropertyVector* rIndex );
53  double getPhotonEnergyDeposit_( const G4ParticleMomentum& p,
54  const G4ThreeVector& x,
55  const G4Step* aStep );
57  bool setPbWO2MaterialProperties_(G4Material* aMaterial );
58 
60  double birk1, birk2, birk3;
61  double slopeLY;
62  DimensionMap xtalLMap; // Store length and width
63 
64  int side;
65 
67  std::unique_ptr<G4PhysicsOrderedFreeVector> chAngleIntegrals_;
68  G4MaterialPropertiesTable* materialPropertiesTable;
69  // Histogramming
70  TTree* ntuple_;
71  int nphotons_;
74 
75 };
76 
77 #endif // DreamSD_h
bool readBothSide_
Definition: DreamSD.h:59
double getPhotonEnergyDeposit_(const G4ParticleMomentum &p, const G4ThreeVector &x, const G4Step *aStep)
Returns energy deposit for a given photon.
Definition: DreamSD.cc:494
std::pair< double, double > Doubles
Definition: DreamSD.h:37
Definition: CaloSD.h:37
float y_[MAXPHOTONS]
Definition: DreamSD.h:73
uint32_t setDetUnitId(const G4Step *) override
Definition: DreamSD.cc:118
double birk1
Definition: DreamSD.h:60
bool setPbWO2MaterialProperties_(G4Material *aMaterial)
Sets material properties at run-time...
Definition: DreamSD.cc:432
double birk3
Definition: DreamSD.h:60
void initRun() override
Definition: DreamSD.cc:98
float px_[MAXPHOTONS]
Definition: DreamSD.h:72
Compact representation of the geometrical detector hierarchy.
Definition: DDCompactView.h:80
TTree * ntuple_
Definition: DreamSD.h:70
std::map< G4LogicalVolume *, Doubles > DimensionMap
Definition: DreamSD.h:38
float z_[MAXPHOTONS]
Definition: DreamSD.h:73
double cherenkovDeposit_(const G4Step *aStep)
Returns the total energy due to Cherenkov radiation.
Definition: DreamSD.cc:228
int nphotons_
Definition: DreamSD.h:71
float pz_[MAXPHOTONS]
Definition: DreamSD.h:72
int side
Definition: DreamSD.h:64
double getAverageNumberOfPhotons_(const double charge, const double beta, const G4Material *aMaterial, const G4MaterialPropertyVector *rIndex)
Returns average number of photons created by track.
Definition: DreamSD.cc:360
~DreamSD() override
Definition: DreamSD.h:25
double birk2
Definition: DreamSD.h:60
const int MAXPHOTONS
Definition: DreamSD.h:13
G4MaterialPropertiesTable * materialPropertiesTable
Definition: DreamSD.h:68
double curve_LY(const G4Step *, int)
Definition: DreamSD.cc:174
const double crystalWidth(G4LogicalVolume *) const
Definition: DreamSD.cc:215
DimensionMap xtalLMap
Definition: DreamSD.h:62
float x_[MAXPHOTONS]
Definition: DreamSD.h:73
std::unique_ptr< G4PhysicsOrderedFreeVector > chAngleIntegrals_
Table of Cherenkov angle integrals vs photon momentum.
Definition: DreamSD.h:67
double getEnergyDeposit(const G4Step *) override
Definition: DreamSD.cc:79
bool doCherenkov_
Definition: DreamSD.h:59
double slopeLY
Definition: DreamSD.h:61
DreamSD(const std::string &, const DDCompactView &, const SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
Definition: DreamSD.cc:29
float py_[MAXPHOTONS]
Definition: DreamSD.h:72
const double crystalLength(G4LogicalVolume *) const
Definition: DreamSD.cc:205
void initMap(const std::string &, const DDCompactView &)
Definition: DreamSD.cc:129
bool useBirk
Definition: DreamSD.h:59