Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:09:41

0001 #include "DQMOffline/Muon/interface/MuonRecoAnalyzer.h"
0002 
0003 #include "DataFormats/Common/interface/Handle.h"
0004 #include "DataFormats/MuonReco/interface/Muon.h"
0005 #include "DataFormats/MuonReco/interface/MuonFwd.h"
0006 #include "DataFormats/MuonReco/interface/MuonEnergy.h"
0007 #include "DataFormats/Scalers/interface/DcsStatus.h"
0008 #include "DataFormats/TrackReco/interface/Track.h"
0009 #include "DataFormats/TrackReco/interface/TrackFwd.h"
0010 
0011 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0012 
0013 #include <string>
0014 #include "TMath.h"
0015 using namespace std;
0016 using namespace edm;
0017 
0018 MuonRecoAnalyzer::MuonRecoAnalyzer(const edm::ParameterSet& pSet) {
0019   parameters = pSet;
0020 
0021   // Input booleans
0022   IsminiAOD = parameters.getParameter<bool>("IsminiAOD");
0023   doMVA = parameters.getParameter<bool>("doMVA");
0024   useGEM = parameters.getUntrackedParameter<bool>("useGEM");
0025   maxGEMhitsSoftMuonMVA = parameters.getUntrackedParameter<int>("maxGEMhitsSoftMuonMVA");
0026 
0027   // the services:
0028   theMuonCollectionLabel_ = consumes<edm::View<reco::Muon> >(parameters.getParameter<edm::InputTag>("MuonCollection"));
0029   theVertexLabel_ = consumes<reco::VertexCollection>(pSet.getParameter<InputTag>("inputTagVertex"));
0030   theBeamSpotLabel_ = consumes<reco::BeamSpot>(pSet.getParameter<InputTag>("inputTagBeamSpot"));
0031   dcsStatusCollection_ =
0032       consumes<DcsStatusCollection>(pSet.getUntrackedParameter<std::string>("dcsStatusCollection", "scalersRawToDigi"));
0033 
0034   ptBin = parameters.getParameter<int>("ptBin");
0035   ptMin = parameters.getParameter<double>("ptMin");
0036   ptMax = parameters.getParameter<double>("ptMax");
0037   pResBin = parameters.getParameter<int>("pResBin");
0038   pResMin = parameters.getParameter<double>("pResMin");
0039   pResMax = parameters.getParameter<double>("pResMax");
0040   rhBin = parameters.getParameter<int>("rhBin");
0041   rhMin = parameters.getParameter<double>("rhMin");
0042   rhMax = parameters.getParameter<double>("rhMax");
0043   pBin = parameters.getParameter<int>("pBin");
0044   pMin = parameters.getParameter<double>("pMin");
0045   pMax = parameters.getParameter<double>("pMax");
0046   chi2Bin = parameters.getParameter<int>("chi2Bin");
0047   chi2Min = parameters.getParameter<double>("chi2Min");
0048   chi2Max = parameters.getParameter<double>("chi2Max");
0049   phiBin = parameters.getParameter<int>("phiBin");
0050   phiMin = parameters.getParameter<double>("phiMin");
0051   phiMax = parameters.getParameter<double>("phiMax");
0052   tunePBin = parameters.getParameter<int>("tunePBin");
0053   tunePMax = parameters.getParameter<double>("tunePMax");
0054   tunePMin = parameters.getParameter<double>("tunePMin");
0055   thetaBin = parameters.getParameter<int>("thetaBin");
0056   thetaMin = parameters.getParameter<double>("thetaMin");
0057   thetaMax = parameters.getParameter<double>("thetaMax");
0058   etaBin = parameters.getParameter<int>("etaBin");
0059   etaMin = parameters.getParameter<double>("etaMin");
0060   etaMax = parameters.getParameter<double>("etaMax");
0061 
0062   theFolder = parameters.getParameter<string>("folder");
0063 }
0064 
0065 MuonRecoAnalyzer::~MuonRecoAnalyzer() {}
0066 void MuonRecoAnalyzer::bookHistograms(DQMStore::IBooker& ibooker,
0067                                       edm::Run const& /*iRun*/,
0068                                       edm::EventSetup const& /* iSetup */) {
0069   ibooker.cd();
0070   ibooker.setCurrentFolder(theFolder);
0071 
0072   muReco = ibooker.book1D("muReco", "muon reconstructed tracks", 6, 1, 7);
0073   muReco->setBinLabel(1, "glb+tk+sta");
0074   muReco->setBinLabel(2, "glb+sta");
0075   muReco->setBinLabel(3, "tk+sta");
0076   muReco->setBinLabel(4, "tk");
0077   muReco->setBinLabel(5, "sta");
0078   muReco->setBinLabel(6, "calo");
0079 
0080   int binFactor = 4;
0081 
0082   /////////////////////////////////////////////////////
0083   // monitoring of eta parameter
0084   /////////////////////////////////////////////////////
0085   std::string histname = "GlbMuon_";
0086   etaGlbTrack.push_back(ibooker.book1D(histname + "Glb_eta", "#eta_{GLB}", etaBin, etaMin, etaMax));
0087   etaGlbTrack.push_back(ibooker.book1D(histname + "Tk_eta", "#eta_{TKfromGLB}", etaBin, etaMin, etaMax));
0088   etaGlbTrack.push_back(ibooker.book1D(histname + "Sta_eta", "#eta_{STAfromGLB}", etaBin, etaMin, etaMax));
0089   etaResolution.push_back(ibooker.book1D(
0090       "Res_TkGlb_eta", "#eta_{TKfromGLB} - #eta_{GLB}", etaBin * binFactor, etaMin / 3000, etaMax / 3000));
0091   etaResolution.push_back(ibooker.book1D(
0092       "Res_GlbSta_eta", "#eta_{GLB} - #eta_{STAfromGLB}", etaBin * binFactor, etaMin / 100, etaMax / 100));
0093   etaResolution.push_back(ibooker.book1D(
0094       "Res_TkSta_eta", "#eta_{TKfromGLB} - #eta_{STAfromGLB}", etaBin * binFactor, etaMin / 100, etaMax / 100));
0095   etaResolution.push_back(ibooker.book2D("ResVsEta_TkGlb_eta",
0096                                          "(#eta_{TKfromGLB} - #eta_{GLB}) vs #eta_{GLB}",
0097                                          etaBin,
0098                                          etaMin,
0099                                          etaMax,
0100                                          etaBin * binFactor,
0101                                          etaMin / 3000,
0102                                          etaMax / 3000));
0103   etaResolution.push_back(ibooker.book2D("ResVsEta_GlbSta_eta",
0104                                          "(#eta_{GLB} - #eta_{STAfromGLB}) vs #eta_{GLB}",
0105                                          etaBin,
0106                                          etaMin,
0107                                          etaMax,
0108                                          etaBin * binFactor,
0109                                          etaMin / 100,
0110                                          etaMax / 100));
0111   etaResolution.push_back(ibooker.book2D("ResVsEta_TkSta_eta",
0112                                          "(#eta_{TKfromGLB} - #eta_{STAfromGLB}) vs #eta_{TKfromGLB}",
0113                                          etaBin,
0114                                          etaMin,
0115                                          etaMax,
0116                                          etaBin * binFactor,
0117                                          etaMin / 100,
0118                                          etaMax / 100));
0119   etaPull = ibooker.book1D("Pull_TkSta_eta", "#eta_{TKfromGLB} - #eta_{GLB} / error", 100, -10, 10);
0120   etaTrack = ibooker.book1D("TkMuon_eta", "#eta_{TK}", etaBin, etaMin, etaMax);
0121   etaStaTrack = ibooker.book1D("StaMuon_eta", "#eta_{STA}", etaBin, etaMin, etaMax);
0122   etaEfficiency.push_back(ibooker.book1D("StaEta", "#eta_{STAfromGLB}", etaBin, etaMin, etaMax));
0123   etaEfficiency.push_back(
0124       ibooker.book1D("StaEta_ifCombinedAlso", "#eta_{STAfromGLB} if isGlb=true", etaBin, etaMin, etaMax));
0125 
0126   //////////////////////////////////////////////////////
0127   // monitoring of theta parameter
0128   /////////////////////////////////////////////////////
0129   thetaGlbTrack.push_back(ibooker.book1D(histname + "Glb_theta", "#theta_{GLB}", thetaBin, thetaMin, thetaMax));
0130   thetaGlbTrack[0]->setAxisTitle("rad");
0131   thetaGlbTrack.push_back(ibooker.book1D(histname + "Tk_theta", "#theta_{TKfromGLB}", thetaBin, thetaMin, thetaMax));
0132   thetaGlbTrack[1]->setAxisTitle("rad");
0133   thetaGlbTrack.push_back(ibooker.book1D(histname + "Sta_theta", "#theta_{STAfromGLB}", thetaBin, thetaMin, thetaMax));
0134   thetaGlbTrack[2]->setAxisTitle("rad");
0135   thetaResolution.push_back(ibooker.book1D("Res_TkGlb_theta",
0136                                            "#theta_{TKfromGLB} - #theta_{GLB}",
0137                                            thetaBin * binFactor,
0138                                            -(thetaMax / 3000),
0139                                            thetaMax / 3000));
0140   thetaResolution[0]->setAxisTitle("rad");
0141   thetaResolution.push_back(ibooker.book1D("Res_GlbSta_theta",
0142                                            "#theta_{GLB} - #theta_{STAfromGLB}",
0143                                            thetaBin * binFactor,
0144                                            -(thetaMax / 100),
0145                                            thetaMax / 100));
0146   thetaResolution[1]->setAxisTitle("rad");
0147   thetaResolution.push_back(ibooker.book1D("Res_TkSta_theta",
0148                                            "#theta_{TKfromGLB} - #theta_{STAfromGLB}",
0149                                            thetaBin * binFactor,
0150                                            -(thetaMax / 100),
0151                                            thetaMax / 100));
0152   thetaResolution[2]->setAxisTitle("rad");
0153   thetaResolution.push_back(ibooker.book2D("ResVsTheta_TkGlb_theta",
0154                                            "(#theta_{TKfromGLB} - #theta_{GLB}) vs #theta_{GLB}",
0155                                            thetaBin,
0156                                            thetaMin,
0157                                            thetaMax,
0158                                            thetaBin * binFactor,
0159                                            -(thetaMax / 3000),
0160                                            thetaMax / 3000));
0161   thetaResolution[3]->setAxisTitle("rad", 1);
0162   thetaResolution[3]->setAxisTitle("rad", 2);
0163   thetaResolution.push_back(ibooker.book2D("ResVsTheta_GlbSta_theta",
0164                                            "(#theta_{GLB} - #theta_{STAfromGLB}) vs #theta_{GLB}",
0165                                            thetaBin,
0166                                            thetaMin,
0167                                            thetaMax,
0168                                            thetaBin * binFactor,
0169                                            -(thetaMax / 100),
0170                                            thetaMax / 100));
0171   thetaResolution[4]->setAxisTitle("rad", 1);
0172   thetaResolution[4]->setAxisTitle("rad", 2);
0173   thetaResolution.push_back(ibooker.book2D("ResVsTheta_TkSta_theta",
0174                                            "(#theta_{TKfromGLB} - #theta_{STAfromGLB}) vs #theta_{TKfromGLB}",
0175                                            thetaBin,
0176                                            thetaMin,
0177                                            thetaMax,
0178                                            thetaBin * binFactor,
0179                                            -(thetaMax / 100),
0180                                            thetaMax / 100));
0181   thetaResolution[5]->setAxisTitle("rad", 1);
0182   thetaResolution[5]->setAxisTitle("rad", 2);
0183   thetaPull = ibooker.book1D("Pull_TkSta_theta", "#theta_{TKfromGLB} - #theta_{STAfromGLB} / error", 100, -10, 10);
0184   thetaTrack = ibooker.book1D("TkMuon_theta", "#theta_{TK}", thetaBin, thetaMin, thetaMax);
0185   thetaTrack->setAxisTitle("rad");
0186   thetaStaTrack = ibooker.book1D("StaMuon_theta", "#theta_{STA}", thetaBin, thetaMin, thetaMax);
0187   thetaStaTrack->setAxisTitle("rad");
0188 
0189   // monitoring tunePMuonBestTrack Pt
0190   tunePResolution = ibooker.book1D(
0191       "Res_TuneP_pt", "Pt_{MuonBestTrack}-Pt_{tunePMuonBestTrack}/Pt_{MuonBestTrack}", tunePBin, tunePMin, tunePMax);
0192 
0193   // monitoring of phi paramater
0194   phiGlbTrack.push_back(ibooker.book1D(histname + "Glb_phi", "#phi_{GLB}", phiBin, phiMin, phiMax));
0195   phiGlbTrack[0]->setAxisTitle("rad");
0196   phiGlbTrack.push_back(ibooker.book1D(histname + "Tk_phi", "#phi_{TKfromGLB}", phiBin, phiMin, phiMax));
0197   phiGlbTrack[1]->setAxisTitle("rad");
0198   phiGlbTrack.push_back(ibooker.book1D(histname + "Sta_phi", "#phi_{STAfromGLB}", phiBin, phiMin, phiMax));
0199   phiGlbTrack[2]->setAxisTitle("rad");
0200   phiResolution.push_back(ibooker.book1D(
0201       "Res_TkGlb_phi", "#phi_{TKfromGLB} - #phi_{GLB}", phiBin * binFactor, phiMin / 3000, phiMax / 3000));
0202   phiResolution[0]->setAxisTitle("rad");
0203   phiResolution.push_back(ibooker.book1D(
0204       "Res_GlbSta_phi", "#phi_{GLB} - #phi_{STAfromGLB}", phiBin * binFactor, phiMin / 100, phiMax / 100));
0205   phiResolution[1]->setAxisTitle("rad");
0206   phiResolution.push_back(ibooker.book1D(
0207       "Res_TkSta_phi", "#phi_{TKfromGLB} - #phi_{STAfromGLB}", phiBin * binFactor, phiMin / 100, phiMax / 100));
0208   phiResolution[2]->setAxisTitle("rad");
0209   phiResolution.push_back(ibooker.book2D("ResVsPhi_TkGlb_phi",
0210                                          "(#phi_{TKfromGLB} - #phi_{GLB}) vs #phi_GLB",
0211                                          phiBin,
0212                                          phiMin,
0213                                          phiMax,
0214                                          phiBin * binFactor,
0215                                          phiMin / 3000,
0216                                          phiMax / 3000));
0217   phiResolution[3]->setAxisTitle("rad", 1);
0218   phiResolution[3]->setAxisTitle("rad", 2);
0219   phiResolution.push_back(ibooker.book2D("ResVsPhi_GlbSta_phi",
0220                                          "(#phi_{GLB} - #phi_{STAfromGLB}) vs #phi_{GLB}",
0221                                          phiBin,
0222                                          phiMin,
0223                                          phiMax,
0224                                          phiBin * binFactor,
0225                                          phiMin / 100,
0226                                          phiMax / 100));
0227   phiResolution[4]->setAxisTitle("rad", 1);
0228   phiResolution[4]->setAxisTitle("rad", 2);
0229   phiResolution.push_back(ibooker.book2D("ResVsPhi_TkSta_phi",
0230                                          "(#phi_{TKfromGLB} - #phi_{STAfromGLB}) vs #phi_{TKfromGLB}",
0231                                          phiBin,
0232                                          phiMin,
0233                                          phiMax,
0234                                          phiBin * binFactor,
0235                                          phiMin / 100,
0236                                          phiMax / 100));
0237   phiResolution[5]->setAxisTitle("rad", 1);
0238   phiResolution[5]->setAxisTitle("rad", 2);
0239   phiPull = ibooker.book1D("Pull_TkSta_phi", "#phi_{TKfromGLB} - #phi_{STAfromGLB} / error", 100, -10, 10);
0240   phiTrack = ibooker.book1D("TkMuon_phi", "#phi_{TK}", phiBin, phiMin, phiMax);
0241   phiTrack->setAxisTitle("rad");
0242   phiStaTrack = ibooker.book1D("StaMuon_phi", "#phi_{STA}", phiBin, phiMin, phiMax);
0243   phiStaTrack->setAxisTitle("rad");
0244   phiEfficiency.push_back(ibooker.book1D("StaPhi", "#phi_{STAfromGLB}", phiBin, phiMin, phiMax));
0245   phiEfficiency[0]->setAxisTitle("rad");
0246   phiEfficiency.push_back(
0247       ibooker.book1D("StaPhi_ifCombinedAlso", "#phi_{STAfromGLB} if the isGlb=true", phiBin, phiMin, phiMax));
0248   phiEfficiency[1]->setAxisTitle("rad");
0249 
0250   // monitoring of the chi2 parameter
0251   chi2OvDFGlbTrack.push_back(
0252       ibooker.book1D(histname + "Glb_chi2OverDf", "#chi_{2}OverDF_{GLB}", chi2Bin, chi2Min, chi2Max));
0253   chi2OvDFGlbTrack.push_back(
0254       ibooker.book1D(histname + "Tk_chi2OverDf", "#chi_{2}OverDF_{TKfromGLB}", phiBin, chi2Min, chi2Max));
0255   chi2OvDFGlbTrack.push_back(
0256       ibooker.book1D(histname + "Sta_chi2OverDf", "#chi_{2}OverDF_{STAfromGLB}", chi2Bin, chi2Min, chi2Max));
0257   chi2OvDFTrack = ibooker.book1D("TkMuon_chi2OverDf", "#chi_{2}OverDF_{TK}", chi2Bin, chi2Min, chi2Max);
0258   chi2OvDFStaTrack = ibooker.book1D("StaMuon_chi2OverDf", "#chi_{2}OverDF_{STA}", chi2Bin, chi2Min, chi2Max);
0259   //--------------------------
0260   probchi2GlbTrack.push_back(ibooker.book1D(histname + "Glb_probchi", "Prob #chi_{GLB}", 120, chi2Min, 1.20));
0261   probchi2GlbTrack.push_back(ibooker.book1D(histname + "Tk_probchi", "Prob #chi_{TKfromGLB}", 120, chi2Min, 1.20));
0262   probchi2GlbTrack.push_back(ibooker.book1D(histname + "Sta_probchi", "Prob #chi_{STAfromGLB}", 120, chi2Min, 1.20));
0263   probchi2Track = ibooker.book1D("TkMuon_probchi", "Prob #chi_{TK}", 120, chi2Min, 1.20);
0264   probchi2StaTrack = ibooker.book1D("StaMuon_probchi", "Prob #chi_{STA}", 120, chi2Min, 1.20);
0265 
0266   // monitoring of the momentum
0267   pGlbTrack.push_back(ibooker.book1D(histname + "Glb_p", "p_{GLB}", pBin, pMin, pMax));
0268   pGlbTrack[0]->setAxisTitle("GeV");
0269   pGlbTrack.push_back(ibooker.book1D(histname + "Tk_p", "p_{TKfromGLB}", pBin, pMin, pMax));
0270   pGlbTrack[1]->setAxisTitle("GeV");
0271   pGlbTrack.push_back(ibooker.book1D(histname + "Sta_p", "p_{STAfromGLB}", pBin, pMin, pMax));
0272   pGlbTrack[2]->setAxisTitle("GeV");
0273   pTrack = ibooker.book1D("TkMuon_p", "p_{TK}", pBin, pMin, pMax);
0274   pTrack->setAxisTitle("GeV");
0275   pStaTrack = ibooker.book1D("StaMuon_p", "p_{STA}", pBin, pMin, pMax);
0276   pStaTrack->setAxisTitle("GeV");
0277 
0278   // monitoring of the transverse momentum
0279   ptGlbTrack.push_back(ibooker.book1D(histname + "Glb_pt", "pt_{GLB}", ptBin, ptMin, ptMax));
0280   ptGlbTrack[0]->setAxisTitle("GeV");
0281   ptGlbTrack.push_back(ibooker.book1D(histname + "Tk_pt", "pt_{TKfromGLB}", ptBin, ptMin, ptMax));
0282   ptGlbTrack[1]->setAxisTitle("GeV");
0283   ptGlbTrack.push_back(ibooker.book1D(histname + "Sta_pt", "pt_{STAfromGLB}", ptBin, ptMin, ptMax));
0284   ptGlbTrack[2]->setAxisTitle("GeV");
0285   ptTrack = ibooker.book1D("TkMuon_pt", "pt_{TK}", ptBin, ptMin, ptMax);
0286   ptTrack->setAxisTitle("GeV");
0287   ptStaTrack = ibooker.book1D("StaMuon_pt", "pt_{STA}", ptBin, ptMin, pMax);
0288   ptStaTrack->setAxisTitle("GeV");
0289 
0290   //monitoring of variables needed by the MVA soft muon
0291 
0292   ptSoftMuonMVA = ibooker.book1D("ptSoftMuonMVA", "pt_{SoftMuon}", 50, 0, 50);
0293   deltaRSoftMuonMVA = ibooker.book1D("deltaRSoftMuonMVA", "#Delta R", 50, 0, 5);
0294   gNchi2SoftMuonMVA = ibooker.book1D("gNchi2SoftMuonMVA", "gNchi2", 50, 0, 3);
0295   vMuHitsSoftMuonMVA = ibooker.book1D("vMuHitsSoftMuonMVA", "vMuHits", 50, 0, 50);
0296   mNuStationsSoftMuonMVA = ibooker.book1D("mNuStationsSoftMuonMVA", "mNuStations", 6, 0, 6);
0297   dxyRefSoftMuonMVA = ibooker.book1D("dxyRefSoftMuonMVA", "dxyRef", 50, -0.1, 0.1);
0298   dzRefSoftMuonMVA = ibooker.book1D("dzRefSoftMuonMVA", "dzRef", 50, -0.1, 0.1);
0299   LWHSoftMuonMVA = ibooker.book1D("LWHSoftMuonMVA", "LWH", 20, 0, 20);
0300   valPixHitsSoftMuonMVA = ibooker.book1D("valPixHitsSoftMuonMVA", "valPixHits", 8, 0, 8);
0301   innerChi2SoftMuonMVA = ibooker.book1D("innerChi2SoftMuonMVA", "innerChi2", 50, 0, 3);
0302   outerChi2SoftMuonMVA = ibooker.book1D("outerChi2SoftMuonMVA", "outerChi2", 50, 0, 4);
0303   iValFracSoftMuonMVA = ibooker.book1D("iValFracSoftMuonMVA", "iValFrac", 50, 0.5, 1.0);
0304   segCompSoftMuonMVA = ibooker.book1D("segCompSoftMuonMVA", "segComp", 50, 0, 1.2);
0305   chi2LocMomSoftMuonMVA = ibooker.book1D("chi2LocMomSoftMuonMVA", "chi2LocMom", 50, 0, 40);
0306   chi2LocPosSoftMuonMVA = ibooker.book1D("chi2LocPosSoftMuonMVA", "chi2LocPos", 0, 0, 8);
0307   glbTrackTailProbSoftMuonMVA = ibooker.book1D("glbTrackTailProbSoftMuonMVA", "glbTrackTailProb", 50, 0, 8);
0308   NTrkVHitsSoftMuonMVA = ibooker.book1D("NTrkVHitsSoftMuonMVA", "NTrkVHits", 50, 0, 35);
0309   kinkFinderSoftMuonMVA = ibooker.book1D("kinkFinderSoftMuonMVA", "kinkFinder", 50, 0, 30);
0310   vRPChitsSoftMuonMVA = ibooker.book1D("vRPChitsSoftMuonMVA", "vRPChits", 50, 0, 50);
0311   glbKinkFinderSoftMuonMVA = ibooker.book1D("glbKinkFinderSoftMuonMVA", "glbKinkFinder", 50, 0, 50);
0312   glbKinkFinderLogSoftMuonMVA = ibooker.book1D("glbKinkFinderLogSoftMuonMVA", "glbKinkFinderLog", 50, 0, 50);
0313   staRelChi2SoftMuonMVA = ibooker.book1D("staRelChi2SoftMuonMVA", "staRelChi2", 50, 0, 2);
0314   glbDeltaEtaPhiSoftMuonMVA = ibooker.book1D("glbDeltaEtaPhiSoftMuonMVA", "glbDeltaEtaPhi", 50, 0, 0.15);
0315   trkRelChi2SoftMuonMVA = ibooker.book1D("trkRelChi2SoftMuonMVA", "trkRelChi2", 50, 0, 1.2);
0316   vDThitsSoftMuonMVA = ibooker.book1D("vDThitsSoftMuonMVA", "vDThits", 50, 0, 50);
0317   vCSChitsSoftMuonMVA = ibooker.book1D("vCSChitsSoftMuonMVA", "vCSChits", 50, 0, 50);
0318   if (useGEM) {
0319     vGEMhitsSoftMuonMVA =
0320         ibooker.book1D("vGEMhitsSoftMuonMVA", "vGEMhits", maxGEMhitsSoftMuonMVA, 0, maxGEMhitsSoftMuonMVA);
0321     vGEMhitsSoftMuonMVA->setXTitle("Number of Valid GEM Hits of Global Muon");
0322   } else {
0323     vGEMhitsSoftMuonMVA = nullptr;
0324   }
0325   timeAtIpInOutSoftMuonMVA = ibooker.book1D("timeAtIpInOutSoftMuonMVA", "timeAtIpInOut", 50, -10.0, 10.0);
0326   timeAtIpInOutErrSoftMuonMVA = ibooker.book1D("timeAtIpInOutErrSoftMuonMVA", "timeAtIpInOutErr", 50, 0, 3.5);
0327   getMuonHitsPerStationSoftMuonMVA =
0328       ibooker.book1D("getMuonHitsPerStationSoftMuonMVA", "getMuonHitsPerStation", 6, 0, 6);
0329   QprodSoftMuonMVA = ibooker.book1D("QprodSoftMuonMVA", "Qprod", 4, -2, 2);
0330 
0331   // monitoring of the muon charge
0332   qGlbTrack.push_back(ibooker.book1D(histname + "Glb_q", "q_{GLB}", 5, -2.5, 2.5));
0333   qGlbTrack.push_back(ibooker.book1D(histname + "Tk_q", "q_{TKfromGLB}", 5, -2.5, 2.5));
0334   qGlbTrack.push_back(ibooker.book1D(histname + "Sta_q", "q_{STAformGLB}", 5, -2.5, 2.5));
0335   qGlbTrack.push_back(ibooker.book1D(
0336       histname + "qComparison", "comparison between q_{GLB} and q_{TKfromGLB}, q_{STAfromGLB}", 8, 0.5, 8.5));
0337   qGlbTrack[3]->setBinLabel(1, "qGlb=qSta");
0338   qGlbTrack[3]->setBinLabel(2, "qGlb!=qSta");
0339   qGlbTrack[3]->setBinLabel(3, "qGlb=qTk");
0340   qGlbTrack[3]->setBinLabel(4, "qGlb!=qTk");
0341   qGlbTrack[3]->setBinLabel(5, "qSta=qTk");
0342   qGlbTrack[3]->setBinLabel(6, "qSta!=qTk");
0343   qGlbTrack[3]->setBinLabel(7, "qGlb!=qSta,qGlb!=Tk");
0344   qGlbTrack[3]->setBinLabel(8, "qGlb=qSta,qGlb=Tk");
0345   qTrack = ibooker.book1D("TkMuon_q", "q_{TK}", 5, -2.5, 2.5);
0346   qStaTrack = ibooker.book1D("StaMuon_q", "q_{STA}", 5, -2.5, 2.5);
0347 
0348   //////////////////////////////////////////////////////////////
0349   // monitoring of the momentum resolution
0350   qOverpResolution.push_back(ibooker.book1D(
0351       "Res_TkGlb_qOverp", "(q/p)_{TKfromGLB} - (q/p)_{GLB}", pResBin * binFactor * 2, pResMin / 10, pResMax / 10));
0352   qOverpResolution[0]->setAxisTitle("GeV^{-1}");
0353   qOverpResolution.push_back(
0354       ibooker.book1D("Res_GlbSta_qOverp", "(q/p)_{GLB} - (q/p)_{STAfromGLB}", pResBin * binFactor, pResMin, pResMax));
0355   qOverpResolution[1]->setAxisTitle("GeV^{-1}");
0356   qOverpResolution.push_back(ibooker.book1D(
0357       "Res_TkSta_qOverp", "(q/p)_{TKfromGLB} - (q/p)_{STAfromGLB}", pResBin * binFactor, pResMin, pResMax));
0358   qOverpResolution[2]->setAxisTitle("GeV^{-1}");
0359   qOverpPull = ibooker.book1D("Pull_TkSta_qOverp", "(q/p)_{TKfromGLB} - (q/p)_{STAfromGLB} / error", 100, -10, 10);
0360 
0361   oneOverpResolution.push_back(ibooker.book1D(
0362       "Res_TkGlb_oneOverp", "(1/p)_{TKfromGLB} - (1/p)_{GLB}", pResBin * binFactor * 2, pResMin / 10, pResMax / 10));
0363   oneOverpResolution[0]->setAxisTitle("GeV^{-1}");
0364   oneOverpResolution.push_back(
0365       ibooker.book1D("Res_GlbSta_oneOverp", "(1/p)_{GLB} - (1/p)_{STAfromGLB}", pResBin * binFactor, pResMin, pResMax));
0366   oneOverpResolution[1]->setAxisTitle("GeV^{-1}");
0367   oneOverpResolution.push_back(ibooker.book1D(
0368       "Res_TkSta_oneOverp", "(q/p)_{TKfromGLB} - (q/p)_{STAfromGLB}", pResBin * binFactor, pResMin, pResMax));
0369   oneOverpResolution[2]->setAxisTitle("GeV^{-1}");
0370   oneOverpPull = ibooker.book1D("Pull_TkSta_oneOverp", "(1/p)_{TKfromGLB} - (1/p)_{STAfromGLB} / error", 100, -10, 10);
0371 
0372   qOverptResolution.push_back(ibooker.book1D("Res_TkGlb_qOverpt",
0373                                              "(q/p_{t})_{TKfromGLB} - (q/p_{t})_{GLB}",
0374                                              pResBin * binFactor * 2,
0375                                              pResMin / 10,
0376                                              pResMax / 10));
0377   qOverptResolution[0]->setAxisTitle("GeV^{-1}");
0378   qOverptResolution.push_back(ibooker.book1D(
0379       "Res_GlbSta_qOverpt", "(q/p_{t})_{GLB} - (q/p_{t})_{STAfromGLB}", pResBin * binFactor, pResMin, pResMax));
0380   qOverptResolution[1]->setAxisTitle("GeV^{-1}");
0381   qOverptResolution.push_back(ibooker.book1D(
0382       "Res_TkSta_qOverpt", "(q/p_{t})_{TKfromGLB} - (q/p_{t})_{STAfromGLB}", pResBin * binFactor, pResMin, pResMax));
0383   qOverptResolution[2]->setAxisTitle("GeV^{-1}");
0384   qOverptPull = ibooker.book1D("Pull_TkSta_qOverpt", "(q/pt)_{TKfromGLB} - (q/pt)_{STAfromGLB} / error", 100, -10, 10);
0385 
0386   oneOverptResolution.push_back(ibooker.book1D("Res_TkGlb_oneOverpt",
0387                                                "(1/p_{t})_{TKfromGLB} - (1/p_{t})_{GLB}",
0388                                                pResBin * binFactor * 2,
0389                                                pResMin / 10,
0390                                                pResMax / 10));
0391   oneOverptResolution[0]->setAxisTitle("GeV^{-1}");
0392   oneOverptResolution.push_back(ibooker.book1D(
0393       "Res_GlbSta_oneOverpt", "(1/p_{t})_{GLB} - (1/p_{t})_{STAfromGLB}", pResBin * binFactor, pResMin, pResMax));
0394   oneOverptResolution[1]->setAxisTitle("GeV^{-1}");
0395   oneOverptResolution.push_back(ibooker.book1D(
0396       "Res_TkSta_oneOverpt", "(1/p_{t})_{TKfromGLB} - (1/p_{t})_{STAfromGLB}", pResBin * binFactor, pResMin, pResMax));
0397   oneOverptResolution[2]->setAxisTitle("GeV^{-1}");
0398   oneOverptResolution.push_back(ibooker.book2D("ResVsEta_TkGlb_oneOverpt",
0399                                                "(#eta_{TKfromGLB} - #eta_{GLB}) vs (1/p_{t})_{GLB}",
0400                                                etaBin,
0401                                                etaMin,
0402                                                etaMax,
0403                                                pResBin * binFactor * 2,
0404                                                pResMin / 10,
0405                                                pResMax / 10));
0406   oneOverptResolution[3]->setAxisTitle("GeV^{-1}", 2);
0407   oneOverptResolution.push_back(ibooker.book2D("ResVsEta_GlbSta_oneOverpt",
0408                                                "(#eta_{GLB} - #eta_{STAfromGLB} vs (1/p_{t})_{GLB}",
0409                                                etaBin,
0410                                                etaMin,
0411                                                etaMax,
0412                                                pResBin * binFactor,
0413                                                pResMin,
0414                                                pResMax));
0415   oneOverptResolution[4]->setAxisTitle("GeV^{-1}", 2);
0416   oneOverptResolution.push_back(ibooker.book2D("ResVsEta_TkSta_oneOverpt",
0417                                                "(#eta_{TKfromGLB} - #eta_{STAfromGLB}) vs (1/p_{t})_{TKfromGLB}",
0418                                                etaBin,
0419                                                etaMin,
0420                                                etaMax,
0421                                                pResBin * binFactor,
0422                                                pResMin,
0423                                                pResMax));
0424   oneOverptResolution[5]->setAxisTitle("GeV^{-1}", 2);
0425   oneOverptResolution.push_back(ibooker.book2D("ResVsPhi_TkGlb_oneOverpt",
0426                                                "(#phi_{TKfromGLB} - #phi_{GLB}) vs (1/p_{t})_{GLB}",
0427                                                phiBin,
0428                                                phiMin,
0429                                                phiMax,
0430                                                pResBin * binFactor * 2,
0431                                                pResMin / 10,
0432                                                pResMax / 10));
0433   oneOverptResolution[6]->setAxisTitle("rad", 1);
0434   oneOverptResolution[6]->setAxisTitle("GeV^{-1}", 2);
0435   oneOverptResolution.push_back(ibooker.book2D("ResVsPhi_GlbSta_oneOverpt",
0436                                                "(#phi_{GLB} - #phi_{STAfromGLB} vs (1/p_{t})_{GLB}",
0437                                                phiBin,
0438                                                phiMin,
0439                                                phiMax,
0440                                                pResBin * binFactor,
0441                                                pResMin,
0442                                                pResMax));
0443   oneOverptResolution[7]->setAxisTitle("rad", 1);
0444   oneOverptResolution[7]->setAxisTitle("GeV^{-1}", 2);
0445   oneOverptResolution.push_back(ibooker.book2D("ResVsPhi_TkSta_oneOverpt",
0446                                                "(#phi_{TKfromGLB} - #phi_{STAfromGLB}) vs (1/p_{t})_{TKfromGLB}",
0447                                                phiBin,
0448                                                phiMin,
0449                                                phiMax,
0450                                                pResBin * binFactor,
0451                                                pResMin,
0452                                                pResMax));
0453   oneOverptResolution[8]->setAxisTitle("rad", 1);
0454   oneOverptResolution[8]->setAxisTitle("GeV^{-1}", 2);
0455   oneOverptResolution.push_back(ibooker.book2D("ResVsPt_TkGlb_oneOverpt",
0456                                                "((1/p_{t})_{TKfromGLB} - (1/p_{t})_{GLB}) vs (1/p_{t})_{GLB}",
0457                                                ptBin / 5,
0458                                                ptMin,
0459                                                ptMax / 100,
0460                                                pResBin * binFactor * 2,
0461                                                pResMin / 10,
0462                                                pResMax / 10));
0463   oneOverptResolution[9]->setAxisTitle("GeV^{-1}", 1);
0464   oneOverptResolution[9]->setAxisTitle("GeV^{-1}", 2);
0465   oneOverptResolution.push_back(ibooker.book2D("ResVsPt_GlbSta_oneOverpt",
0466                                                "((1/p_{t})_{GLB} - (1/p_{t})_{STAfromGLB} vs (1/p_{t})_{GLB}",
0467                                                ptBin / 5,
0468                                                ptMin,
0469                                                ptMax / 100,
0470                                                pResBin * binFactor,
0471                                                pResMin,
0472                                                pResMax));
0473   oneOverptResolution[10]->setAxisTitle("GeV^{-1}", 1);
0474   oneOverptResolution[10]->setAxisTitle("GeV^{-1}", 2);
0475   oneOverptResolution.push_back(
0476       ibooker.book2D("ResVsPt_TkSta_oneOverpt",
0477                      "((1/p_{t})_{TKfromGLB} - (1/p_{t})_{STAfromGLB}) vs (1/p_{t})_{TKfromGLB}",
0478                      ptBin / 5,
0479                      ptMin,
0480                      ptMax / 100,
0481                      pResBin * binFactor,
0482                      pResMin,
0483                      pResMax));
0484   oneOverptResolution[11]->setAxisTitle("GeV^{-1}", 1);
0485   oneOverptResolution[11]->setAxisTitle("GeV^{-1}", 2);
0486   oneOverptPull =
0487       ibooker.book1D("Pull_TkSta_oneOverpt", "(1/pt)_{TKfromGLB} - (1/pt)_{STAfromGLB} / error", 100, -10, 10);
0488 
0489   //////////////////////////////////////////////////////////////
0490   // monitoring of the phi-eta
0491   phiVsetaGlbTrack.push_back(ibooker.book2D(
0492       histname + "Glb_phiVSeta", "#phi vs #eta (GLB)", etaBin / 2, etaMin, etaMax, phiBin / 2, phiMin, phiMax));
0493   phiVsetaGlbTrack.push_back(ibooker.book2D(
0494       histname + "Tk_phiVSeta", "#phi vs #eta (TKfromGLB)", etaBin / 2, etaMin, etaMax, phiBin / 2, phiMin, phiMax));
0495   phiVsetaGlbTrack.push_back(ibooker.book2D(
0496       histname + "Sta_phiVseta", "#phi vs #eta (STAfromGLB)", etaBin / 2, etaMin, etaMax, phiBin / 2, phiMin, phiMax));
0497 
0498   phiVsetaGlbTrack_badlumi.push_back(ibooker.book2D(
0499       histname + "Glb_phiVSeta_badlumi", "#phi vs #eta (GLB)", etaBin / 2, etaMin, etaMax, phiBin / 2, phiMin, phiMax));
0500   phiVsetaGlbTrack_badlumi.push_back(ibooker.book2D(histname + "Tk_phiVSeta_badlumi",
0501                                                     "#phi vs #eta (TKfromGLB)",
0502                                                     etaBin / 2,
0503                                                     etaMin,
0504                                                     etaMax,
0505                                                     phiBin / 2,
0506                                                     phiMin,
0507                                                     phiMax));
0508   phiVsetaGlbTrack_badlumi.push_back(ibooker.book2D(histname + "Sta_phiVseta_badlumi",
0509                                                     "#phi vs #eta (STAfromGLB)",
0510                                                     etaBin / 2,
0511                                                     etaMin,
0512                                                     etaMax,
0513                                                     phiBin / 2,
0514                                                     phiMin,
0515                                                     phiMax));
0516 
0517   //////////////////////////////////////////////////////////////
0518   // monitoring of the recHits provenance
0519   rhAnalysis.push_back(ibooker.book1D("StaRh_Frac_inGlb", "recHits_{STAinGLB} / recHits_{GLB}", rhBin, rhMin, rhMax));
0520   rhAnalysis.push_back(ibooker.book1D("TkRh_Frac_inGlb", "recHits_{TKinGLB} / recHits_{GLB}", rhBin, rhMin, rhMax));
0521   rhAnalysis.push_back(
0522       ibooker.book1D("StaRh_inGlb_Div_RhAssoSta", "recHits_{STAinGLB} / recHits_{STAfromGLB}", rhBin, rhMin, rhMax));
0523   rhAnalysis.push_back(
0524       ibooker.book1D("TkRh_inGlb_Div_RhAssoTk", "recHits_{TKinGLB} / recHits_{TKfromGLB}", rhBin, rhMin, rhMax));
0525   rhAnalysis.push_back(ibooker.book1D(
0526       "GlbRh_Div_RhAssoStaTk", "recHits_{GLB} / (recHits_{TKfromGLB}+recHits_{STAfromGLB})", rhBin, rhMin, rhMax));
0527   rhAnalysis.push_back(ibooker.book1D("invalidRh_Frac_inTk", "Invalid recHits / rechits_{GLB}", rhBin, rhMin, rhMax));
0528 
0529   //////////////////////////////////////////////////////////////
0530   // monitoring of the muon system rotation w.r.t. tracker
0531   muVStkSytemRotation.push_back(ibooker.book2D(
0532       "muVStkSytemRotation_posMu", "pT_{TK} / pT_{GLB} vs pT_{GLB} for #mu^{+}", 50, 0, 200, 100, 0.8, 1.2));
0533   muVStkSytemRotation.push_back(ibooker.book2D(
0534       "muVStkSytemRotation_negMu", "pT_{TK} / pT_{GLB} vs pT_{GLB} for #mu^{-}", 50, 0, 200, 100, 0.8, 1.2));
0535 }
0536 
0537 void MuonRecoAnalyzer::GetRes(reco::TrackRef t1, reco::TrackRef t2, string par, float& res, float& pull) {
0538   float p1 = 0, p2 = 0, p1e = 1, p2e = 1;
0539 
0540   if (par == "eta") {
0541     p1 = t1->eta();
0542     p1e = t1->etaError();
0543     p2 = t2->eta();
0544     p2e = t2->etaError();
0545   } else if (par == "theta") {
0546     p1 = t1->theta();
0547     p1e = t1->thetaError();
0548     p2 = t2->theta();
0549     p2e = t2->thetaError();
0550   } else if (par == "phi") {
0551     p1 = t1->phi();
0552     p1e = t1->phiError();
0553     p2 = t2->phi();
0554     p2e = t2->phiError();
0555   } else if (par == "qOverp") {
0556     p1 = t1->charge() / t1->p();
0557     p1e = t1->qoverpError();
0558     p2 = t2->charge() / t2->p();
0559     p2e = t2->qoverpError();
0560   } else if (par == "oneOverp") {
0561     p1 = 1. / t1->p();
0562     p1e = t1->qoverpError();
0563     p2 = 1. / t2->p();
0564     p2e = t2->qoverpError();
0565   } else if (par == "qOverpt") {
0566     p1 = t1->charge() / t1->pt();
0567     p1e = t1->ptError() * p1 * p1;
0568     p2 = t2->charge() / t2->pt();
0569     p2e = t2->ptError() * p2 * p2;
0570   } else if (par == "oneOverpt") {
0571     p1 = 1. / t1->pt();
0572     p1e = t1->ptError() * p1 * p1;
0573     p2 = 1. / t2->pt();
0574     p2e = t2->ptError() * p2 * p2;
0575   }
0576 
0577   res = p1 - p2;
0578   if (p1e != 0 || p2e != 0)
0579     pull = res / sqrt(p1e * p1e + p2e * p2e);
0580   else
0581     pull = -99;
0582   return;
0583 }
0584 
0585 void MuonRecoAnalyzer::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0586   LogTrace(metname) << "[MuonRecoAnalyzer] Analyze the mu";
0587 
0588   // Take the muon container
0589   edm::Handle<edm::View<reco::Muon> > muons;
0590   iEvent.getByToken(theMuonCollectionLabel_, muons);
0591 
0592   Handle<reco::BeamSpot> beamSpot;
0593   Handle<reco::VertexCollection> vertex;
0594   if (doMVA) {
0595     iEvent.getByToken(theBeamSpotLabel_, beamSpot);
0596     if (!beamSpot.isValid()) {
0597       edm::LogInfo("MuonRecoAnalyzer") << "Error: Can't get the beamspot" << endl;
0598       doMVA = false;
0599     }
0600     iEvent.getByToken(theVertexLabel_, vertex);
0601     if (!vertex.isValid()) {
0602       edm::LogInfo("MuonRecoAnalyzer") << "Error: Can't get the vertex collection" << endl;
0603       doMVA = false;
0604     }
0605   }
0606 
0607   //In this part we determine if we want to fill the plots for events where the DCS flag was set to bad
0608   edm::Handle<DcsStatusCollection> dcsStatus;
0609   bool fillBadLumi = false;
0610   if (iEvent.getByToken(dcsStatusCollection_, dcsStatus) && dcsStatus.isValid()) {
0611     for (auto const& dcsStatusItr : *dcsStatus) {
0612       if (!dcsStatusItr.ready(DcsStatus::CSCp))
0613         fillBadLumi = true;
0614       if (!dcsStatusItr.ready(DcsStatus::CSCm))
0615         fillBadLumi = true;
0616       if (!dcsStatusItr.ready(DcsStatus::DT0))
0617         fillBadLumi = true;
0618       if (!dcsStatusItr.ready(DcsStatus::DTp))
0619         fillBadLumi = true;
0620       if (!dcsStatusItr.ready(DcsStatus::DTm))
0621         fillBadLumi = true;
0622       if (!dcsStatusItr.ready(DcsStatus::EBp))
0623         fillBadLumi = true;
0624       if (!dcsStatusItr.ready(DcsStatus::EBm))
0625         fillBadLumi = true;
0626       if (!dcsStatusItr.ready(DcsStatus::EEp))
0627         fillBadLumi = true;
0628       if (!dcsStatusItr.ready(DcsStatus::EEm))
0629         fillBadLumi = true;
0630       if (!dcsStatusItr.ready(DcsStatus::ESp))
0631         fillBadLumi = true;
0632       if (!dcsStatusItr.ready(DcsStatus::ESm))
0633         fillBadLumi = true;
0634       if (!dcsStatusItr.ready(DcsStatus::HBHEa))
0635         fillBadLumi = true;
0636       if (!dcsStatusItr.ready(DcsStatus::HBHEb))
0637         fillBadLumi = true;
0638       if (!dcsStatusItr.ready(DcsStatus::HBHEc))
0639         fillBadLumi = true;
0640       if (!dcsStatusItr.ready(DcsStatus::HF))
0641         fillBadLumi = true;
0642       if (!dcsStatusItr.ready(DcsStatus::HO))
0643         fillBadLumi = true;
0644       if (!dcsStatusItr.ready(DcsStatus::BPIX))
0645         fillBadLumi = true;
0646       if (!dcsStatusItr.ready(DcsStatus::FPIX))
0647         fillBadLumi = true;
0648       if (!dcsStatusItr.ready(DcsStatus::RPC))
0649         fillBadLumi = true;
0650       if (!dcsStatusItr.ready(DcsStatus::TIBTID))
0651         fillBadLumi = true;
0652       if (!dcsStatusItr.ready(DcsStatus::TOB))
0653         fillBadLumi = true;
0654       if (!dcsStatusItr.ready(DcsStatus::TECp))
0655         fillBadLumi = true;
0656       if (!dcsStatusItr.ready(DcsStatus::TECm))
0657         fillBadLumi = true;
0658       //if (!dcsStatusItr.ready(DcsStatus::CASTOR)) fillBadLumi = true;
0659     }
0660   }
0661 
0662   float res = 0, pull = 0;
0663   if (!muons.isValid())
0664     return;
0665 
0666   for (edm::View<reco::Muon>::const_iterator muon = muons->begin(); muon != muons->end(); ++muon) {
0667     //Needed for MVA soft muon
0668 
0669     reco::TrackRef gTrack = muon->globalTrack();
0670     reco::TrackRef iTrack = muon->innerTrack();
0671     reco::TrackRef oTrack = muon->outerTrack();
0672     if (iTrack.isNonnull() && oTrack.isNonnull() && gTrack.isNonnull()) {
0673       const reco::HitPattern gHits = gTrack->hitPattern();
0674       const reco::HitPattern iHits = iTrack->hitPattern();
0675       const reco::MuonQuality muonQuality = muon->combinedQuality();
0676       int pvIndex = 0;
0677       math::XYZPoint refPoint;
0678       if (doMVA) {
0679         pvIndex = getPv(iTrack.index(), &(*vertex));  //HFDumpUtitilies
0680         if (pvIndex > -1) {
0681           refPoint = vertex->at(pvIndex).position();
0682         } else {
0683           if (beamSpot.isValid()) {
0684             refPoint = beamSpot->position();
0685           } else {
0686             edm::LogInfo("MuonRecoAnalyzer") << "ERROR: No beam sport found!" << endl;
0687           }
0688         }
0689       }
0690       ptSoftMuonMVA->Fill(iTrack->eta());
0691       deltaRSoftMuonMVA->Fill(getDeltaR(*iTrack, *oTrack));
0692       gNchi2SoftMuonMVA->Fill(gTrack->normalizedChi2());
0693       vMuHitsSoftMuonMVA->Fill(gHits.numberOfValidMuonHits());
0694       mNuStationsSoftMuonMVA->Fill(muon->numberOfMatchedStations());
0695       if (doMVA) {
0696         dxyRefSoftMuonMVA->Fill(iTrack->dxy(refPoint));
0697         dzRefSoftMuonMVA->Fill(iTrack->dz(refPoint));
0698       }
0699       LWHSoftMuonMVA->Fill(iHits.trackerLayersWithMeasurement());
0700       valPixHitsSoftMuonMVA->Fill(iHits.numberOfValidPixelHits());
0701       innerChi2SoftMuonMVA->Fill(iTrack->normalizedChi2());
0702       outerChi2SoftMuonMVA->Fill(oTrack->normalizedChi2());
0703       iValFracSoftMuonMVA->Fill(iTrack->validFraction());
0704       //segCompSoftMuonMVA->Fill(reco::Muon::segmentCompatibility(*muon));
0705       chi2LocMomSoftMuonMVA->Fill(muonQuality.chi2LocalMomentum);
0706       chi2LocPosSoftMuonMVA->Fill(muonQuality.chi2LocalPosition);
0707       glbTrackTailProbSoftMuonMVA->Fill(muonQuality.glbTrackProbability);
0708       NTrkVHitsSoftMuonMVA->Fill(iHits.numberOfValidTrackerHits());
0709       kinkFinderSoftMuonMVA->Fill(muonQuality.trkKink);
0710       vRPChitsSoftMuonMVA->Fill(gHits.numberOfValidMuonRPCHits());
0711       glbKinkFinderSoftMuonMVA->Fill(muonQuality.glbKink);
0712       glbKinkFinderLogSoftMuonMVA->Fill(TMath::Log(2 + muonQuality.glbKink));
0713       staRelChi2SoftMuonMVA->Fill(muonQuality.staRelChi2);
0714       glbDeltaEtaPhiSoftMuonMVA->Fill(muonQuality.globalDeltaEtaPhi);
0715       trkRelChi2SoftMuonMVA->Fill(muonQuality.trkRelChi2);
0716       vDThitsSoftMuonMVA->Fill(gHits.numberOfValidMuonDTHits());
0717       vCSChitsSoftMuonMVA->Fill(gHits.numberOfValidMuonCSCHits());
0718       if (useGEM) {
0719         vGEMhitsSoftMuonMVA->Fill(gHits.numberOfValidMuonGEMHits());
0720       }
0721       timeAtIpInOutSoftMuonMVA->Fill(muon->time().timeAtIpInOut);
0722       timeAtIpInOutErrSoftMuonMVA->Fill(muon->time().timeAtIpInOutErr);
0723       //getMuonHitsPerStationSoftMuonMVA->Fill(gTrack);
0724       QprodSoftMuonMVA->Fill((iTrack->charge() * oTrack->charge()));
0725     }
0726 
0727     if (muon->isGlobalMuon()) {
0728       LogTrace(metname) << "[MuonRecoAnalyzer] The mu is global - filling the histos";
0729       if (muon->isTrackerMuon() && muon->isStandAloneMuon())
0730         muReco->Fill(1);
0731       if (!(muon->isTrackerMuon()) && muon->isStandAloneMuon())
0732         muReco->Fill(2);
0733       if (!muon->isStandAloneMuon())
0734         LogTrace(metname) << "[MuonRecoAnalyzer] ERROR: the mu is global but not standalone!";
0735       // get the track combinig the information from both the Tracker and the Spectrometer
0736       reco::TrackRef recoCombinedGlbTrack = muon->combinedMuon();
0737 
0738       // get the track using only the tracker data
0739       reco::TrackRef recoTkGlbTrack = muon->track();
0740       // get the track using only the mu spectrometer data
0741       reco::TrackRef recoStaGlbTrack = muon->standAloneMuon();
0742       etaGlbTrack[0]->Fill(recoCombinedGlbTrack->eta());
0743       etaGlbTrack[1]->Fill(recoTkGlbTrack->eta());
0744       etaGlbTrack[2]->Fill(recoStaGlbTrack->eta());
0745 
0746       phiVsetaGlbTrack[0]->Fill(recoCombinedGlbTrack->eta(), recoCombinedGlbTrack->phi());
0747       phiVsetaGlbTrack[1]->Fill(recoTkGlbTrack->eta(), recoTkGlbTrack->phi());
0748       phiVsetaGlbTrack[2]->Fill(recoStaGlbTrack->eta(), recoStaGlbTrack->phi());
0749 
0750       if (fillBadLumi) {
0751         phiVsetaGlbTrack_badlumi[0]->Fill(recoCombinedGlbTrack->eta(), recoCombinedGlbTrack->phi());
0752         phiVsetaGlbTrack_badlumi[1]->Fill(recoTkGlbTrack->eta(), recoTkGlbTrack->phi());
0753         phiVsetaGlbTrack_badlumi[2]->Fill(recoStaGlbTrack->eta(), recoStaGlbTrack->phi());
0754       }
0755 
0756       GetRes(recoTkGlbTrack, recoCombinedGlbTrack, "eta", res, pull);
0757       etaResolution[0]->Fill(res);
0758       GetRes(recoCombinedGlbTrack, recoStaGlbTrack, "eta", res, pull);
0759       etaResolution[1]->Fill(res);
0760       GetRes(recoTkGlbTrack, recoStaGlbTrack, "eta", res, pull);
0761       etaResolution[2]->Fill(res);
0762       etaPull->Fill(pull);
0763       etaResolution[3]->Fill(recoCombinedGlbTrack->eta(), recoTkGlbTrack->eta() - recoCombinedGlbTrack->eta());
0764       etaResolution[4]->Fill(recoCombinedGlbTrack->eta(), -recoStaGlbTrack->eta() + recoCombinedGlbTrack->eta());
0765       etaResolution[5]->Fill(recoCombinedGlbTrack->eta(), recoTkGlbTrack->eta() - recoStaGlbTrack->eta());
0766 
0767       thetaGlbTrack[0]->Fill(recoCombinedGlbTrack->theta());
0768       thetaGlbTrack[1]->Fill(recoTkGlbTrack->theta());
0769       thetaGlbTrack[2]->Fill(recoStaGlbTrack->theta());
0770       GetRes(recoTkGlbTrack, recoCombinedGlbTrack, "theta", res, pull);
0771       thetaResolution[0]->Fill(res);
0772       GetRes(recoCombinedGlbTrack, recoStaGlbTrack, "theta", res, pull);
0773       thetaResolution[1]->Fill(res);
0774 
0775       GetRes(recoTkGlbTrack, recoStaGlbTrack, "theta", res, pull);
0776       thetaResolution[2]->Fill(res);
0777       thetaPull->Fill(pull);
0778       thetaResolution[3]->Fill(recoCombinedGlbTrack->theta(), recoTkGlbTrack->theta() - recoCombinedGlbTrack->theta());
0779       thetaResolution[4]->Fill(recoCombinedGlbTrack->theta(),
0780                                -recoStaGlbTrack->theta() + recoCombinedGlbTrack->theta());
0781       thetaResolution[5]->Fill(recoCombinedGlbTrack->theta(), recoTkGlbTrack->theta() - recoStaGlbTrack->theta());
0782 
0783       phiGlbTrack[0]->Fill(recoCombinedGlbTrack->phi());
0784       phiGlbTrack[1]->Fill(recoTkGlbTrack->phi());
0785       phiGlbTrack[2]->Fill(recoStaGlbTrack->phi());
0786       GetRes(recoTkGlbTrack, recoCombinedGlbTrack, "phi", res, pull);
0787       phiResolution[0]->Fill(res);
0788       GetRes(recoCombinedGlbTrack, recoStaGlbTrack, "phi", res, pull);
0789       phiResolution[1]->Fill(res);
0790       GetRes(recoTkGlbTrack, recoStaGlbTrack, "phi", res, pull);
0791       phiResolution[2]->Fill(res);
0792       phiPull->Fill(pull);
0793       phiResolution[3]->Fill(recoCombinedGlbTrack->phi(), recoTkGlbTrack->phi() - recoCombinedGlbTrack->phi());
0794       phiResolution[4]->Fill(recoCombinedGlbTrack->phi(), -recoStaGlbTrack->phi() + recoCombinedGlbTrack->phi());
0795       phiResolution[5]->Fill(recoCombinedGlbTrack->phi(), recoTkGlbTrack->phi() - recoStaGlbTrack->phi());
0796 
0797       chi2OvDFGlbTrack[0]->Fill(recoCombinedGlbTrack->normalizedChi2());
0798       chi2OvDFGlbTrack[1]->Fill(recoTkGlbTrack->normalizedChi2());
0799       chi2OvDFGlbTrack[2]->Fill(recoStaGlbTrack->normalizedChi2());
0800       //-------------------------
0801       //    double probchi = TMath::Prob(recoCombinedGlbTrack->normalizedChi2(),recoCombinedGlbTrack->ndof());
0802       //    cout << "rellenando histos."<<endl;
0803       probchi2GlbTrack[0]->Fill(TMath::Prob(recoCombinedGlbTrack->chi2(), recoCombinedGlbTrack->ndof()));
0804       probchi2GlbTrack[1]->Fill(TMath::Prob(recoTkGlbTrack->chi2(), recoTkGlbTrack->ndof()));
0805       probchi2GlbTrack[2]->Fill(TMath::Prob(recoStaGlbTrack->chi2(), recoStaGlbTrack->ndof()));
0806       //    cout << "rellenados histos."<<endl;
0807       //-------------------------
0808 
0809       pGlbTrack[0]->Fill(recoCombinedGlbTrack->p());
0810       pGlbTrack[1]->Fill(recoTkGlbTrack->p());
0811       pGlbTrack[2]->Fill(recoStaGlbTrack->p());
0812 
0813       ptGlbTrack[0]->Fill(recoCombinedGlbTrack->pt());
0814       ptGlbTrack[1]->Fill(recoTkGlbTrack->pt());
0815       ptGlbTrack[2]->Fill(recoStaGlbTrack->pt());
0816 
0817       qGlbTrack[0]->Fill(recoCombinedGlbTrack->charge());
0818       qGlbTrack[1]->Fill(recoTkGlbTrack->charge());
0819       qGlbTrack[2]->Fill(recoStaGlbTrack->charge());
0820       if (recoCombinedGlbTrack->charge() == recoStaGlbTrack->charge())
0821         qGlbTrack[3]->Fill(1);
0822       else
0823         qGlbTrack[3]->Fill(2);
0824       if (recoCombinedGlbTrack->charge() == recoTkGlbTrack->charge())
0825         qGlbTrack[3]->Fill(3);
0826       else
0827         qGlbTrack[3]->Fill(4);
0828       if (recoStaGlbTrack->charge() == recoTkGlbTrack->charge())
0829         qGlbTrack[3]->Fill(5);
0830       else
0831         qGlbTrack[3]->Fill(6);
0832       if (recoCombinedGlbTrack->charge() != recoStaGlbTrack->charge() &&
0833           recoCombinedGlbTrack->charge() != recoTkGlbTrack->charge())
0834         qGlbTrack[3]->Fill(7);
0835       if (recoCombinedGlbTrack->charge() == recoStaGlbTrack->charge() &&
0836           recoCombinedGlbTrack->charge() == recoTkGlbTrack->charge())
0837         qGlbTrack[3]->Fill(8);
0838 
0839       GetRes(recoTkGlbTrack, recoCombinedGlbTrack, "qOverp", res, pull);
0840       qOverpResolution[0]->Fill(res);
0841       GetRes(recoCombinedGlbTrack, recoStaGlbTrack, "qOverp", res, pull);
0842       qOverpResolution[1]->Fill(res);
0843       GetRes(recoTkGlbTrack, recoStaGlbTrack, "qOverp", res, pull);
0844       qOverpResolution[2]->Fill(res);
0845       qOverpPull->Fill(pull);
0846 
0847       GetRes(recoTkGlbTrack, recoCombinedGlbTrack, "oneOverp", res, pull);
0848       oneOverpResolution[0]->Fill(res);
0849       GetRes(recoCombinedGlbTrack, recoStaGlbTrack, "oneOverp", res, pull);
0850       oneOverpResolution[1]->Fill(res);
0851       GetRes(recoTkGlbTrack, recoStaGlbTrack, "oneOverp", res, pull);
0852       oneOverpResolution[2]->Fill(res);
0853       oneOverpPull->Fill(pull);
0854 
0855       GetRes(recoTkGlbTrack, recoCombinedGlbTrack, "qOverpt", res, pull);
0856       qOverptResolution[0]->Fill(res);
0857       GetRes(recoCombinedGlbTrack, recoStaGlbTrack, "qOverpt", res, pull);
0858       qOverptResolution[1]->Fill(res);
0859       GetRes(recoTkGlbTrack, recoStaGlbTrack, "qOverpt", res, pull);
0860       qOverptResolution[2]->Fill(res);
0861       qOverptPull->Fill(pull);
0862 
0863       GetRes(recoTkGlbTrack, recoCombinedGlbTrack, "oneOverpt", res, pull);
0864       oneOverptResolution[0]->Fill(res);
0865       GetRes(recoCombinedGlbTrack, recoStaGlbTrack, "oneOverpt", res, pull);
0866       oneOverptResolution[1]->Fill(res);
0867       GetRes(recoTkGlbTrack, recoStaGlbTrack, "oneOverpt", res, pull);
0868       oneOverptResolution[2]->Fill(res);
0869       oneOverptPull->Fill(pull);
0870 
0871       // //--- Test new tunePMuonBestTrack() method from Muon.h
0872 
0873       reco::TrackRef recoBestTrack = muon->muonBestTrack();
0874 
0875       reco::TrackRef recoTunePBestTrack = muon->tunePMuonBestTrack();
0876 
0877       double bestTrackPt = recoBestTrack->pt();
0878 
0879       double tunePBestTrackPt = recoTunePBestTrack->pt();
0880 
0881       double tunePBestTrackRes = (bestTrackPt - tunePBestTrackPt) / bestTrackPt;
0882 
0883       tunePResolution->Fill(tunePBestTrackRes);
0884 
0885       oneOverptResolution[3]->Fill(recoCombinedGlbTrack->eta(),
0886                                    (1 / recoTkGlbTrack->pt()) - (1 / recoCombinedGlbTrack->pt()));
0887       oneOverptResolution[4]->Fill(recoCombinedGlbTrack->eta(),
0888                                    -(1 / recoStaGlbTrack->pt()) + (1 / recoCombinedGlbTrack->pt()));
0889       oneOverptResolution[5]->Fill(recoCombinedGlbTrack->eta(),
0890                                    (1 / recoTkGlbTrack->pt()) - (1 / recoStaGlbTrack->pt()));
0891       oneOverptResolution[6]->Fill(recoCombinedGlbTrack->phi(),
0892                                    (1 / recoTkGlbTrack->pt()) - (1 / recoCombinedGlbTrack->pt()));
0893       oneOverptResolution[7]->Fill(recoCombinedGlbTrack->phi(),
0894                                    -(1 / recoStaGlbTrack->pt()) + (1 / recoCombinedGlbTrack->pt()));
0895       oneOverptResolution[8]->Fill(recoCombinedGlbTrack->phi(),
0896                                    (1 / recoTkGlbTrack->pt()) - (1 / recoStaGlbTrack->pt()));
0897       oneOverptResolution[9]->Fill(recoCombinedGlbTrack->pt(),
0898                                    (1 / recoTkGlbTrack->pt()) - (1 / recoCombinedGlbTrack->pt()));
0899       oneOverptResolution[10]->Fill(recoCombinedGlbTrack->pt(),
0900                                     -(1 / recoStaGlbTrack->pt()) + (1 / recoCombinedGlbTrack->pt()));
0901       oneOverptResolution[11]->Fill(recoCombinedGlbTrack->pt(),
0902                                     (1 / recoTkGlbTrack->pt()) - (1 / recoStaGlbTrack->pt()));
0903 
0904       if (!IsminiAOD) {
0905         // valid hits Glb track
0906         double rhGlb = recoCombinedGlbTrack->found();
0907         // valid hits Glb track from Tracker
0908         double rhGlb_StaProvenance = 0;
0909         // valid hits Glb track from Sta system
0910         double rhGlb_TkProvenance = 0;
0911 
0912         for (trackingRecHit_iterator recHit = recoCombinedGlbTrack->recHitsBegin();
0913              recHit != recoCombinedGlbTrack->recHitsEnd();
0914              ++recHit) {
0915           if ((*recHit)->isValid()) {
0916             DetId id = (*recHit)->geographicalId();
0917             if (id.det() == DetId::Muon)
0918               rhGlb_StaProvenance++;
0919             if (id.det() == DetId::Tracker)
0920               rhGlb_TkProvenance++;
0921           }
0922         }
0923         // valid hits Sta track associated to Glb track
0924         double rhStaGlb = recoStaGlbTrack->recHitsSize();
0925         // valid hits Traker track associated to Glb track
0926         double rhTkGlb = recoTkGlbTrack->found();
0927         // invalid hits Traker track associated to Glb track
0928         double rhTkGlb_notValid = recoTkGlbTrack->lost();
0929 
0930         // fill the histos
0931         rhAnalysis[0]->Fill(rhGlb_StaProvenance / rhGlb);
0932         rhAnalysis[1]->Fill(rhGlb_TkProvenance / rhGlb);
0933         rhAnalysis[2]->Fill(rhGlb_StaProvenance / rhStaGlb);
0934         rhAnalysis[3]->Fill(rhGlb_TkProvenance / rhTkGlb);
0935         rhAnalysis[4]->Fill(rhGlb / (rhStaGlb + rhTkGlb));
0936         rhAnalysis[5]->Fill(rhTkGlb_notValid / rhGlb);
0937       }
0938       // aligment plots (mu system w.r.t. tracker rotation)
0939       if (recoCombinedGlbTrack->charge() > 0)
0940         muVStkSytemRotation[0]->Fill(recoCombinedGlbTrack->pt(), recoTkGlbTrack->pt() / recoCombinedGlbTrack->pt());
0941       else
0942         muVStkSytemRotation[1]->Fill(recoCombinedGlbTrack->pt(), recoTkGlbTrack->pt() / recoCombinedGlbTrack->pt());
0943     }
0944 
0945     if (muon->isTrackerMuon() && !(muon->isGlobalMuon())) {
0946       LogTrace(metname) << "[MuonRecoAnalyzer] The mu is tracker only - filling the histos";
0947       if (muon->isStandAloneMuon())
0948         muReco->Fill(3);
0949       if (!(muon->isStandAloneMuon()))
0950         muReco->Fill(4);
0951 
0952       // get the track using only the tracker data
0953       reco::TrackRef recoTrack = muon->track();
0954 
0955       etaTrack->Fill(recoTrack->eta());
0956       thetaTrack->Fill(recoTrack->theta());
0957       phiTrack->Fill(recoTrack->phi());
0958       chi2OvDFTrack->Fill(recoTrack->normalizedChi2());
0959       probchi2Track->Fill(TMath::Prob(recoTrack->chi2(), recoTrack->ndof()));
0960       pTrack->Fill(recoTrack->p());
0961       ptTrack->Fill(recoTrack->pt());
0962       qTrack->Fill(recoTrack->charge());
0963     }
0964 
0965     if (muon->isStandAloneMuon() && !(muon->isGlobalMuon())) {
0966       LogTrace(metname) << "[MuonRecoAnalyzer] The mu is STA only - filling the histos";
0967       if (!(muon->isTrackerMuon()))
0968         muReco->Fill(5);
0969 
0970       // get the track using only the mu spectrometer data
0971       reco::TrackRef recoStaTrack = muon->standAloneMuon();
0972 
0973       etaStaTrack->Fill(recoStaTrack->eta());
0974       thetaStaTrack->Fill(recoStaTrack->theta());
0975       phiStaTrack->Fill(recoStaTrack->phi());
0976       chi2OvDFStaTrack->Fill(recoStaTrack->normalizedChi2());
0977       probchi2StaTrack->Fill(TMath::Prob(recoStaTrack->chi2(), recoStaTrack->ndof()));
0978       pStaTrack->Fill(recoStaTrack->p());
0979       ptStaTrack->Fill(recoStaTrack->pt());
0980       qStaTrack->Fill(recoStaTrack->charge());
0981     }
0982 
0983     if (muon->isCaloMuon() && !(muon->isGlobalMuon()) && !(muon->isTrackerMuon()) && !(muon->isStandAloneMuon()))
0984       muReco->Fill(6);
0985 
0986     //efficiency plots
0987 
0988     // get the track using only the mu spectrometer data
0989     reco::TrackRef recoStaGlbTrack = muon->standAloneMuon();
0990 
0991     if (muon->isStandAloneMuon()) {
0992       etaEfficiency[0]->Fill(recoStaGlbTrack->eta());
0993       phiEfficiency[0]->Fill(recoStaGlbTrack->phi());
0994     }
0995     if (muon->isStandAloneMuon() && muon->isGlobalMuon()) {
0996       etaEfficiency[1]->Fill(recoStaGlbTrack->eta());
0997       phiEfficiency[1]->Fill(recoStaGlbTrack->phi());
0998     }
0999   }
1000 }
1001 
1002 //Needed by MVA Soft Muon
1003 double MuonRecoAnalyzer::getDeltaR(reco::Track track1, reco::Track track2) {
1004   double dphi = acos(cos(track1.phi() - track2.phi()));
1005   double deta = track1.eta() - track2.eta();
1006   return sqrt(dphi * dphi + deta * deta);
1007 }
1008 
1009 // ----------------------------------------------------------------------
1010 int MuonRecoAnalyzer::getPv(int tidx, const reco::VertexCollection* vc) {
1011   if (vc) {
1012     for (unsigned int i = 0; i < vc->size(); ++i) {
1013       reco::Vertex::trackRef_iterator v1TrackIter;
1014       reco::Vertex::trackRef_iterator v1TrackBegin = vc->at(i).tracks_begin();
1015       reco::Vertex::trackRef_iterator v1TrackEnd = vc->at(i).tracks_end();
1016       for (v1TrackIter = v1TrackBegin; v1TrackIter != v1TrackEnd; v1TrackIter++) {
1017         if (static_cast<unsigned int>(tidx) == v1TrackIter->key())
1018           return i;
1019       }
1020     }
1021   }
1022   return -1;
1023 }