CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
DTCalibValidation.h
Go to the documentation of this file.
1 
2 #ifndef DTCalibValidation_H
3 #define DTCalibValidation_H
4 
17 
20 
26 
28 
29 #include <string>
30 #include <map>
31 #include <vector>
32 
33 // To remove into CMSSW versions before 20X
34 // To add into CMSSW versions before 20X
35 //class DaqMonitorBEInterface;
36 
37 class DTGeometry;
38 class DTChamber;
39 
40 // FR class DTCalibValidation: public edm::EDAnalyzer{
42 public:
45 
47  ~DTCalibValidation() override;
48 
50  void dqmBeginRun(const edm::Run& r, const edm::EventSetup& c) override;
51 
52  // Operations
53  void analyze(const edm::Event& event, const edm::EventSetup& setup) override;
54 
55 protected:
56  void bookHistograms(DQMStore::IBooker&, edm::Run const&, edm::EventSetup const&) override;
57 
58 private:
59  // Switch for verbosity
60  //bool debug;
64  int nevent;
65  // the analysis type
67  // the geometry
70 
71  // Lable of 1D rechits in the event
73  // Lable of 2D segments in the event
75  // Lable of 4D segments in the event
77 
78  // Return a map between DTRecHit1DPair and wireId
79  std::map<DTWireId, std::vector<DTRecHit1DPair> > map1DRecHitsPerWire(const DTRecHitCollection* dt1DRecHitPairs);
80 
81  // Return a map between DTRecHit1D and wireId
82  std::map<DTWireId, std::vector<DTRecHit1D> > map1DRecHitsPerWire(const DTRecSegment2DCollection* segment2Ds);
83 
84  // Return a map between DTRecHit1D and wireId
85  std::map<DTWireId, std::vector<DTRecHit1D> > map1DRecHitsPerWire(const DTRecSegment4DCollection* segment4Ds);
86 
87  template <typename type>
88  const type* findBestRecHit(const DTLayer* layer,
89  DTWireId wireId,
90  const std::vector<type>& recHits,
91  const float simHitDist);
92 
93  // Compute the distance from wire (cm) of a hits in a DTRecHit1DPair
94  float recHitDistFromWire(const DTRecHit1DPair& hitPair, const DTLayer* layer);
95  // Compute the distance from wire (cm) of a hits in a DTRecHit1D
96  float recHitDistFromWire(const DTRecHit1D& recHit, const DTLayer* layer);
97  // Compute the position with respect to the wire (cm) of a hits in a DTRecHit1DPair
98  float recHitPosition(
99  const DTRecHit1DPair& hitPair, const DTLayer* layer, const DTChamber* chamber, float segmPos, int sl);
100  // Compute the position with respect to the wire (cm) of a hits in a DTRecHit1D
101  float recHitPosition(const DTRecHit1D& recHit, const DTLayer* layer, const DTChamber* chamber, float segmPos, int sl);
102 
103  // Does the real job
104  template <typename type>
105  void compute(const DTGeometry* dtGeom,
106  const DTRecSegment4D& segment,
107  const std::map<DTWireId, std::vector<type> >& recHitsPerWire,
108  int step);
109 
110  // Book a set of histograms for a give chamber
111  void bookHistos(DTSuperLayerId slId, int step);
112  // Fill a set of histograms for a give chamber
113  void fillHistos(DTSuperLayerId slId,
114  float distance,
115  float residualOnDistance,
116  float position,
117  float residualOnPosition,
118  int step);
119 
120  std::map<std::pair<DTSuperLayerId, int>, std::vector<MonitorElement*> > histosPerSL;
121 };
122 #endif
123 
124 /* Local Variables: */
125 /* show-trailing-whitespace: t */
126 /* truncate-lines: t */
127 /* End: */
const edm::EventSetup & c
const DTGeometry * dtGeom
float recHitDistFromWire(const DTRecHit1DPair &hitPair, const DTLayer *layer)
constexpr std::array< uint8_t, layerIndexSize > layer
void compute(const DTGeometry *dtGeom, const DTRecSegment4D &segment, const std::map< DTWireId, std::vector< type > > &recHitsPerWire, int step)
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
edm::EDGetTokenT< DTRecSegment4DCollection > segment4DToken_
std::map< DTWireId, std::vector< DTRecHit1DPair > > map1DRecHitsPerWire(const DTRecHitCollection *dt1DRecHitPairs)
void fillHistos(DTSuperLayerId slId, float distance, float residualOnDistance, float position, float residualOnPosition, int step)
~DTCalibValidation() override
Destructor.
void analyze(const edm::Event &event, const edm::EventSetup &setup) override
edm::EDGetTokenT< DTRecHitCollection > recHits1DToken_
edm::ESGetToken< DTGeometry, MuonGeometryRecord > muonGeomToken_
void bookHistos(DTSuperLayerId slId, int step)
void dqmBeginRun(const edm::Run &r, const edm::EventSetup &c) override
BeginRun.
static int position[264][3]
Definition: ReadPGInfo.cc:289
std::map< std::pair< DTSuperLayerId, int >, std::vector< MonitorElement * > > histosPerSL
DTCalibValidation(const edm::ParameterSet &pset)
Constructor.
edm::ParameterSet parameters
step
Definition: StallMonitor.cc:94
float recHitPosition(const DTRecHit1DPair &hitPair, const DTLayer *layer, const DTChamber *chamber, float segmPos, int sl)
Definition: Run.h:45
const type * findBestRecHit(const DTLayer *layer, DTWireId wireId, const std::vector< type > &recHits, const float simHitDist)
edm::EDGetTokenT< DTRecSegment2DCollection > segment2DToken_