Line Code
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169

/*
 *  See header file for a description of this class.
 *
 *  \author G. Mila - INFN Torino
 */

#include <DQMOffline/Muon/interface/MuonRecoTest.h>

// Framework
#include <FWCore/Framework/interface/Event.h>
#include "DataFormats/Common/interface/Handle.h"
#include <FWCore/Framework/interface/ESHandle.h>
#include <FWCore/Framework/interface/MakerMacros.h>
#include <FWCore/Framework/interface/EventSetup.h>
#include <FWCore/ParameterSet/interface/ParameterSet.h>

#include "DQMServices/Core/interface/DQMStore.h"
#include "FWCore/MessageLogger/interface/MessageLogger.h"
#include "FWCore/Framework/interface/Run.h"

#include <iostream>
#include <cstdio>
#include <string>
#include <cmath>
#include "TF1.h"

using namespace edm;
using namespace std;

MuonRecoTest::MuonRecoTest(const edm::ParameterSet& ps) {
  parameters = ps;

  prescaleFactor = parameters.getUntrackedParameter<int>("diagnosticPrescale", 1);

  // Parameters
  etaBin = parameters.getParameter<int>("etaBin");
  etaMin = parameters.getParameter<double>("etaMin");
  etaMax = parameters.getParameter<double>("etaMax");

  phiBin = parameters.getParameter<int>("phiBin");
  phiMin = parameters.getParameter<double>("phiMin");
  phiMax = parameters.getParameter<double>("phiMax");

  EfficiencyCriterionName = parameters.getUntrackedParameter<string>("efficiencyTestName", "EfficiencyInRange");
}

void MuonRecoTest::dqmEndRun(DQMStore::IBooker& ibooker,
                             DQMStore::IGetter& igetter,
                             edm::Run const&,
                             edm::EventSetup const&) {
  // efficiency plot
  ibooker.setCurrentFolder("Muons/Tests/muonRecoTest");
  etaEfficiency = ibooker.book1D("etaEfficiency_staMuon", "#eta_{STA} efficiency", etaBin, etaMin, etaMax);
  phiEfficiency = ibooker.book1D("phiEfficiency_staMuon", "#phi_{STA} efficiency", phiBin, phiMin, phiMax);

  /// GETTING PREVIOUS HISTOS AND DO SOME OPERATIONS
  string path = "Muons/MuonRecoAnalyzer/StaEta_ifCombinedAlso";
  MonitorElement* staEtaIfComb_histo = igetter.get(path);
  path = "Muons/MuonRecoAnalyzer/StaEta";
  MonitorElement* staEta_histo = igetter.get(path);

  if (staEtaIfComb_histo && staEta_histo) {
    TH1F* staEtaIfComb_root = staEtaIfComb_histo->getTH1F();
    TH1F* staEta_root = staEta_histo->getTH1F();

    if (staEtaIfComb_root->GetXaxis()->GetNbins() != etaBin || staEtaIfComb_root->GetXaxis()->GetXmax() != etaMax ||
        staEtaIfComb_root->GetXaxis()->GetXmin() != etaMin) {
      LogTrace(metname) << "[MuonRecoTest] wrong histo binning on eta histograms";
      return;
    }

    for (int i = 1; i <= etaBin; i++) {
      if (staEta_root->GetBinContent(i) != 0) {
        double efficiency = double(staEtaIfComb_root->GetBinContent(i)) / double(staEta_root->GetBinContent(i));
        etaEfficiency->setBinContent(i, efficiency);
      }
    }
  }

  path = "Muons/MuonRecoAnalyzer/StaPhi_ifCombinedAlso";
  MonitorElement* staPhiIfComb_histo = igetter.get(path);
  path = "Muons/MuonRecoAnalyzer/StaPhi";
  MonitorElement* staPhi_histo = igetter.get(path);

  if (staPhiIfComb_histo && staPhi_histo) {
    TH1F* staPhiIfComb_root = staPhiIfComb_histo->getTH1F();
    TH1F* staPhi_root = staPhi_histo->getTH1F();

    if (staPhiIfComb_root->GetXaxis()->GetNbins() != phiBin || staPhiIfComb_root->GetXaxis()->GetXmax() != phiMax ||
        staPhiIfComb_root->GetXaxis()->GetXmin() != phiMin) {
      LogTrace(metname) << "[MuonRecoTest] wrong histo binning on phi histograms";
      return;
    }

    for (int i = 1; i <= etaBin; i++) {
      if (staPhi_root->GetBinContent(i) != 0) {
        double efficiency = double(staPhiIfComb_root->GetBinContent(i)) / double(staPhi_root->GetBinContent(i));
        phiEfficiency->setBinContent(i, efficiency);
      }
    }
  }
}

