CMS 3D CMS Logo

CaloTrkProcessing.h
Go to the documentation of this file.
1 #ifndef SimG4CMS_CaloTrkProcessing_H
2 #define SimG4CMS_CaloTrkProcessing_H
3 
7 
9 
10 #include "G4VTouchable.hh"
11 
12 #include <map>
13 #include <vector>
14 #include <string>
15 #include <iostream>
16 
17 class SimTrackManager;
18 class BeginOfEvent;
19 class G4LogicalVolume;
20 class G4Step;
21 class SimTrackManager;
22 
24  public Observer<const BeginOfEvent *>,
25  public Observer<const G4Step *> {
26 
27 public:
28 
29  CaloTrkProcessing(const std::string& aSDname, const DDCompactView & cpv,
30  const SensitiveDetectorCatalog & clg,
31  edm::ParameterSet const & p, const SimTrackManager*);
32  ~CaloTrkProcessing() override;
33  void Initialize(G4HCofThisEvent * ) override {}
34  void clearHits() override {}
35  bool ProcessHits(G4Step * , G4TouchableHistory * ) override {
36  return true;
37  }
38  uint32_t setDetUnitId(const G4Step * step) override {return 0;}
39  void EndOfEvent(G4HCofThisEvent * ) override {}
40  void fillHits(edm::PCaloHitContainer&, const std::string& ) override {}
41 
42 private:
43 
44  void update(const BeginOfEvent * evt) override;
45  void update(const G4Step *) override;
46  std::vector<std::string> getNames(G4String, const DDsvalues_type&);
47  std::vector<double> getNumbers(G4String, const DDsvalues_type&);
48  int isItCalo(const G4VTouchable*);
49  int isItInside(const G4VTouchable*, int, int);
50 
51 
52  struct Detector {
53  Detector() {}
55  G4LogicalVolume* lv;
56  int level;
57  std::vector<std::string> fromDets;
58  std::vector<G4LogicalVolume*> fromDetL;
59  std::vector<int> fromLevels;
60  };
61 
62  // Utilities to get detector levels during a step
63  int detLevels(const G4VTouchable*) const;
64  G4LogicalVolume* detLV(const G4VTouchable*, int) const;
65  void detectorLevel(const G4VTouchable*, int&, int*,
66  G4String*) const;
67 
69  double eMin;
71  std::vector<Detector> detectors;
73 };
74 
75 #endif
76 
77 
78 
79 
80 
81 
82 
83 
void update(const BeginOfEvent *evt) override
This routine will be called when the appropriate signal arrives.
void EndOfEvent(G4HCofThisEvent *) override
G4LogicalVolume * detLV(const G4VTouchable *, int) const
std::vector< PCaloHit > PCaloHitContainer
CaloTrkProcessing(const std::string &aSDname, const DDCompactView &cpv, const SensitiveDetectorCatalog &clg, edm::ParameterSet const &p, const SimTrackManager *)
int detLevels(const G4VTouchable *) const
int isItCalo(const G4VTouchable *)
Compact representation of the geometrical detector hierarchy.
Definition: DDCompactView.h:83
void detectorLevel(const G4VTouchable *, int &, int *, G4String *) const
int isItInside(const G4VTouchable *, int, int)
bool ProcessHits(G4Step *, G4TouchableHistory *) override
const std::vector< std::string > & getNames() const
std::vector< int > fromLevels
std::vector< Detector > detectors
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:20
std::vector< double > getNumbers(G4String, const DDsvalues_type &)
std::vector< G4LogicalVolume * > fromDetL
void fillHits(edm::PCaloHitContainer &, const std::string &) override
std::vector< std::string > fromDets
uint32_t setDetUnitId(const G4Step *step) override
~CaloTrkProcessing() override
void Initialize(G4HCofThisEvent *) override
void clearHits() override
const SimTrackManager * m_trackManager
step