CMS 3D CMS Logo

CaloTrkProcessing.h
Go to the documentation of this file.
1 #ifndef SimG4CMS_CaloTrkProcessing_H
2 #define SimG4CMS_CaloTrkProcessing_H
3 
6 
8 
9 #include "G4VTouchable.hh"
10 
11 #include <map>
12 #include <vector>
13 #include <string>
14 #include <iostream>
15 
16 class SimTrackManager;
17 class BeginOfEvent;
18 class G4LogicalVolume;
19 class G4Step;
20 class SimTrackManager;
21 
23  public Observer<const BeginOfEvent*>,
24  public Observer<const G4Step*> {
25 public:
26  CaloTrkProcessing(const std::string& aSDname,
27  const edm::EventSetup& es,
28  const SensitiveDetectorCatalog& clg,
29  edm::ParameterSet const& p,
30  const SimTrackManager*);
31  ~CaloTrkProcessing() override;
32  void Initialize(G4HCofThisEvent*) override {}
33  void clearHits() override {}
34  bool ProcessHits(G4Step*, G4TouchableHistory*) override { return true; }
35  uint32_t setDetUnitId(const G4Step* step) override { return 0; }
36  void EndOfEvent(G4HCofThisEvent*) override {}
37  void fillHits(edm::PCaloHitContainer&, const std::string&) override {}
38 
39 private:
40  struct Detector {
41  Detector() {}
43  G4LogicalVolume* lv;
44  int level;
45  std::vector<std::string> fromDets;
46  std::vector<G4LogicalVolume*> fromDetL;
47  std::vector<int> fromLevels;
48  };
49 
50  void update(const BeginOfEvent* evt) override;
51  void update(const G4Step*) override;
52  int isItCalo(const G4VTouchable*, const std::vector<Detector>&);
53  int isItInside(const G4VTouchable*, int, int);
54 
55  // Utilities to get detector levels during a step
56  int detLevels(const G4VTouchable*) const;
57  G4LogicalVolume* detLV(const G4VTouchable*, int) const;
58  void detectorLevel(const G4VTouchable*, int&, int*, G4String*) const;
59 
63  std::vector<Detector> detectors_, fineDetectors_;
64 };
65 
66 #endif
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
int detLevels(const G4VTouchable *) const
void detectorLevel(const G4VTouchable *, int &, int *, G4String *) const
int isItInside(const G4VTouchable *, int, int)
CaloTrkProcessing(const std::string &aSDname, const edm::EventSetup &es, const SensitiveDetectorCatalog &clg, edm::ParameterSet const &p, const SimTrackManager *)
bool ProcessHits(G4Step *, G4TouchableHistory *) override
std::vector< int > fromLevels
std::vector< Detector > fineDetectors_
std::vector< G4LogicalVolume * > fromDetL
void fillHits(edm::PCaloHitContainer &, const std::string &) override
std::vector< std::string > fromDets
uint32_t setDetUnitId(const G4Step *step) override
int isItCalo(const G4VTouchable *, const std::vector< Detector > &)
~CaloTrkProcessing() override
void Initialize(G4HCofThisEvent *) override
void clearHits() override
step
Definition: StallMonitor.cc:94
std::vector< Detector > detectors_