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 
12 
13 #include "G4VTouchable.hh"
14 #include "G4Track.hh"
16 
17 #include <map>
18 #include <vector>
19 #include <string>
20 #include <iostream>
21 
22 class BeginOfEvent;
23 class G4LogicalVolume;
24 class G4Step;
25 class SimTrackManager;
26 
28  public Observer<const BeginOfEvent*>,
29  public Observer<const G4Step*> {
30 public:
31  CaloTrkProcessing(const std::string& aSDname,
32  const CaloSimulationParameters& csps,
33  const SensitiveDetectorCatalog& clg,
34  bool testBeam,
35  double eMin,
36  bool putHistory,
37  bool doFineCalo,
38  double eMinFine,
39  int addlevel,
40  const std::vector<std::string>& fineNames,
41  const std::vector<int>& fineLevels,
42  const std::vector<int>& useFines,
43  const SimTrackManager*);
44  ~CaloTrkProcessing() override;
45  void Initialize(G4HCofThisEvent*) override {}
46  void clearHits() override {}
47  bool ProcessHits(G4Step*, G4TouchableHistory*) override { return true; }
48  uint32_t setDetUnitId(const G4Step* step) override { return 0; }
49  void EndOfEvent(G4HCofThisEvent*) override {}
50  void fillHits(edm::PCaloHitContainer&, const std::string&) override {}
51 
52 private:
53  struct Detector {
54  Detector() {}
56  G4LogicalVolume* lv;
57  int level;
58  std::vector<std::string> fromDets;
59  std::vector<G4LogicalVolume*> fromDetL;
60  std::vector<int> fromLevels;
61  };
62 
63  void update(const BeginOfEvent* evt) override;
64  void update(const G4Step*) override;
65  int isItCalo(const G4VTouchable*, const std::vector<Detector>&);
66  int isItInside(const G4VTouchable*, int, int);
67 
68  // Utilities to get detector levels during a step
69  int detLevels(const G4VTouchable*) const;
70  G4LogicalVolume* detLV(const G4VTouchable*, int) const;
71  void detectorLevel(const G4VTouchable*, int&, int*, G4String*) const;
72 
73  const bool testBeam_;
74  const double eMin_;
75  const bool putHistory_;
77  const double eMinFine_;
78  const int addlevel_;
80  std::vector<Detector> detectors_, fineDetectors_;
81 };
82 
83 #endif
SimTrackManager
Definition: SimTrackManager.h:35
CaloTrkProcessing::Detector::Detector
Detector()
Definition: CaloTrkProcessing.h:54
CaloSimulationParameters
Definition: CaloSimulationParameters.h:6
Observer
Definition: Observer.h:23
MessageLogger.h
SensitiveCaloDetector.h
step
step
Definition: StallMonitor.cc:94
CaloTrkProcessing::Detector::fromLevels
std::vector< int > fromLevels
Definition: CaloTrkProcessing.h:60
CaloTrkProcessing::detectorLevel
void detectorLevel(const G4VTouchable *, int &, int *, G4String *) const
Definition: CaloTrkProcessing.cc:404
CaloTrkProcessing::Detector::fromDets
std::vector< std::string > fromDets
Definition: CaloTrkProcessing.h:58
CaloTrkProcessing::Detector::level
int level
Definition: CaloTrkProcessing.h:57
CaloTrkProcessing::lastTrackID_
int lastTrackID_
Definition: CaloTrkProcessing.h:79
Observer.h
CaloTrkProcessing::doFineCalo_
bool doFineCalo_
Definition: CaloTrkProcessing.h:76
CaloTrkProcessing::isItInside
int isItInside(const G4VTouchable *, int, int)
Definition: CaloTrkProcessing.cc:322
SensitiveCaloDetector
Definition: SensitiveCaloDetector.h:10
CaloSimulationParameters.h
CaloTrkProcessing::detLevels
int detLevels(const G4VTouchable *) const
Definition: CaloTrkProcessing.cc:383
CaloTrkProcessing::detLV
G4LogicalVolume * detLV(const G4VTouchable *, int) const
Definition: CaloTrkProcessing.cc:391
photonAnalyzer_cfi.eMin
eMin
Definition: photonAnalyzer_cfi.py:50
CaloTrkProcessing::Detector::lv
G4LogicalVolume * lv
Definition: CaloTrkProcessing.h:56
CaloTrkProcessing::update
void update(const BeginOfEvent *evt) override
This routine will be called when the appropriate signal arrives.
Definition: CaloTrkProcessing.cc:170
CaloTrkProcessing::Detector::fromDetL
std::vector< G4LogicalVolume * > fromDetL
Definition: CaloTrkProcessing.h:59
SensitiveDetectorCatalog
Definition: SensitiveDetectorCatalog.h:10
CaloTrkProcessing::fineDetectors_
std::vector< Detector > fineDetectors_
Definition: CaloTrkProcessing.h:80
CaloTrkProcessing::fillHits
void fillHits(edm::PCaloHitContainer &, const std::string &) override
Definition: CaloTrkProcessing.h:50
CaloTrkProcessing
Definition: CaloTrkProcessing.h:27
CaloTrkProcessing::EndOfEvent
void EndOfEvent(G4HCofThisEvent *) override
Definition: CaloTrkProcessing.h:49
CaloTrkProcessing::eMinFine_
const double eMinFine_
Definition: CaloTrkProcessing.h:77
CaloTrkProcessing::addlevel_
const int addlevel_
Definition: CaloTrkProcessing.h:78
LorentzVector.h
CaloTrkProcessing::setDetUnitId
uint32_t setDetUnitId(const G4Step *step) override
Definition: CaloTrkProcessing.h:48
BeginOfEvent
Definition: BeginOfEvent.h:6
CaloTrkProcessing::Detector::name
std::string name
Definition: CaloTrkProcessing.h:55
CaloTrkProcessing::testBeam_
const bool testBeam_
Definition: CaloTrkProcessing.h:73
AlCaHLTBitMon_QueryRunRegistry.string
string string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
CaloTrkProcessing::Detector
Definition: CaloTrkProcessing.h:53
CaloTrkProcessing::ProcessHits
bool ProcessHits(G4Step *, G4TouchableHistory *) override
Definition: CaloTrkProcessing.h:47
Frameworkfwd.h
CaloTrkProcessing::Initialize
void Initialize(G4HCofThisEvent *) override
Definition: CaloTrkProcessing.h:45
CaloTrkProcessing::clearHits
void clearHits() override
Definition: CaloTrkProcessing.h:46
edm::PCaloHitContainer
std::vector< PCaloHit > PCaloHitContainer
Definition: PCaloHitContainer.h:8
ParameterSetfwd.h
CaloTrkProcessing::isItCalo
int isItCalo(const G4VTouchable *, const std::vector< Detector > &)
Definition: CaloTrkProcessing.cc:293
CaloTrkProcessing::detectors_
std::vector< Detector > detectors_
Definition: CaloTrkProcessing.h:80
CaloTrkProcessing::putHistory_
const bool putHistory_
Definition: CaloTrkProcessing.h:75
CaloTrkProcessing::~CaloTrkProcessing
~CaloTrkProcessing() override
Definition: CaloTrkProcessing.cc:168
CaloTrkProcessing::CaloTrkProcessing
CaloTrkProcessing(const std::string &aSDname, const CaloSimulationParameters &csps, const SensitiveDetectorCatalog &clg, bool testBeam, double eMin, bool putHistory, bool doFineCalo, double eMinFine, int addlevel, const std::vector< std::string > &fineNames, const std::vector< int > &fineLevels, const std::vector< int > &useFines, const SimTrackManager *)
Definition: CaloTrkProcessing.cc:21
CaloTrkProcessing::eMin_
const double eMin_
Definition: CaloTrkProcessing.h:74