CMS 3D CMS Logo

DreamSD.h
Go to the documentation of this file.
1 #ifndef SimG4CMS_DreamSD_h
2 #define SimG4CMS_DreamSD_h
3 
6 
7 #include "G4PhysicsOrderedFreeVector.hh"
8 
9 #include <map>
10 
11 const int MAXPHOTONS = 500; // Maximum number of photons we can store
12 
13 class G4LogicalVolume;
14 
15 class DreamSD : public CaloSD {
16 public:
17  DreamSD(const std::string &,
18  const DDCompactView &,
20  edm::ParameterSet const &,
21  const SimTrackManager *);
22  ~DreamSD() override {}
23 
24  uint32_t setDetUnitId(const G4Step *) override;
25 
26 protected:
27  double getEnergyDeposit(const G4Step *) override;
28  void initRun() override;
29 
30 private:
31  typedef std::pair<double, double> Doubles;
32  typedef std::map<G4LogicalVolume *, Doubles> DimensionMap;
33 
34  void initMap(const std::string &, const DDCompactView &);
35  double curve_LY(const G4Step *, int);
36  double crystalLength(G4LogicalVolume *) const;
37  double crystalWidth(G4LogicalVolume *) const;
38 
40  double cherenkovDeposit_(const G4Step *aStep);
42  double getAverageNumberOfPhotons_(const double charge,
43  const double beta,
44  const G4Material *aMaterial,
45  const G4MaterialPropertyVector *rIndex);
47  double getPhotonEnergyDeposit_(const G4ParticleMomentum &p, const G4ThreeVector &x, const G4Step *aStep);
49  bool setPbWO2MaterialProperties_(G4Material *aMaterial);
50 
52  double birk1, birk2, birk3;
53  double slopeLY;
54  DimensionMap xtalLMap; // Store length and width
55 
56  int side;
57 
59  std::unique_ptr<G4PhysicsOrderedFreeVector> chAngleIntegrals_;
60  G4MaterialPropertiesTable *materialPropertiesTable;
61 
62  int nphotons_;
63 };
64 
65 #endif // DreamSD_h
bool readBothSide_
Definition: DreamSD.h:51
double getPhotonEnergyDeposit_(const G4ParticleMomentum &p, const G4ThreeVector &x, const G4Step *aStep)
Returns energy deposit for a given photon.
Definition: DreamSD.cc:451
Definition: CaloSD.h:38
uint32_t setDetUnitId(const G4Step *) override
Definition: DreamSD.cc:92
double birk1
Definition: DreamSD.h:52
bool setPbWO2MaterialProperties_(G4Material *aMaterial)
Sets material properties at run-time...
Definition: DreamSD.cc:371
double birk3
Definition: DreamSD.h:52
void initRun() override
Definition: DreamSD.cc:74
Compact representation of the geometrical detector hierarchy.
Definition: DDCompactView.h:80
double cherenkovDeposit_(const G4Step *aStep)
Returns the total energy due to Cherenkov radiation.
Definition: DreamSD.cc:191
int nphotons_
Definition: DreamSD.h:62
int side
Definition: DreamSD.h:56
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:304
~DreamSD() override
Definition: DreamSD.h:22
double crystalWidth(G4LogicalVolume *) const
Definition: DreamSD.cc:180
double birk2
Definition: DreamSD.h:52
std::pair< double, double > Doubles
Definition: DreamSD.h:31
const int MAXPHOTONS
Definition: DreamSD.h:11
std::map< G4LogicalVolume *, Doubles > DimensionMap
Definition: DreamSD.h:32
G4MaterialPropertiesTable * materialPropertiesTable
Definition: DreamSD.h:60
double curve_LY(const G4Step *, int)
Definition: DreamSD.cc:146
DimensionMap xtalLMap
Definition: DreamSD.h:54
double crystalLength(G4LogicalVolume *) const
Definition: DreamSD.cc:171
std::unique_ptr< G4PhysicsOrderedFreeVector > chAngleIntegrals_
Table of Cherenkov angle integrals vs photon momentum.
Definition: DreamSD.h:59
double getEnergyDeposit(const G4Step *) override
Definition: DreamSD.cc:55
bool doCherenkov_
Definition: DreamSD.h:51
double slopeLY
Definition: DreamSD.h:53
DreamSD(const std::string &, const DDCompactView &, const SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
Definition: DreamSD.cc:27
void initMap(const std::string &, const DDCompactView &)
Definition: DreamSD.cc:104
bool useBirk
Definition: DreamSD.h:51