Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #ifndef DTCHAMBEREFFICIENCY_H
0002 #define DTCHAMBEREFFICIENCY_H
0003 
0004 /** \class DTChamberEfficiency
0005  *
0006  * Description:
0007  *
0008  * This class provides the histograms for the calculation of the
0009  * efficiency of muons reconstruction in the DTs. It is applicable
0010  * both in presence or absence of a magnetic field.
0011  * Histos are 2D Sector vs Chamber plots for each wheel
0012  *
0013  * \author : Mario Pelliccioni - INFN Torino <pellicci@cern.ch>
0014  * $date   : 05/12/2008 16:51:04 CET $
0015  *
0016  * Modification:
0017  *
0018  */
0019 
0020 #include "DataFormats/TrackReco/interface/TrackFwd.h"
0021 #include "DataFormats/TrackReco/interface/Track.h"
0022 
0023 #include "FWCore/Framework/interface/ESHandle.h"
0024 #include "FWCore/Utilities/interface/InputTag.h"
0025 
0026 #include "FWCore/Framework/interface/FrameworkfwdMostUsed.h"
0027 
0028 #include "DQMServices/Core/interface/DQMStore.h"
0029 
0030 #include "DQMServices/Core/interface/DQMEDAnalyzer.h"
0031 
0032 #include "RecoMuon/MeasurementDet/interface/MuonDetLayerMeasurements.h"
0033 
0034 #include <string>
0035 #include <vector>
0036 
0037 namespace reco {
0038   class TransientTrack;
0039 }
0040 
0041 class Chi2MeasurementEstimator;
0042 class MuonServiceProxy;
0043 
0044 class FreeTrajectoryState;
0045 class DetLayer;
0046 class DetId;
0047 class NavigationSchool;
0048 
0049 class DTChamberEfficiency : public DQMEDAnalyzer {
0050 public:
0051   //Constructor
0052   DTChamberEfficiency(const edm::ParameterSet& pset);
0053 
0054   //Destructor
0055   ~DTChamberEfficiency() override;
0056 
0057   //Operations
0058   void analyze(const edm::Event& event, const edm::EventSetup& eventSetup) override;
0059   void dqmBeginRun(const edm::Run&, const edm::EventSetup&) override;
0060 
0061 protected:
0062   // Book the histograms
0063   void bookHistograms(DQMStore::IBooker&, edm::Run const&, edm::EventSetup const&) override;
0064 
0065 private:
0066   //functions
0067   std::vector<const DetLayer*> compatibleLayers(const NavigationSchool& navigationSchool,
0068                                                 const DetLayer* initialLayer,
0069                                                 const FreeTrajectoryState& fts,
0070                                                 PropagationDirection propDir);
0071 
0072   MeasurementContainer segQualityCut(const MeasurementContainer& seg_list) const;
0073   bool chamberSelection(const DetId& idDetLay, reco::TransientTrack& trans_track) const;
0074   inline edm::ESHandle<Propagator> propagator() const;
0075 
0076   //data members
0077   bool debug;
0078 
0079   edm::InputTag theTracksLabel_;
0080   edm::EDGetTokenT<reco::TrackCollection> theTracksToken_;
0081 
0082   edm::InputTag labelRPCRecHits;
0083   edm::InputTag thedt4DSegments;
0084   edm::InputTag thecscSegments;
0085 
0086   double theMaxChi2;
0087   double theNSigma;
0088   int theMinNrec;
0089 
0090   std::string theNavigationType;
0091 
0092   MuonServiceProxy* theService;
0093   MuonDetLayerMeasurements* theMeasurementExtractor;
0094   Chi2MeasurementEstimator* theEstimator;
0095 
0096   std::vector<std::vector<MonitorElement*> > histosPerW;
0097 };
0098 
0099 #endif  // DTANALYZER_H
0100 
0101 /* Local Variables: */
0102 /* show-trailing-whitespace: t */
0103 /* truncate-lines: t */
0104 /* End: */