CMS 3D CMS Logo

CaloSteppingAction.h
Go to the documentation of this file.
1 #ifndef SimG4CMS_CaloSteppingAction_H
2 #define SimG4CMS_CaloSteppingAction_H
3 //#define HcalNumberingTest
4 
12 
20 
23 
27 
31 #ifdef HcalNumberingTest
34 #endif
35 
36 #include "G4LogicalVolume.hh"
37 #include "G4Region.hh"
38 #include "G4Step.hh"
39 #include "G4UserSteppingAction.hh"
40 #include "G4VPhysicalVolume.hh"
41 #include "G4VTouchable.hh"
42 #include "G4Track.hh"
43 
44 #include <algorithm>
45 #include <iostream>
46 #include <memory>
47 #include <string>
48 #include <vector>
49 
51  public Observer<const BeginOfJob *>,
52  public Observer<const BeginOfRun *>,
53  public Observer<const BeginOfEvent *>,
54  public Observer<const EndOfEvent *>,
55  public Observer<const G4Step *> {
56 public:
58  ~CaloSteppingAction() override;
59 
60  void produce(edm::Event &, const edm::EventSetup &) override;
61 
62 private:
63  void fillHits(edm::PCaloHitContainer &cc, int type);
64  // observer classes
65  void update(const BeginOfJob *job) override;
66  void update(const BeginOfRun *run) override;
67  void update(const BeginOfEvent *evt) override;
68  void update(const G4Step *step) override;
69  void update(const EndOfEvent *evt) override;
70 
71  void NaNTrap(const G4Step *) const;
72  uint32_t getDetIDHC(int det, int lay, int depth, const math::XYZVectorD &pos) const;
73  void fillHit(uint32_t id, double dE, double time, int primID, uint16_t depth, double em, int flag);
74  uint16_t getDepth(bool flag, double crystalDepth, double radl) const;
75  double curve_LY(double crystalLength, double crystalDepth) const;
76  double getBirkL3(double dE, double step, double chg, double dens) const;
77  double getBirkHC(double dE, double step, double chg, double dens) const;
78  void saveHits(int flag);
79 
80  static const int nSD_ = 3;
81  std::unique_ptr<EcalBarrelNumberingScheme> ebNumberingScheme_;
82  std::unique_ptr<EcalEndcapNumberingScheme> eeNumberingScheme_;
83  std::unique_ptr<HcalNumberingFromPS> hcNumberingPS_;
84 #ifdef HcalNumberingTest
85  std::unique_ptr<HcalNumberingFromDDD> hcNumbering_;
86 #endif
87  std::unique_ptr<HcalNumberingScheme> hcNumberingScheme_;
88  std::unique_ptr<CaloSlaveSD> slave_[nSD_];
89 
90  std::vector<std::string> nameEBSD_, nameEESD_, nameHCSD_;
91  std::vector<std::string> nameHitC_;
92  std::vector<const G4LogicalVolume *> volEBSD_, volEESD_, volHCSD_;
93  std::map<const G4LogicalVolume *, double> xtalMap_;
97  double birkC3HC_;
98  std::map<std::pair<int, CaloHitID>, CaloGVHit> hitMap_[nSD_];
99 };
100 
101 #endif
std::map< std::pair< int, CaloHitID >, CaloGVHit > hitMap_[nSD_]
std::vector< const G4LogicalVolume * > volEBSD_
type
Definition: HCALResponse.h:21
std::vector< PCaloHit > PCaloHitContainer
double getBirkL3(double dE, double step, double chg, double dens) const
const float chg[109]
Definition: CoreSimTrack.cc:5
void produce(edm::Event &, const edm::EventSetup &) override
uint32_t getDetIDHC(int det, int lay, int depth, const math::XYZVectorD &pos) const
double getBirkHC(double dE, double step, double chg, double dens) const
double curve_LY(double crystalLength, double crystalDepth) const
void fillHits(edm::PCaloHitContainer &cc, int type)
std::vector< const G4LogicalVolume * > volEESD_
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< double > > XYZVectorD
spatial vector with cartesian internal representation
Definition: Vector3D.h:8
std::vector< std::string > nameHitC_
CaloSteppingAction(const edm::ParameterSet &p)
std::vector< const G4LogicalVolume * > volHCSD_
std::vector< std::string > nameHCSD_
std::map< const G4LogicalVolume *, double > xtalMap_
std::unique_ptr< CaloSlaveSD > slave_[nSD_]
std::vector< std::string > nameEBSD_
void fillHit(uint32_t id, double dE, double time, int primID, uint16_t depth, double em, int flag)
void NaNTrap(const G4Step *) const
std::unique_ptr< HcalNumberingFromPS > hcNumberingPS_
std::unique_ptr< HcalNumberingScheme > hcNumberingScheme_
void update(const BeginOfJob *job) override
This routine will be called when the appropriate signal arrives.
std::unique_ptr< EcalBarrelNumberingScheme > ebNumberingScheme_
std::unique_ptr< EcalEndcapNumberingScheme > eeNumberingScheme_
std::vector< std::string > nameEESD_
step
Definition: StallMonitor.cc:94
static const int nSD_
uint16_t getDepth(bool flag, double crystalDepth, double radl) const
~CaloSteppingAction() override