Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:33:46

0001 #include "Validation/TrackerRecHits/interface/SiStripRecHitsValid.h"
0002 
0003 //needed for the geometry:
0004 #include "DataFormats/DetId/interface/DetId.h"
0005 #include "DataFormats/SiStripDetId/interface/StripSubdetector.h"
0006 #include "DataFormats/TrackerCommon/interface/SiStripSubStructure.h"
0007 #include "DQM/SiStripCommon/interface/SiStripFolderOrganizer.h"
0008 #include "DataFormats/GeometryVector/interface/LocalPoint.h"
0009 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
0010 #include "DQM/SiStripCommon/interface/SiStripHistoId.h"
0011 
0012 //--- for RecHit
0013 #include "DataFormats/SiStripCluster/interface/SiStripCluster.h"
0014 #include "DataFormats/SiStripCluster/interface/SiStripClusterCollection.h"
0015 #include "DataFormats/Common/interface/OwnVector.h"
0016 #include "DQMServices/Core/interface/DQMStore.h"
0017 
0018 #include "CommonTools/TriggerUtils/interface/GenericTriggerEventFlag.h"
0019 
0020 using namespace std;
0021 using namespace edm;
0022 
0023 //Constructor
0024 SiStripRecHitsValid::SiStripRecHitsValid(const ParameterSet& ps)
0025     : m_geomToken(esConsumes()),
0026       m_topoToken(esConsumes()),
0027       m_topoTokenBR(esConsumes<edm::Transition::BeginRun>()),
0028       m_SiStripDetCablingToken(esConsumes<edm::Transition::BeginRun>()),
0029       conf_(ps),
0030       trackerHitAssociatorConfig_(ps, consumesCollector())
0031 // matchedRecHits_( ps.getParameter<edm::InputTag>("matchedRecHits") ),
0032 // rphiRecHits_( ps.getParameter<edm::InputTag>("rphiRecHits") ),
0033 // stereoRecHits_( ps.getParameter<edm::InputTag>("stereoRecHits") )
0034 {
0035   matchedRecHitsToken_ = consumes<SiStripMatchedRecHit2DCollection>(ps.getParameter<edm::InputTag>("matchedRecHits"));
0036 
0037   rphiRecHitsToken_ = consumes<SiStripRecHit2DCollection>(ps.getParameter<edm::InputTag>("rphiRecHits"));
0038 
0039   stereoRecHitsToken_ = consumes<SiStripRecHit2DCollection>(ps.getParameter<edm::InputTag>("stereoRecHits"));
0040 
0041   topFolderName_ = conf_.getParameter<std::string>("TopFolderName");
0042 
0043   SubDetList_ = conf_.getParameter<std::vector<std::string> >("SubDetList");
0044 
0045   edm::ParameterSet ParametersNumTotrphi = conf_.getParameter<edm::ParameterSet>("TH1NumTotrphi");
0046   switchNumTotrphi = ParametersNumTotrphi.getParameter<bool>("switchon");
0047 
0048   edm::ParameterSet ParametersNumTotStereo = conf_.getParameter<edm::ParameterSet>("TH1NumTotStereo");
0049   switchNumTotStereo = ParametersNumTotStereo.getParameter<bool>("switchon");
0050 
0051   edm::ParameterSet ParametersNumTotMatched = conf_.getParameter<edm::ParameterSet>("TH1NumTotMatched");
0052   switchNumTotMatched = ParametersNumTotMatched.getParameter<bool>("switchon");
0053 
0054   edm::ParameterSet ParametersNumrphi = conf_.getParameter<edm::ParameterSet>("TH1Numrphi");
0055   switchNumrphi = ParametersNumrphi.getParameter<bool>("switchon");
0056 
0057   edm::ParameterSet ParametersBunchrphi = conf_.getParameter<edm::ParameterSet>("TH1Bunchrphi");
0058   switchBunchrphi = ParametersBunchrphi.getParameter<bool>("switchon");
0059 
0060   edm::ParameterSet ParametersEventrphi = conf_.getParameter<edm::ParameterSet>("TH1Eventrphi");
0061   switchEventrphi = ParametersEventrphi.getParameter<bool>("switchon");
0062 
0063   edm::ParameterSet ParametersNumStereo = conf_.getParameter<edm::ParameterSet>("TH1NumStereo");
0064   switchNumStereo = ParametersNumStereo.getParameter<bool>("switchon");
0065 
0066   edm::ParameterSet ParametersBunchStereo = conf_.getParameter<edm::ParameterSet>("TH1BunchStereo");
0067   switchBunchStereo = ParametersBunchStereo.getParameter<bool>("switchon");
0068 
0069   edm::ParameterSet ParametersEventStereo = conf_.getParameter<edm::ParameterSet>("TH1EventStereo");
0070   switchEventStereo = ParametersEventStereo.getParameter<bool>("switchon");
0071 
0072   edm::ParameterSet ParametersNumMatched = conf_.getParameter<edm::ParameterSet>("TH1NumMatched");
0073   switchNumMatched = ParametersNumMatched.getParameter<bool>("switchon");
0074 
0075   edm::ParameterSet ParametersBunchMatched = conf_.getParameter<edm::ParameterSet>("TH1BunchMatched");
0076   switchBunchMatched = ParametersBunchMatched.getParameter<bool>("switchon");
0077 
0078   edm::ParameterSet ParametersEventMatched = conf_.getParameter<edm::ParameterSet>("TH1EventMatched");
0079   switchEventMatched = ParametersEventMatched.getParameter<bool>("switchon");
0080 
0081   edm::ParameterSet ParametersWclusrphi = conf_.getParameter<edm::ParameterSet>("TH1Wclusrphi");
0082   switchWclusrphi = ParametersWclusrphi.getParameter<bool>("switchon");
0083 
0084   edm::ParameterSet ParametersAdcrphi = conf_.getParameter<edm::ParameterSet>("TH1Adcrphi");
0085   switchAdcrphi = ParametersAdcrphi.getParameter<bool>("switchon");
0086 
0087   edm::ParameterSet ParametersPosxrphi = conf_.getParameter<edm::ParameterSet>("TH1Posxrphi");
0088   switchPosxrphi = ParametersPosxrphi.getParameter<bool>("switchon");
0089 
0090   edm::ParameterSet ParametersResolxrphi = conf_.getParameter<edm::ParameterSet>("TH1Resolxrphi");
0091   switchResolxrphi = ParametersResolxrphi.getParameter<bool>("switchon");
0092 
0093   edm::ParameterSet ParametersResrphi = conf_.getParameter<edm::ParameterSet>("TH1Resrphi");
0094   switchResrphi = ParametersResrphi.getParameter<bool>("switchon");
0095 
0096   edm::ParameterSet ParametersPullLFrphi = conf_.getParameter<edm::ParameterSet>("TH1PullLFrphi");
0097   switchPullLFrphi = ParametersPullLFrphi.getParameter<bool>("switchon");
0098 
0099   edm::ParameterSet ParametersPullMFrphi = conf_.getParameter<edm::ParameterSet>("TH1PullMFrphi");
0100   switchPullMFrphi = ParametersPullMFrphi.getParameter<bool>("switchon");
0101 
0102   edm::ParameterSet ParametersChi2rphi = conf_.getParameter<edm::ParameterSet>("TH1Chi2rphi");
0103   switchChi2rphi = ParametersChi2rphi.getParameter<bool>("switchon");
0104 
0105   edm::ParameterSet ParametersNsimHitrphi = conf_.getParameter<edm::ParameterSet>("TH1NsimHitrphi");
0106   switchNsimHitrphi = ParametersNsimHitrphi.getParameter<bool>("switchon");
0107 
0108   edm::ParameterSet ParametersWclusStereo = conf_.getParameter<edm::ParameterSet>("TH1WclusStereo");
0109   switchWclusStereo = ParametersWclusStereo.getParameter<bool>("switchon");
0110 
0111   edm::ParameterSet ParametersAdcStereo = conf_.getParameter<edm::ParameterSet>("TH1AdcStereo");
0112   switchAdcStereo = ParametersAdcStereo.getParameter<bool>("switchon");
0113 
0114   edm::ParameterSet ParametersPosxStereo = conf_.getParameter<edm::ParameterSet>("TH1PosxStereo");
0115   switchPosxStereo = ParametersPosxStereo.getParameter<bool>("switchon");
0116 
0117   edm::ParameterSet ParametersResolxStereo = conf_.getParameter<edm::ParameterSet>("TH1ResolxStereo");
0118   switchResolxStereo = ParametersResolxStereo.getParameter<bool>("switchon");
0119 
0120   edm::ParameterSet ParametersResStereo = conf_.getParameter<edm::ParameterSet>("TH1ResStereo");
0121   switchResStereo = ParametersResStereo.getParameter<bool>("switchon");
0122 
0123   edm::ParameterSet ParametersPullLFStereo = conf_.getParameter<edm::ParameterSet>("TH1PullLFStereo");
0124   switchPullLFStereo = ParametersPullLFStereo.getParameter<bool>("switchon");
0125 
0126   edm::ParameterSet ParametersPullMFStereo = conf_.getParameter<edm::ParameterSet>("TH1PullMFStereo");
0127   switchPullMFStereo = ParametersPullMFStereo.getParameter<bool>("switchon");
0128 
0129   edm::ParameterSet ParametersChi2Stereo = conf_.getParameter<edm::ParameterSet>("TH1Chi2Stereo");
0130   switchChi2Stereo = ParametersChi2Stereo.getParameter<bool>("switchon");
0131 
0132   edm::ParameterSet ParametersNsimHitStereo = conf_.getParameter<edm::ParameterSet>("TH1NsimHitStereo");
0133   switchNsimHitStereo = ParametersNsimHitStereo.getParameter<bool>("switchon");
0134 
0135   edm::ParameterSet ParametersPosxMatched = conf_.getParameter<edm::ParameterSet>("TH1PosxMatched");
0136   switchPosxMatched = ParametersPosxMatched.getParameter<bool>("switchon");
0137 
0138   edm::ParameterSet ParametersPosyMatched = conf_.getParameter<edm::ParameterSet>("TH1PosyMatched");
0139   switchPosyMatched = ParametersPosyMatched.getParameter<bool>("switchon");
0140 
0141   edm::ParameterSet ParametersResolxMatched = conf_.getParameter<edm::ParameterSet>("TH1ResolxMatched");
0142   switchResolxMatched = ParametersResolxMatched.getParameter<bool>("switchon");
0143 
0144   edm::ParameterSet ParametersResolyMatched = conf_.getParameter<edm::ParameterSet>("TH1ResolyMatched");
0145   switchResolyMatched = ParametersResolyMatched.getParameter<bool>("switchon");
0146 
0147   edm::ParameterSet ParametersResxMatched = conf_.getParameter<edm::ParameterSet>("TH1ResxMatched");
0148   switchResxMatched = ParametersResxMatched.getParameter<bool>("switchon");
0149 
0150   edm::ParameterSet ParametersResyMatched = conf_.getParameter<edm::ParameterSet>("TH1ResyMatched");
0151   switchResyMatched = ParametersResyMatched.getParameter<bool>("switchon");
0152 
0153   edm::ParameterSet ParametersChi2Matched = conf_.getParameter<edm::ParameterSet>("TH1Chi2Matched");
0154   switchChi2Matched = ParametersChi2Matched.getParameter<bool>("switchon");
0155 
0156   edm::ParameterSet ParametersNsimHitMatched = conf_.getParameter<edm::ParameterSet>("TH1NsimHitMatched");
0157   switchNsimHitMatched = ParametersNsimHitMatched.getParameter<bool>("switchon");
0158 }
0159 
0160 SiStripRecHitsValid::~SiStripRecHitsValid() {}
0161 
0162 //--------------------------------------------------------------------------------------------
0163 void SiStripRecHitsValid::bookHistograms(DQMStore::IBooker& ibooker, const edm::Run& run, const edm::EventSetup& es) {
0164   if (watchSiStripDetCablingRcd_.check(es)) {
0165     edm::LogInfo("SiStripRecHitsValid") << "SiStripRecHitsValid::beginRun: "
0166                                         << " Creating MEs for new Cabling ";
0167 
0168     createMEs(ibooker, es);
0169   }
0170 }
0171 
0172 void SiStripRecHitsValid::analyze(const edm::Event& e, const edm::EventSetup& es) {
0173   LogInfo("EventInfo") << " Run = " << e.id().run() << " Event = " << e.id().event();
0174 
0175   //Retrieve tracker topology from geometry
0176   const TrackerTopology* const tTopo = &es.getData(m_topoToken);
0177 
0178   // Step A: Get Inputs
0179   edm::Handle<SiStripMatchedRecHit2DCollection> rechitsmatched;
0180   edm::Handle<SiStripRecHit2DCollection> rechitsrphi;
0181   edm::Handle<SiStripRecHit2DCollection> rechitsstereo;
0182   e.getByToken(matchedRecHitsToken_, rechitsmatched);
0183   e.getByToken(rphiRecHitsToken_, rechitsrphi);
0184   e.getByToken(stereoRecHitsToken_, rechitsstereo);
0185 
0186   //Variables in order to count total num of rechitrphi,rechitstereo, rechitmatched in subdetectors
0187   std::map<std::string, int> totnumrechitrphi;
0188   std::map<std::string, int> totnumrechitstereo;
0189   std::map<std::string, int> totnumrechitmatched;
0190   int totrechitrphi = 0;
0191   int totrechitstereo = 0;
0192   int totrechitmatched = 0;
0193 
0194   TrackerHitAssociator associate(e, trackerHitAssociatorConfig_);
0195 
0196   const TrackerGeometry& tracker = es.getData(m_geomToken);
0197 
0198   SiStripHistoId hidmanager;
0199   SiStripFolderOrganizer fold_organ;
0200   for (auto const& theDetSet : *rechitsrphi) {
0201     DetId detid = theDetSet.detId();
0202     uint32_t myid = detid.rawId();
0203     totrechitrphi += theDetSet.size();
0204 
0205     std::string label = hidmanager.getSubdetid(myid, tTopo, true);
0206     std::map<std::string, LayerMEs>::iterator iLayerME = LayerMEsMap.find(label);
0207     std::pair<std::string, int32_t> det_lay_pair = fold_organ.GetSubDetAndLayer(myid, tTopo, true);
0208 
0209     totnumrechitrphi[det_lay_pair.first] = totnumrechitrphi[det_lay_pair.first] + theDetSet.size();
0210     //loop over rechits-rphi in the same subdetector
0211     if (iLayerME != LayerMEsMap.end()) {
0212       for (auto const& rechit : theDetSet) {
0213         const GeomDetUnit* det = tracker.idToDetUnit(detid);
0214         const StripTopology& topol = static_cast<const StripGeomDetUnit*>(det)->specificTopology();
0215         //analyze RecHits
0216         rechitanalysis(rechit, topol, associate);
0217         // fill the result in a histogram
0218         fillME(iLayerME->second.meWclusrphi, rechitpro.clusiz);
0219         fillME(iLayerME->second.meAdcrphi, rechitpro.cluchg);
0220         fillME(iLayerME->second.mePosxrphi, rechitpro.x);
0221         fillME(iLayerME->second.meResolxrphi, rechitpro.resolxx);
0222         fillME(iLayerME->second.meNsimHitrphi, rechitpro.NsimHit);
0223         if (rechitpro.NsimHit > 0) {
0224           std::map<std::string, SubDetMEs>::iterator iSubDetME = SubDetMEsMap.find(det_lay_pair.first);
0225           fillME(iSubDetME->second.meBunchrphi, rechitpro.bunch);
0226           if (rechitpro.bunch == 0)
0227             fillME(iSubDetME->second.meEventrphi, rechitpro.event);
0228         }
0229         if (rechitpro.resx != -999999. || rechitpro.pullMF != -999999. || rechitpro.chi2 != -999999.) {
0230           fillME(iLayerME->second.meResrphi, rechitpro.resx);
0231           fillME(iLayerME->second.mePullLFrphi, rechitpro.resx / sqrt(rechitpro.resolxx));
0232           fillME(iLayerME->second.mePullMFrphi, rechitpro.pullMF);
0233           fillME(iLayerME->second.meChi2rphi, rechitpro.chi2);
0234         }
0235       }
0236     }
0237   }
0238 
0239   // start loops over detectors with detected rechitsstereo
0240   for (auto const& theDetSet : *rechitsstereo) {
0241     DetId detid = theDetSet.detId();
0242     uint32_t myid = detid.rawId();
0243     totrechitstereo += theDetSet.size();
0244 
0245     std::string label = hidmanager.getSubdetid(myid, tTopo, true);
0246     std::map<std::string, StereoAndMatchedMEs>::iterator iStereoAndMatchedME = StereoAndMatchedMEsMap.find(label);
0247     std::pair<std::string, int32_t> det_lay_pair = fold_organ.GetSubDetAndLayer(myid, tTopo, true);
0248 
0249     totnumrechitstereo[det_lay_pair.first] = totnumrechitstereo[det_lay_pair.first] + theDetSet.size();
0250     //loop over rechits-stereo in the same subdetector
0251     if (iStereoAndMatchedME != StereoAndMatchedMEsMap.end()) {
0252       for (auto const& rechit : theDetSet) {
0253         const GeomDetUnit* det = tracker.idToDetUnit(detid);
0254         const StripTopology& topol = static_cast<const StripGeomDetUnit*>(det)->specificTopology();
0255         //analyze RecHits
0256         rechitanalysis(rechit, topol, associate);
0257         // fill the result in a histogram
0258         fillME(iStereoAndMatchedME->second.meWclusStereo, rechitpro.clusiz);
0259         fillME(iStereoAndMatchedME->second.meAdcStereo, rechitpro.cluchg);
0260         fillME(iStereoAndMatchedME->second.mePosxStereo, rechitpro.x);
0261         fillME(iStereoAndMatchedME->second.meResolxStereo, sqrt(rechitpro.resolxx));
0262         fillME(iStereoAndMatchedME->second.meNsimHitStereo, rechitpro.NsimHit);
0263         if (rechitpro.NsimHit > 0) {
0264           std::map<std::string, SubDetMEs>::iterator iSubDetME = SubDetMEsMap.find(det_lay_pair.first);
0265           fillME(iSubDetME->second.meBunchStereo, rechitpro.bunch);
0266           if (rechitpro.bunch == 0)
0267             fillME(iSubDetME->second.meEventStereo, rechitpro.event);
0268         }
0269         if (rechitpro.resx != -999999. || rechitpro.pullMF != -999999. || rechitpro.chi2 != -999999.) {
0270           fillME(iStereoAndMatchedME->second.meResStereo, rechitpro.resx);
0271           fillME(iStereoAndMatchedME->second.mePullLFStereo, rechitpro.resx / sqrt(rechitpro.resolxx));
0272           fillME(iStereoAndMatchedME->second.mePullMFStereo, rechitpro.pullMF);
0273           fillME(iStereoAndMatchedME->second.meChi2Stereo, rechitpro.chi2);
0274         }
0275       }
0276     }
0277   }
0278 
0279   // start loops over detectors with detected rechitmatched
0280   for (auto const& theDetSet : *rechitsmatched) {
0281     DetId detid = theDetSet.detId();
0282     uint32_t myid = detid.rawId();
0283     totrechitmatched += theDetSet.size();
0284 
0285     std::string label = hidmanager.getSubdetid(myid, tTopo, true);
0286     std::map<std::string, StereoAndMatchedMEs>::iterator iStereoAndMatchedME = StereoAndMatchedMEsMap.find(label);
0287     std::pair<std::string, int32_t> det_lay_pair = fold_organ.GetSubDetAndLayer(myid, tTopo, true);
0288 
0289     totnumrechitmatched[det_lay_pair.first] = totnumrechitmatched[det_lay_pair.first] + theDetSet.size();
0290     //loop over rechits-matched in the same subdetector
0291     if (iStereoAndMatchedME != StereoAndMatchedMEsMap.end()) {
0292       for (auto const& rechit : theDetSet) {
0293         const GluedGeomDet* gluedDet = static_cast<const GluedGeomDet*>(tracker.idToDet(rechit.geographicalId()));
0294         //analyze RecHits
0295         rechitanalysis_matched(rechit, gluedDet, associate);
0296         // fill the result in a histogram
0297         fillME(iStereoAndMatchedME->second.mePosxMatched, rechitpro.x);
0298         fillME(iStereoAndMatchedME->second.mePosyMatched, rechitpro.y);
0299         fillME(iStereoAndMatchedME->second.meResolxMatched, sqrt(rechitpro.resolxx));
0300         fillME(iStereoAndMatchedME->second.meResolyMatched, sqrt(rechitpro.resolyy));
0301         fillME(iStereoAndMatchedME->second.meNsimHitMatched, rechitpro.NsimHit);
0302         if (rechitpro.NsimHit > 0) {
0303           std::map<std::string, SubDetMEs>::iterator iSubDetME = SubDetMEsMap.find(det_lay_pair.first);
0304           fillME(iSubDetME->second.meBunchMatched, rechitpro.bunch);
0305           if (rechitpro.bunch == 0)
0306             fillME(iSubDetME->second.meEventMatched, rechitpro.event);
0307         }
0308         if (rechitpro.resx != -999999. || rechitpro.resy != -999999. || rechitpro.chi2 != -999999.) {
0309           fillME(iStereoAndMatchedME->second.meResxMatched, rechitpro.resx);
0310           fillME(iStereoAndMatchedME->second.meResyMatched, rechitpro.resy);
0311           fillME(iStereoAndMatchedME->second.meChi2Matched, rechitpro.chi2);
0312         }
0313       }
0314     }
0315   }  //End of loops over detectors
0316 
0317   //now fill the cumulative histograms of the hits
0318   for (std::vector<std::string>::iterator iSubdet = SubDetList_.begin(); iSubdet != SubDetList_.end(); ++iSubdet) {
0319     std::map<std::string, SubDetMEs>::iterator iSubDetME = SubDetMEsMap.find((*iSubdet));
0320     fillME(iSubDetME->second.meNumrphi, totnumrechitrphi[(*iSubdet)]);
0321     fillME(iSubDetME->second.meNumStereo, totnumrechitstereo[(*iSubdet)]);
0322     fillME(iSubDetME->second.meNumMatched, totnumrechitmatched[(*iSubdet)]);
0323   }
0324 
0325   fillME(totalMEs.meNumTotrphi, totrechitrphi);
0326   fillME(totalMEs.meNumTotStereo, totrechitstereo);
0327   fillME(totalMEs.meNumTotMatched, totrechitmatched);
0328 }
0329 
0330 //needed by to do the residual for matched hits
0331 std::pair<LocalPoint, LocalVector> SiStripRecHitsValid::projectHit(const PSimHit& hit,
0332                                                                    const StripGeomDetUnit* stripDet,
0333                                                                    const BoundPlane& plane) {
0334   //  const StripGeomDetUnit* stripDet = dynamic_cast<const StripGeomDetUnit*>(hit.det());
0335   //if (stripDet == nullptr) throw MeasurementDetException("HitMatcher hit is not on StripGeomDetUnit");
0336 
0337   const StripTopology& topol = stripDet->specificTopology();
0338   GlobalPoint globalpos = stripDet->surface().toGlobal(hit.localPosition());
0339   LocalPoint localHit = plane.toLocal(globalpos);
0340   //track direction
0341   LocalVector locdir = hit.localDirection();
0342   //rotate track in new frame
0343 
0344   GlobalVector globaldir = stripDet->surface().toGlobal(locdir);
0345   LocalVector dir = plane.toLocal(globaldir);
0346   double scale = -localHit.z() / dir.z();
0347 
0348   LocalPoint projectedPos = localHit + scale * dir;
0349 
0350   double selfAngle = topol.stripAngle(topol.strip(hit.localPosition()));
0351 
0352   LocalVector stripDir(sin(selfAngle), cos(selfAngle), 0);  // vector along strip in hit frame
0353 
0354   LocalVector localStripDir(plane.toLocal(stripDet->surface().toGlobal(stripDir)));
0355 
0356   return std::pair<LocalPoint, LocalVector>(projectedPos, localStripDir);
0357 }
0358 //--------------------------------------------------------------------------------------------
0359 void SiStripRecHitsValid::rechitanalysis(SiStripRecHit2D const rechit,
0360                                          const StripTopology& topol,
0361                                          TrackerHitAssociator& associate) {
0362   rechitpro.resx = -999999.;
0363   rechitpro.resy = -999999.;
0364   rechitpro.pullMF = -999999.;
0365   rechitpro.chi2 = -999999.;
0366   rechitpro.bunch = -999999.;
0367   rechitpro.event = -999999.;
0368 
0369   LocalPoint position = rechit.localPosition();
0370   LocalError error = rechit.localPositionError();
0371   MeasurementPoint Mposition = topol.measurementPosition(position);
0372   MeasurementError Merror = topol.measurementError(position, error);
0373   const auto& amplitudes = (rechit.cluster())->amplitudes();
0374   int totcharge = 0;
0375   for (auto ia : amplitudes) {
0376     totcharge += ia;
0377   }
0378   rechitpro.x = position.x();
0379   rechitpro.y = position.y();
0380   //rechitpro.z = position.z();
0381   rechitpro.resolxx = error.xx();
0382   //rechitpro.resolxy = error.xy();
0383   rechitpro.resolyy = error.yy();
0384   rechitpro.clusiz = amplitudes.size();
0385   rechitpro.cluchg = totcharge;
0386 
0387   auto const& matched = associate.associateHit(rechit);
0388   rechitpro.NsimHit = matched.size();
0389 
0390   if (!matched.empty()) {
0391     float mindist = std::numeric_limits<float>::max();
0392     float dist = std::numeric_limits<float>::max();
0393     PSimHit const* closest = nullptr;
0394 
0395     for (auto const& m : matched) {
0396       dist = fabs(rechitpro.x - m.localPosition().x());
0397       if (dist < mindist) {
0398         mindist = dist;
0399         closest = &m;
0400       }
0401     }
0402 
0403     if (!closest)
0404       return;
0405 
0406     rechitpro.bunch = closest->eventId().bunchCrossing();
0407     rechitpro.event = closest->eventId().event();
0408     rechitpro.resx = rechitpro.x - closest->localPosition().x();
0409     rechitpro.pullMF = (Mposition.x() - (topol.measurementPosition(closest->localPosition())).x()) / sqrt(Merror.uu());
0410 
0411     //chi2test compare rechit errors with the simhit position ( using null matrix for the simhit).
0412     //Can spot problems in the geometry better than a simple residual. (thanks to BorisM)
0413     AlgebraicVector rhparameters(2);  //= rechit.parameters();
0414     rhparameters[0] = position.x();
0415     rhparameters[1] = position.y();
0416     AlgebraicVector shparameters(2);
0417     shparameters[0] = closest->localPosition().x();
0418     shparameters[1] = closest->localPosition().y();
0419     AlgebraicVector r(rhparameters - shparameters);
0420     AlgebraicSymMatrix R(2);  //  = rechit.parametersError();
0421     R[0][0] = error.xx();
0422     R[0][1] = error.xy();
0423     R[1][1] = error.yy();
0424     int ierr;
0425     R.invert(ierr);  // if (ierr != 0) throw exception;
0426     float est = R.similarity(r);
0427     rechitpro.chi2 = est;
0428   }
0429 }
0430 
0431 //--------------------------------------------------------------------------------------------
0432 void SiStripRecHitsValid::rechitanalysis_matched(SiStripMatchedRecHit2D const rechit,
0433                                                  const GluedGeomDet* gluedDet,
0434                                                  TrackerHitAssociator& associate) {
0435   rechitpro.resx = -999999.;
0436   rechitpro.resy = -999999.;
0437   rechitpro.pullMF = -999999.;
0438   rechitpro.chi2 = -999999.;
0439   rechitpro.bunch = -999999.;
0440   rechitpro.event = -999999.;
0441   rechitpro.clusiz = -999999.;
0442   rechitpro.cluchg = -999999.;
0443 
0444   LocalPoint position = rechit.localPosition();
0445   LocalError error = rechit.localPositionError();
0446 
0447   rechitpro.x = position.x();
0448   rechitpro.y = position.y();
0449   //rechitpro.z = position.z();
0450   rechitpro.resolxx = error.xx();
0451   //rechitpro.resolxy = error.xy();
0452   rechitpro.resolyy = error.yy();
0453 
0454   auto const& matched = associate.associateHit(rechit);
0455   rechitpro.NsimHit = matched.size();
0456 
0457   if (matched.empty())
0458     return;
0459   float mindist = std::numeric_limits<float>::max();
0460   float dist = std::numeric_limits<float>::max();
0461   float dist2 = std::numeric_limits<float>::max();
0462   float distx = std::numeric_limits<float>::max();
0463   float disty = std::numeric_limits<float>::max();
0464   PSimHit const* closest = nullptr;
0465   std::pair<LocalPoint, LocalVector> closestPair;
0466 
0467   const StripGeomDetUnit* partnerstripdet = static_cast<const StripGeomDetUnit*>(gluedDet->stereoDet());
0468   std::pair<LocalPoint, LocalVector> hitPair;
0469 
0470   for (auto const& m : matched) {
0471     SiStripDetId hitDetId(m.detUnitId());
0472     if (hitDetId.stereo()) {  // project from the stereo sensor
0473                               //project simhit;
0474       hitPair = projectHit(m, partnerstripdet, gluedDet->surface());
0475       distx = rechitpro.x - hitPair.first.x();
0476       disty = rechitpro.y - hitPair.first.y();
0477       dist2 = distx * distx + disty * disty;
0478       dist = sqrt(dist2);
0479       if (dist < mindist) {
0480         mindist = dist;
0481         closestPair = hitPair;
0482         closest = &m;
0483       }
0484     }
0485   }
0486   if (!closest)
0487     return;
0488   rechitpro.bunch = closest->eventId().bunchCrossing();
0489   rechitpro.event = closest->eventId().event();
0490   rechitpro.resx = rechitpro.x - closestPair.first.x();
0491   rechitpro.resy = rechitpro.y - closestPair.first.y();
0492 
0493   //chi2test compare rechit errors with the simhit position ( using null matrix for the simhit).
0494   //Can spot problems in the geometry better than a simple residual. (thanks to BorisM)
0495   AlgebraicVector rhparameters(2);  //= rechit.parameters();
0496   rhparameters[0] = position.x();
0497   rhparameters[1] = position.y();
0498   LocalPoint sh = closestPair.first;
0499   AlgebraicVector shparameters(2);
0500   shparameters[0] = sh.x();
0501   shparameters[1] = sh.y();
0502   AlgebraicVector r(rhparameters - shparameters);
0503   AlgebraicSymMatrix R(2);  //  = rechit.parametersError();
0504   R[0][0] = error.xx();
0505   R[0][1] = error.xy();
0506   R[1][1] = error.yy();
0507   int ierr;
0508   R.invert(ierr);  // if (ierr != 0) throw exception;
0509   float est = R.similarity(r);
0510   rechitpro.chi2 = est;
0511 }
0512 
0513 //--------------------------------------------------------------------------------------------
0514 void SiStripRecHitsValid::createMEs(DQMStore::IBooker& ibooker, const edm::EventSetup& es) {
0515   //Retrieve tracker topology from geometry
0516   const TrackerTopology* const tTopo = &es.getData(m_topoTokenBR);
0517 
0518   // take from eventSetup the SiStripDetCabling object - here will use SiStripDetControl later on
0519   const auto& SiStripDetCabling_ = &es.getData(m_SiStripDetCablingToken);
0520 
0521   // get list of active detectors from SiStripDetCabling
0522   std::vector<uint32_t> activeDets;
0523   SiStripDetCabling_->addActiveDetectorsRawIds(activeDets);
0524 
0525   SiStripFolderOrganizer folder_organizer;
0526   // folder_organizer.setSiStripFolderName(topFolderName_);
0527   std::string curfold = topFolderName_;
0528   folder_organizer.setSiStripFolderName(curfold);
0529   folder_organizer.setSiStripFolder();
0530 
0531   createTotalMEs(ibooker);
0532   // loop over detectors and book MEs
0533   edm::LogInfo("SiStripTkRecHits|SiStripRecHitsValid") << "nr. of activeDets:  " << activeDets.size();
0534   const std::string &tec = "TEC", tid = "TID", tob = "TOB", tib = "TIB";
0535   for (auto detid_iterator = activeDets.begin(), detid_end = activeDets.end(); detid_iterator != detid_end;
0536        ++detid_iterator) {
0537     uint32_t detid = (*detid_iterator);
0538     // remove any eventual zero elements - there should be none, but just in case
0539     if (detid == 0) {
0540       activeDets.erase(detid_iterator);
0541       continue;
0542     }
0543 
0544     // Create Layer Level MEs
0545     std::pair<std::string, int32_t> det_layer_pair = folder_organizer.GetSubDetAndLayer(detid, tTopo, true);
0546     SiStripHistoId hidmanager;
0547     std::string label = hidmanager.getSubdetid(detid, tTopo, true);
0548 
0549     if (LayerMEsMap.find(label) == LayerMEsMap.end()) {
0550       // get detids for the layer
0551       // Keep in mind that when we are on the TID or TEC we deal with rings not wheel
0552       int32_t lnumber = det_layer_pair.second;
0553       const std::string& lname = det_layer_pair.first;
0554       std::vector<uint32_t> layerDetIds;
0555       if (lname == tec) {
0556         if (lnumber > 0) {
0557           SiStripSubStructure::getTECDetectors(activeDets, layerDetIds, tTopo, 2, 0, 0, 0, abs(lnumber), 0);
0558         } else if (lnumber < 0) {
0559           SiStripSubStructure::getTECDetectors(activeDets, layerDetIds, tTopo, 1, 0, 0, 0, abs(lnumber), 0);
0560         }
0561       } else if (lname == tid) {
0562         if (lnumber > 0) {
0563           SiStripSubStructure::getTIDDetectors(activeDets, layerDetIds, tTopo, 2, 0, abs(lnumber), 0);
0564         } else if (lnumber < 0) {
0565           SiStripSubStructure::getTIDDetectors(activeDets, layerDetIds, tTopo, 1, 0, abs(lnumber), 0);
0566         }
0567       } else if (lname == tob) {
0568         SiStripSubStructure::getTOBDetectors(activeDets, layerDetIds, tTopo, lnumber, 0, 0);
0569       } else if (lname == tib) {
0570         SiStripSubStructure::getTIBDetectors(activeDets, layerDetIds, tTopo, lnumber, 0, 0, 0);
0571       }
0572       LayerDetMap[label] = layerDetIds;
0573 
0574       // book Layer MEs
0575       folder_organizer.setLayerFolder(detid, tTopo, det_layer_pair.second, true);
0576       createLayerMEs(ibooker, label);
0577     }
0578     // book sub-detector plots
0579     if (SubDetMEsMap.find(det_layer_pair.first) == SubDetMEsMap.end()) {
0580       auto sdet_pair = folder_organizer.getSubDetFolderAndTag(detid, tTopo);
0581       ibooker.setCurrentFolder(sdet_pair.first);
0582       createSubDetMEs(ibooker, det_layer_pair.first);
0583     }
0584     //Create StereoAndMatchedMEs
0585     if (StereoAndMatchedMEsMap.find(label) == StereoAndMatchedMEsMap.end()) {
0586       // get detids for the stereo and matched layer. We are going to need a bool for these layers
0587       bool isStereo = false;
0588       // Keep in mind that when we are on the TID or TEC we deal with rings not wheel
0589       std::vector<uint32_t> stereoandmatchedDetIds;
0590       int32_t stereolnumber = det_layer_pair.second;
0591       const std::string& stereolname = det_layer_pair.first;
0592       if (stereolname == tec && (tTopo->tecIsStereo(detid))) {
0593         if (stereolnumber > 0) {
0594           SiStripSubStructure::getTECDetectors(
0595               activeDets, stereoandmatchedDetIds, tTopo, 2, 0, 0, 0, abs(stereolnumber), 1);
0596           isStereo = true;
0597         } else if (stereolnumber < 0) {
0598           SiStripSubStructure::getTECDetectors(
0599               activeDets, stereoandmatchedDetIds, tTopo, 1, 0, 0, 0, abs(stereolnumber), 1);
0600           isStereo = true;
0601         }
0602       } else if (stereolname == tid && (tTopo->tidIsStereo(detid))) {
0603         if (stereolnumber > 0) {
0604           SiStripSubStructure::getTIDDetectors(activeDets, stereoandmatchedDetIds, tTopo, 2, 0, abs(stereolnumber), 1);
0605           isStereo = true;
0606         } else if (stereolnumber < 0) {
0607           SiStripSubStructure::getTIDDetectors(activeDets, stereoandmatchedDetIds, tTopo, 1, 0, abs(stereolnumber), 1);
0608           isStereo = true;
0609         }
0610       } else if (stereolname == tob && (tTopo->tobIsStereo(detid))) {
0611         SiStripSubStructure::getTOBDetectors(activeDets, stereoandmatchedDetIds, tTopo, stereolnumber, 0, 0);
0612         isStereo = true;
0613       } else if (stereolname == tib && (tTopo->tibIsStereo(detid))) {
0614         SiStripSubStructure::getTIBDetectors(activeDets, stereoandmatchedDetIds, tTopo, stereolnumber, 0, 0, 0);
0615         isStereo = true;
0616       }
0617 
0618       StereoAndMatchedDetMap[label] = stereoandmatchedDetIds;
0619 
0620       if (isStereo) {
0621         //book StereoAndMatched MEs
0622         folder_organizer.setLayerFolder(detid, tTopo, det_layer_pair.second, true);
0623         //Create the Monitor Elements only when we have a stereo module
0624         createStereoAndMatchedMEs(ibooker, label);
0625       }
0626     }
0627   }  //end of loop over detectors
0628 }
0629 //------------------------------------------------------------------------------------------
0630 void SiStripRecHitsValid::createTotalMEs(DQMStore::IBooker& ibooker) {
0631   totalMEs.meNumTotrphi = nullptr;
0632   totalMEs.meNumTotStereo = nullptr;
0633   totalMEs.meNumTotMatched = nullptr;
0634 
0635   //NumTotrphi
0636   if (switchNumTotrphi) {
0637     totalMEs.meNumTotrphi = bookME1D(ibooker, "TH1NumTotrphi", "TH1NumTotrphi", "Num of RecHits rphi");
0638     totalMEs.meNumTotrphi->setAxisTitle("Total number of RecHits");
0639   }
0640   //NumTotStereo
0641   if (switchNumTotStereo) {
0642     totalMEs.meNumTotStereo = bookME1D(ibooker, "TH1NumTotStereo", "TH1NumTotStereo", "Num of RecHits stereo");
0643     totalMEs.meNumTotStereo->setAxisTitle("Total number of RecHits (stereo)");
0644   }
0645   //NumTotMatched
0646   if (switchNumTotMatched) {
0647     totalMEs.meNumTotMatched = bookME1D(ibooker, "TH1NumTotMatched", "TH1NumTotMatched", "Num of RecHits matched");
0648     totalMEs.meNumTotMatched->setAxisTitle("Total number of matched RecHits");
0649   }
0650 }
0651 //------------------------------------------------------------------------------------------
0652 void SiStripRecHitsValid::createLayerMEs(DQMStore::IBooker& ibooker, std::string label) {
0653   SiStripHistoId hidmanager;
0654   LayerMEs layerMEs;
0655 
0656   layerMEs.meWclusrphi = nullptr;
0657   layerMEs.meAdcrphi = nullptr;
0658   layerMEs.mePosxrphi = nullptr;
0659   layerMEs.meResolxrphi = nullptr;
0660   layerMEs.meResrphi = nullptr;
0661   layerMEs.mePullLFrphi = nullptr;
0662   layerMEs.mePullMFrphi = nullptr;
0663   layerMEs.meChi2rphi = nullptr;
0664   layerMEs.meNsimHitrphi = nullptr;
0665 
0666   //Wclusrphi
0667   if (switchWclusrphi) {
0668     layerMEs.meWclusrphi = bookME1D(ibooker,
0669                                     "TH1Wclusrphi",
0670                                     hidmanager.createHistoLayer("Wclus_rphi", "layer", label, "").c_str(),
0671                                     "Cluster Width - Number of strips that belong to the RecHit cluster");
0672     layerMEs.meWclusrphi->setAxisTitle("Cluster Width [nr strips] in " + label);
0673   }
0674   //Adcrphi
0675   if (switchAdcrphi) {
0676     layerMEs.meAdcrphi = bookME1D(ibooker,
0677                                   "TH1Adcrphi",
0678                                   hidmanager.createHistoLayer("Adc_rphi", "layer", label, "").c_str(),
0679                                   "RecHit Cluster Charge");
0680     layerMEs.meAdcrphi->setAxisTitle("cluster charge [ADC] in " + label);
0681   }
0682   //Posxrphi
0683   if (switchPosxrphi) {
0684     layerMEs.mePosxrphi = bookME1D(ibooker,
0685                                    "TH1Posxrphi",
0686                                    hidmanager.createHistoLayer("Posx_rphi", "layer", label, "").c_str(),
0687                                    "RecHit x coord.");
0688     layerMEs.mePosxrphi->setAxisTitle("x RecHit coord. (local frame) in " + label);
0689   }
0690   //Resolxrphi
0691   if (switchResolxrphi) {
0692     layerMEs.meResolxrphi = bookME1D(ibooker,
0693                                      "TH1Resolxrphi",
0694                                      hidmanager.createHistoLayer("Resolx_rphi", "layer", label, "").c_str(),
0695                                      "RecHit resol(x) coord.");  //<resolor>~20micron
0696     layerMEs.meResolxrphi->setAxisTitle("resol(x) RecHit coord. (local frame) in " + label);
0697   }
0698   //Resrphi
0699   if (switchResrphi) {
0700     layerMEs.meResrphi = bookME1D(ibooker,
0701                                   "TH1Resrphi",
0702                                   hidmanager.createHistoLayer("Res_rphi", "layer", label, "").c_str(),
0703                                   "Residuals of the hit x coordinate");
0704     layerMEs.meResrphi->setAxisTitle("RecHit Res(x) in " + label);
0705   }
0706   //PullLFrphi
0707   if (switchPullLFrphi) {
0708     layerMEs.mePullLFrphi = bookME1D(ibooker,
0709                                      "TH1PullLFrphi",
0710                                      hidmanager.createHistoLayer("Pull_LF_rphi", "layer", label, "").c_str(),
0711                                      "Pull distribution");
0712     layerMEs.mePullLFrphi->setAxisTitle("Pull distribution (local frame) in " + label);
0713   }
0714   //PullMFrphi
0715   if (switchPullMFrphi) {
0716     layerMEs.mePullMFrphi = bookME1D(ibooker,
0717                                      "TH1PullMFrphi",
0718                                      hidmanager.createHistoLayer("Pull_MF_rphi", "layer", label, "").c_str(),
0719                                      "Pull distribution");
0720     layerMEs.mePullMFrphi->setAxisTitle("Pull distribution (measurement frame) in " + label);
0721   }
0722   //Chi2rphi
0723   if (switchChi2rphi) {
0724     layerMEs.meChi2rphi = bookME1D(ibooker,
0725                                    "TH1Chi2rphi",
0726                                    hidmanager.createHistoLayer("Chi2_rphi", "layer", label, "").c_str(),
0727                                    "RecHit Chi2 test");
0728     layerMEs.meChi2rphi->setAxisTitle("RecHit Chi2 test in " + label);
0729   }
0730   //NsimHitrphi
0731   if (switchNsimHitrphi) {
0732     layerMEs.meNsimHitrphi = bookME1D(ibooker,
0733                                       "TH1NsimHitrphi",
0734                                       hidmanager.createHistoLayer("NsimHit_rphi", "layer", label, "").c_str(),
0735                                       "No. of assoc. simHits");
0736     layerMEs.meNsimHitrphi->setAxisTitle("Number of assoc. simHits in " + label);
0737   }
0738 
0739   LayerMEsMap[label] = layerMEs;
0740 }
0741 //------------------------------------------------------------------------------------------
0742 void SiStripRecHitsValid::createStereoAndMatchedMEs(DQMStore::IBooker& ibooker, std::string label) {
0743   SiStripHistoId hidmanager;
0744   StereoAndMatchedMEs stereoandmatchedMEs;
0745 
0746   stereoandmatchedMEs.meWclusStereo = nullptr;
0747   stereoandmatchedMEs.meAdcStereo = nullptr;
0748   stereoandmatchedMEs.mePosxStereo = nullptr;
0749   stereoandmatchedMEs.meResolxStereo = nullptr;
0750   stereoandmatchedMEs.meResStereo = nullptr;
0751   stereoandmatchedMEs.mePullLFStereo = nullptr;
0752   stereoandmatchedMEs.mePullMFStereo = nullptr;
0753   stereoandmatchedMEs.meChi2Stereo = nullptr;
0754   stereoandmatchedMEs.meNsimHitStereo = nullptr;
0755   stereoandmatchedMEs.mePosxMatched = nullptr;
0756   stereoandmatchedMEs.mePosyMatched = nullptr;
0757   stereoandmatchedMEs.meResolxMatched = nullptr;
0758   stereoandmatchedMEs.meResolyMatched = nullptr;
0759   stereoandmatchedMEs.meResxMatched = nullptr;
0760   stereoandmatchedMEs.meResyMatched = nullptr;
0761   stereoandmatchedMEs.meChi2Matched = nullptr;
0762   stereoandmatchedMEs.meNsimHitMatched = nullptr;
0763 
0764   //WclusStereo
0765   if (switchWclusStereo) {
0766     stereoandmatchedMEs.meWclusStereo =
0767         bookME1D(ibooker,
0768                  "TH1WclusStereo",
0769                  hidmanager.createHistoLayer("Wclus_stereo", "layer", label, "").c_str(),
0770                  "Cluster Width - Number of strips that belong to the RecHit cluster");
0771     stereoandmatchedMEs.meWclusStereo->setAxisTitle("Cluster Width [nr strips] in stereo modules in " + label);
0772   }
0773   //AdcStereo
0774   if (switchAdcStereo) {
0775     stereoandmatchedMEs.meAdcStereo = bookME1D(ibooker,
0776                                                "TH1AdcStereo",
0777                                                hidmanager.createHistoLayer("Adc_stereo", "layer", label, "").c_str(),
0778                                                "RecHit Cluster Charge");
0779     stereoandmatchedMEs.meAdcStereo->setAxisTitle("cluster charge [ADC] in stereo modules in " + label);
0780   }
0781   //PosxStereo
0782   if (switchPosxStereo) {
0783     stereoandmatchedMEs.mePosxStereo = bookME1D(ibooker,
0784                                                 "TH1PosxStereo",
0785                                                 hidmanager.createHistoLayer("Posx_stereo", "layer", label, "").c_str(),
0786                                                 "RecHit x coord.");
0787     stereoandmatchedMEs.mePosxStereo->setAxisTitle("x RecHit coord. (local frame) in stereo modules in " + label);
0788   }
0789   //ResolxStereo
0790   if (switchResolxStereo) {
0791     stereoandmatchedMEs.meResolxStereo =
0792         bookME1D(ibooker,
0793                  "TH1ResolxStereo",
0794                  hidmanager.createHistoLayer("Resolx_stereo", "layer", label, "").c_str(),
0795                  "RecHit resol(x) coord.");
0796     stereoandmatchedMEs.meResolxStereo->setAxisTitle("resol(x) RecHit coord. (local frame) in stereo modules in " +
0797                                                      label);
0798   }
0799   //ResStereo
0800   if (switchResStereo) {
0801     stereoandmatchedMEs.meResStereo = bookME1D(ibooker,
0802                                                "TH1ResStereo",
0803                                                hidmanager.createHistoLayer("Res_stereo", "layer", label, "").c_str(),
0804                                                "Residuals of the hit x coordinate");
0805     stereoandmatchedMEs.meResStereo->setAxisTitle("RecHit Res(x) in stereo modules in " + label);
0806   }
0807   //PullLFStereo
0808   if (switchPullLFStereo) {
0809     stereoandmatchedMEs.mePullLFStereo =
0810         bookME1D(ibooker,
0811                  "TH1PullLFStereo",
0812                  hidmanager.createHistoLayer("Pull_LF_stereo", "layer", label, "").c_str(),
0813                  "Pull distribution");
0814     stereoandmatchedMEs.mePullLFStereo->setAxisTitle("Pull distribution (local frame) in stereo modules in " + label);
0815   }
0816   //PullMFStereo
0817   if (switchPullMFStereo) {
0818     stereoandmatchedMEs.mePullMFStereo =
0819         bookME1D(ibooker,
0820                  "TH1PullMFStereo",
0821                  hidmanager.createHistoLayer("Pull_MF_stereo", "layer", label, "").c_str(),
0822                  "Pull distribution");
0823     stereoandmatchedMEs.mePullMFStereo->setAxisTitle("Pull distribution (measurement frame) in stereo modules in " +
0824                                                      label);
0825   }
0826   //Chi2Stereo
0827   if (switchChi2Stereo) {
0828     stereoandmatchedMEs.meChi2Stereo = bookME1D(ibooker,
0829                                                 "TH1Chi2Stereo",
0830                                                 hidmanager.createHistoLayer("Chi2_stereo", "layer", label, "").c_str(),
0831                                                 "RecHit Chi2 test");
0832     stereoandmatchedMEs.meChi2Stereo->setAxisTitle("RecHit Chi2 test in stereo modules in " + label);
0833   }
0834   //NsimHitStereo
0835   if (switchNsimHitStereo) {
0836     stereoandmatchedMEs.meNsimHitStereo =
0837         bookME1D(ibooker,
0838                  "TH1NsimHitStereo",
0839                  hidmanager.createHistoLayer("NsimHit_stereo", "layer", label, "").c_str(),
0840                  "No. of assoc. simHits");
0841     stereoandmatchedMEs.meNsimHitStereo->setAxisTitle("Number of assoc. simHits in stereo modules in " + label);
0842   }
0843   //PosxMatched
0844   if (switchPosxMatched) {
0845     stereoandmatchedMEs.mePosxMatched =
0846         bookME1D(ibooker,
0847                  "TH1PosxMatched",
0848                  hidmanager.createHistoLayer("Posx_matched", "layer", label, "").c_str(),
0849                  "RecHit x coord.");
0850     stereoandmatchedMEs.mePosxMatched->setAxisTitle("x coord. matched RecHit (local frame) in " + label);
0851   }
0852   //PosyMatched
0853   if (switchPosyMatched) {
0854     stereoandmatchedMEs.mePosyMatched =
0855         bookME1D(ibooker,
0856                  "TH1PosyMatched",
0857                  hidmanager.createHistoLayer("Posy_matched", "layer", label, "").c_str(),
0858                  "RecHit y coord.");
0859     stereoandmatchedMEs.mePosyMatched->setAxisTitle("y coord. matched RecHit (local frame) in " + label);
0860   }
0861   //ResolxMatched
0862   if (switchResolxMatched) {
0863     stereoandmatchedMEs.meResolxMatched =
0864         bookME1D(ibooker,
0865                  "TH1ResolxMatched",
0866                  hidmanager.createHistoLayer("Resolx_matched", "layer", label, "").c_str(),
0867                  "RecHit resol(x) coord.");
0868     stereoandmatchedMEs.meResolxMatched->setAxisTitle("resol(x) coord. matched RecHit (local frame) in " + label);
0869   }
0870   //ResolyMatched
0871   if (switchResolyMatched) {
0872     stereoandmatchedMEs.meResolyMatched =
0873         bookME1D(ibooker,
0874                  "TH1ResolyMatched",
0875                  hidmanager.createHistoLayer("Resoly_matched", "layer", label, "").c_str(),
0876                  "RecHit resol(y) coord.");
0877     stereoandmatchedMEs.meResolyMatched->setAxisTitle("resol(y) coord. matched RecHit (local frame) in " + label);
0878   }
0879   //ResxMatched
0880   if (switchResxMatched) {
0881     stereoandmatchedMEs.meResxMatched =
0882         bookME1D(ibooker,
0883                  "TH1ResxMatched",
0884                  hidmanager.createHistoLayer("Resx_matched", "layer", label, "").c_str(),
0885                  "Residuals of the hit x coord.");
0886     stereoandmatchedMEs.meResxMatched->setAxisTitle("Res(x) in matched RecHit in " + label);
0887   }
0888   //ResyMatched
0889   if (switchResyMatched) {
0890     stereoandmatchedMEs.meResyMatched =
0891         bookME1D(ibooker,
0892                  "TH1ResyMatched",
0893                  hidmanager.createHistoLayer("Resy_matched", "layer", label, "").c_str(),
0894                  "Residuals of the hit y coord.");
0895     stereoandmatchedMEs.meResyMatched->setAxisTitle("Res(y) in matched RecHit in " + label);
0896   }
0897   //Chi2Matched
0898   if (switchChi2Matched) {
0899     stereoandmatchedMEs.meChi2Matched =
0900         bookME1D(ibooker,
0901                  "TH1Chi2Matched",
0902                  hidmanager.createHistoLayer("Chi2_matched", "layer", label, "").c_str(),
0903                  "RecHit Chi2 test");
0904     stereoandmatchedMEs.meChi2Matched->setAxisTitle("Matched RecHit Chi2 test in " + label);
0905   }
0906   //NsimHitMatched
0907   if (switchNsimHitMatched) {
0908     stereoandmatchedMEs.meNsimHitMatched =
0909         bookME1D(ibooker,
0910                  "TH1NsimHitMatched",
0911                  hidmanager.createHistoLayer("NsimHit_matched", "layer", label, "").c_str(),
0912                  "No. of assoc. simHits");
0913     stereoandmatchedMEs.meNsimHitMatched->setAxisTitle("Number of assoc. simHits in " + label);
0914   }
0915 
0916   StereoAndMatchedMEsMap[label] = stereoandmatchedMEs;
0917 }
0918 //------------------------------------------------------------------------------------------
0919 void SiStripRecHitsValid::createSubDetMEs(DQMStore::IBooker& ibooker, std::string label) {
0920   SubDetMEs subdetMEs;
0921   subdetMEs.meNumrphi = nullptr;
0922   subdetMEs.meBunchrphi = nullptr;
0923   subdetMEs.meEventrphi = nullptr;
0924   subdetMEs.meNumStereo = nullptr;
0925   subdetMEs.meBunchStereo = nullptr;
0926   subdetMEs.meEventStereo = nullptr;
0927   subdetMEs.meNumMatched = nullptr;
0928   subdetMEs.meBunchMatched = nullptr;
0929   subdetMEs.meEventMatched = nullptr;
0930 
0931   std::string HistoName;
0932   //Numrphi
0933   if (switchNumrphi) {
0934     HistoName = "TH1Numrphi__" + label;
0935     subdetMEs.meNumrphi = bookME1D(ibooker, "TH1Numrphi", HistoName.c_str(), "Num of RecHits");
0936     subdetMEs.meNumrphi->setAxisTitle("Total number of RecHits in " + label);
0937   }
0938   //Bunchrphi
0939   if (switchBunchrphi) {
0940     HistoName = "TH1Bunchrphi__" + label;
0941     subdetMEs.meBunchrphi = bookME1D(ibooker, "TH1Bunchrphi", HistoName.c_str(), "Bunch Crossing");
0942     subdetMEs.meBunchrphi->setAxisTitle("Bunch crossing in " + label);
0943   }
0944   //Eventrphi
0945   if (switchEventrphi) {
0946     HistoName = "TH1Eventrphi__" + label;
0947     subdetMEs.meEventrphi = bookME1D(ibooker, "TH1Eventrphi", HistoName.c_str(), "Event (in-time bunch)");
0948     subdetMEs.meEventrphi->setAxisTitle("Event (in-time bunch) in " + label);
0949   }
0950   //NumStereo
0951   if (switchNumStereo) {
0952     HistoName = "TH1NumStereo__" + label;
0953     subdetMEs.meNumStereo = bookME1D(ibooker, "TH1NumStereo", HistoName.c_str(), "Num of RecHits in stereo modules");
0954     subdetMEs.meNumStereo->setAxisTitle("Total number of RecHits, stereo modules in " + label);
0955   }
0956   //BunchStereo
0957   if (switchBunchStereo) {
0958     HistoName = "TH1BunchStereo__" + label;
0959     subdetMEs.meBunchStereo = bookME1D(ibooker, "TH1BunchStereo", HistoName.c_str(), "Bunch Crossing");
0960     subdetMEs.meBunchStereo->setAxisTitle("Bunch crossing, stereo modules in " + label);
0961   }
0962   //EventStereo
0963   if (switchEventStereo) {
0964     HistoName = "TH1EventStereo__" + label;
0965     subdetMEs.meEventStereo = bookME1D(ibooker, "TH1EventStereo", HistoName.c_str(), "Event (in-time bunch)");
0966     subdetMEs.meEventStereo->setAxisTitle("Event (in-time bunch), stereo modules in " + label);
0967   }
0968   //NumMatched
0969   if (switchNumMatched) {
0970     HistoName = "TH1NumMatched__" + label;
0971     subdetMEs.meNumMatched = bookME1D(ibooker, "TH1NumMatched", HistoName.c_str(), "Num of matched RecHits");
0972     subdetMEs.meNumMatched->setAxisTitle("Total number of matched RecHits in " + label);
0973   }
0974   //BunchMatched
0975   if (switchBunchMatched) {
0976     HistoName = "TH1BunchMatched__" + label;
0977     subdetMEs.meBunchMatched = bookME1D(ibooker, "TH1BunchMatched", HistoName.c_str(), "Bunch Crossing");
0978     subdetMEs.meBunchMatched->setAxisTitle("Bunch crossing, matched RecHits in " + label);
0979   }
0980   //EventMatched
0981   if (switchEventMatched) {
0982     HistoName = "TH1EventMatched__" + label;
0983     subdetMEs.meEventMatched = bookME1D(ibooker, "TH1EventMatched", HistoName.c_str(), "Event (in-time bunch)");
0984     subdetMEs.meEventMatched->setAxisTitle("Event (in-time bunch), matched RecHits in " + label);
0985   }
0986 
0987   SubDetMEsMap[label] = subdetMEs;
0988 }
0989 //------------------------------------------------------------------------------------------
0990 inline SiStripRecHitsValid::MonitorElement* SiStripRecHitsValid::bookME1D(DQMStore::IBooker& ibooker,
0991                                                                           const char* ParameterSetLabel,
0992                                                                           const char* HistoName,
0993                                                                           const char* HistoTitle) {
0994   edm::ParameterSet parameters = conf_.getParameter<edm::ParameterSet>(ParameterSetLabel);
0995   return ibooker.book1D(HistoName,
0996                         HistoTitle,
0997                         parameters.getParameter<int32_t>("Nbinx"),
0998                         parameters.getParameter<double>("xmin"),
0999                         parameters.getParameter<double>("xmax"));
1000 }
1001 
1002 //define this as a plug-in
1003 DEFINE_FWK_MODULE(SiStripRecHitsValid);