void MuonRecoTest::dqmEndJob(DQMStore::IBooker& ibooker, DQMStore::IGetter& igetter) {
  /// BOOKING NEW HISTOGRAMS
  ibooker.setCurrentFolder("Muons/Tests/muonRecoTest");
  // alignment plots
  globalRotation.push_back(ibooker.book1D(
      "muVStkSytemRotation_posMu_profile", "pT_{TK} / pT_{GLB} vs pT_{GLB} profile for #mu^{+}", 50, 0, 200));
  globalRotation.push_back(ibooker.book1D(
      "muVStkSytemRotation_negMu_profile", "pT_{TK} / pT_{GLB} vs pT_{GLB} profile for #mu^{-}", 50, 0, 200));
  globalRotation.push_back(ibooker.book1D(
      "muVStkSytemRotation_profile", "pT_{TK} / pT_{GLB} vs pT_{GLB} profile for #mu^{+}-#mu^{-}", 50, 0, 200));

  // efficiency test

  // eta efficiency
  const QReport* theEtaQReport = etaEfficiency->getQReport(EfficiencyCriterionName);
  if (theEtaQReport) {
    vector<dqm::me_util::Channel> badChannels = theEtaQReport->getBadChannels();
    for (vector<dqm::me_util::Channel>::iterator channel = badChannels.begin(); channel != badChannels.end();
         channel++) {
      LogTrace(metname) << "[etaEfficiency test] bad ranges: " << (*channel).getBin()
                        << "  Contents : " << (*channel).getContents() << endl;
    }
    LogTrace(metname) << "-------- type: [etaEfficiency]  " << theEtaQReport->getMessage() << " ------- "
                      << theEtaQReport->getStatus() << endl;
  }
  // phi efficiency
  const QReport* thePhiQReport = phiEfficiency->getQReport(EfficiencyCriterionName);
  if (thePhiQReport) {
    vector<dqm::me_util::Channel> badChannels = thePhiQReport->getBadChannels();
    for (vector<dqm::me_util::Channel>::iterator channel = badChannels.begin(); channel != badChannels.end();
         channel++) {
      LogTrace(metname) << "[phiEfficiency test] bad ranges: " << (*channel).getBin()
                        << "  Contents : " << (*channel).getContents() << endl;
    }
    LogTrace(metname) << "-------- type: [phiEfficiency]  " << thePhiQReport->getMessage() << " ------- "
                      << thePhiQReport->getStatus() << endl;
  }

  //alignment plot
  string pathPos = "Muons/MuonRecoAnalyzer/muVStkSytemRotation_posMu";
  MonitorElement* muVStkSytemRotation_posMu_histo = igetter.get(pathPos);
  string pathNeg = "Muons/MuonRecoAnalyzer/muVStkSytemRotation_negMu";
  MonitorElement* muVStkSytemRotation_negMu_histo = igetter.get(pathNeg);
  if (muVStkSytemRotation_posMu_histo && muVStkSytemRotation_negMu_histo) {
    TH2F* muVStkSytemRotation_posMu_root = muVStkSytemRotation_posMu_histo->getTH2F();
    TProfile* muVStkSytemRotation_posMu_profile = muVStkSytemRotation_posMu_root->ProfileX("", 1, 100);
    TH2F* muVStkSytemRotation_negMu_root = muVStkSytemRotation_negMu_histo->getTH2F();
    TProfile* muVStkSytemRotation_negMu_profile = muVStkSytemRotation_negMu_root->ProfileX("", 1, 100);

    for (int x = 1; x < 50; x++) {
      globalRotation[0]->Fill((x * 4) - 1, muVStkSytemRotation_posMu_profile->GetBinContent(x));
      globalRotation[0]->setBinError(x, muVStkSytemRotation_posMu_profile->GetBinError(x));
      globalRotation[1]->Fill((x * 4) - 1, muVStkSytemRotation_negMu_profile->GetBinContent(x));
      globalRotation[1]->setBinError(x, muVStkSytemRotation_negMu_profile->GetBinError(x));
      globalRotation[2]->Fill(
          (x * 4) - 1,
          muVStkSytemRotation_posMu_profile->GetBinContent(x) - muVStkSytemRotation_negMu_profile->GetBinContent(x));
      globalRotation[2]->setBinError(
          x,
          sqrt(
              (muVStkSytemRotation_posMu_profile->GetBinError(x) * muVStkSytemRotation_posMu_profile->GetBinError(x)) +
              (muVStkSytemRotation_negMu_profile->GetBinError(x) * muVStkSytemRotation_negMu_profile->GetBinError(x))));
    }
  }
}