Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:07:06

0001 
0002 #ifndef DTCalibValidation_H
0003 #define DTCalibValidation_H
0004 
0005 /** \class DTCalibValidation
0006  *  Analysis on DT residuals to validate the kFactor
0007  *
0008  *
0009  *  \author G. Mila - INFN Torino
0010  */
0011 
0012 #include "FWCore/Framework/interface/Frameworkfwd.h"
0013 #include "DataFormats/MuonDetId/interface/DTSuperLayerId.h"
0014 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0015 #include "FWCore/Framework/interface/LuminosityBlock.h"
0016 
0017 #include "DQMServices/Core/interface/DQMStore.h"
0018 #include "DQMServices/Core/interface/DQMEDAnalyzer.h"
0019 
0020 #include "DataFormats/MuonDetId/interface/DTWireId.h"
0021 #include "DataFormats/DTRecHit/interface/DTRecHitCollection.h"
0022 #include "DataFormats/DTRecHit/interface/DTRecSegment2DCollection.h"
0023 #include "DataFormats/DTRecHit/interface/DTRecSegment4DCollection.h"
0024 #include "Geometry/Records/interface/MuonGeometryRecord.h"
0025 
0026 #include "FWCore/Framework/interface/ESHandle.h"
0027 
0028 #include <string>
0029 #include <map>
0030 #include <vector>
0031 
0032 // To remove into CMSSW versions before 20X
0033 // To add into CMSSW versions before 20X
0034 //class DaqMonitorBEInterface;
0035 
0036 class DTGeometry;
0037 class DTChamber;
0038 
0039 // FR class DTCalibValidation: public edm::EDAnalyzer{
0040 class DTCalibValidation : public DQMEDAnalyzer {
0041 public:
0042   /// Constructor
0043   DTCalibValidation(const edm::ParameterSet& pset);
0044 
0045   /// Destructor
0046   ~DTCalibValidation() override;
0047 
0048   /// BeginRun
0049   void dqmBeginRun(const edm::Run& r, const edm::EventSetup& c) override;
0050 
0051   // Operations
0052   void analyze(const edm::Event& event, const edm::EventSetup& setup) override;
0053 
0054 protected:
0055   void bookHistograms(DQMStore::IBooker&, edm::Run const&, edm::EventSetup const&) override;
0056 
0057 private:
0058   // Switch for verbosity
0059   //bool debug;
0060   edm::ParameterSet parameters;
0061   int wrongSegment;
0062   int rightSegment;
0063   int nevent;
0064   // the analysis type
0065   bool detailedAnalysis;
0066   // the geometry
0067   edm::ESGetToken<DTGeometry, MuonGeometryRecord> muonGeomToken_;
0068   const DTGeometry* dtGeom;
0069 
0070   // Lable of 1D rechits in the event
0071   edm::EDGetTokenT<DTRecHitCollection> recHits1DToken_;
0072   // Lable of 2D segments in the event
0073   edm::EDGetTokenT<DTRecSegment2DCollection> segment2DToken_;
0074   // Lable of 4D segments in the event
0075   edm::EDGetTokenT<DTRecSegment4DCollection> segment4DToken_;
0076 
0077   // Return a map between DTRecHit1DPair and wireId
0078   std::map<DTWireId, std::vector<DTRecHit1DPair> > map1DRecHitsPerWire(const DTRecHitCollection* dt1DRecHitPairs);
0079 
0080   // Return a map between DTRecHit1D and wireId
0081   std::map<DTWireId, std::vector<DTRecHit1D> > map1DRecHitsPerWire(const DTRecSegment2DCollection* segment2Ds);
0082 
0083   // Return a map between DTRecHit1D and wireId
0084   std::map<DTWireId, std::vector<DTRecHit1D> > map1DRecHitsPerWire(const DTRecSegment4DCollection* segment4Ds);
0085 
0086   template <typename type>
0087   const type* findBestRecHit(const DTLayer* layer,
0088                              DTWireId wireId,
0089                              const std::vector<type>& recHits,
0090                              const float simHitDist);
0091 
0092   // Compute the distance from wire (cm) of a hits in a DTRecHit1DPair
0093   float recHitDistFromWire(const DTRecHit1DPair& hitPair, const DTLayer* layer);
0094   // Compute the distance from wire (cm) of a hits in a DTRecHit1D
0095   float recHitDistFromWire(const DTRecHit1D& recHit, const DTLayer* layer);
0096   // Compute the position with respect to the wire (cm) of a hits in a DTRecHit1DPair
0097   float recHitPosition(
0098       const DTRecHit1DPair& hitPair, const DTLayer* layer, const DTChamber* chamber, float segmPos, int sl);
0099   // Compute the position with respect to the wire (cm) of a hits in a DTRecHit1D
0100   float recHitPosition(const DTRecHit1D& recHit, const DTLayer* layer, const DTChamber* chamber, float segmPos, int sl);
0101 
0102   // Does the real job
0103   template <typename type>
0104   void compute(const DTGeometry* dtGeom,
0105                const DTRecSegment4D& segment,
0106                const std::map<DTWireId, std::vector<type> >& recHitsPerWire,
0107                int step);
0108 
0109   // Book a set of histograms for a give chamber
0110   void bookHistos(DTSuperLayerId slId, int step);
0111   // Fill a set of histograms for a give chamber
0112   void fillHistos(DTSuperLayerId slId,
0113                   float distance,
0114                   float residualOnDistance,
0115                   float position,
0116                   float residualOnPosition,
0117                   int step);
0118 
0119   std::map<std::pair<DTSuperLayerId, int>, std::vector<MonitorElement*> > histosPerSL;
0120 };
0121 #endif
0122 
0123 /* Local Variables: */
0124 /* show-trailing-whitespace: t */
0125 /* truncate-lines: t */
0126 /* End: */