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
#include "DataFormats/HLTReco/interface/TriggerObject.h"
#include "FWCore/Common/interface/TriggerNames.h"
#include "FWCore/Framework/interface/Frameworkfwd.h"
#include "FWCore/Framework/interface/MakerMacros.h"
#include "FWCore/MessageLogger/interface/MessageLogger.h"
#include "HLTriggerOffline/SUSYBSM/interface/SUSY_HLT_ElecFakes.h"
#include <iostream>

SUSY_HLT_ElecFakes::SUSY_HLT_ElecFakes(const edm::ParameterSet &ps) {
  edm::LogInfo("SUSY_HLT_ElecFakes") << "Constructor SUSY_HLT_ElecFakes::SUSY_HLT_ElecFakes " << std::endl;
  // Get parameters from configuration file
  theTrigSummary_ = consumes<trigger::TriggerEvent>(ps.getParameter<edm::InputTag>("trigSummary"));
  triggerResults_ = consumes<edm::TriggerResults>(ps.getParameter<edm::InputTag>("TriggerResults"));
  HLTProcess_ = ps.getParameter<std::string>("HLTProcess");
  triggerPath_ = ps.getParameter<std::string>("TriggerPath");
  triggerFilter_ = ps.getParameter<edm::InputTag>("TriggerFilter");
  triggerJetFilter_ = ps.getParameter<edm::InputTag>("TriggerJetFilter");
}

SUSY_HLT_ElecFakes::~SUSY_HLT_ElecFakes() {
  edm::LogInfo("SUSY_HLT_ElecFakes") << "Destructor SUSY_HLT_ElecFakes::~SUSY_HLT_ElecFakes " << std::endl;
}

void SUSY_HLT_ElecFakes::dqmBeginRun(edm::Run const &run, edm::EventSetup const &e) {
  bool changed;

  if (!fHltConfig.init(run, e, HLTProcess_, changed)) {
    edm::LogError("SUSY_HLT_ElecFakes") << "Initialization of HLTConfigProvider failed!!";
    return;
  }

  bool pathFound = false;
  const std::vector<std::string> allTrigNames = fHltConfig.triggerNames();
  for (size_t j = 0; j < allTrigNames.size(); ++j) {
    if (allTrigNames[j].find(triggerPath_) != std::string::npos) {
      pathFound = true;
    }
  }

  if (!pathFound) {
    edm::LogInfo("SUSY_HLT_ElecFakes") << "Path not found"
                                       << "\n";
    return;
  }

  edm::LogInfo("SUSY_HLT_ElecFakes") << "SUSY_HLT_ElecFakes::beginRun" << std::endl;
}

void SUSY_HLT_ElecFakes::bookHistograms(DQMStore::IBooker &ibooker_, edm::Run const &, edm::EventSetup const &) {
  edm::LogInfo("SUSY_HLT_ElecFakes") << "SUSY_HLT_ElecFakes::bookHistograms" << std::endl;
  // book at beginRun
  bookHistos(ibooker_);
}

void SUSY_HLT_ElecFakes::analyze(edm::Event const &e, edm::EventSetup const &eSetup) {
  edm::LogInfo("SUSY_HLT_ElecFakes") << "SUSY_HLT_ElecFakes::analyze" << std::endl;

  //-------------------------------
  //--- Trigger
  //-------------------------------
  edm::Handle<edm::TriggerResults> hltresults;
  e.getByToken(triggerResults_, hltresults);
  if (!hltresults.isValid()) {
    edm::LogError("SUSY_HLT_ElecFakes") << "invalid collection: TriggerResults"
                                        << "\n";
    return;
  }
  edm::Handle<trigger::TriggerEvent> triggerSummary;
  e.getByToken(theTrigSummary_, triggerSummary);
  if (!triggerSummary.isValid()) {
    edm::LogError("SUSY_HLT_ElecFakes") << "invalid collection: TriggerSummary"
                                        << "\n";
    return;
  }

  // get online objects
  size_t filterIndex = triggerSummary->filterIndex(triggerFilter_);
  trigger::TriggerObjectCollection triggerObjects = triggerSummary->getObjects();
  if (!(filterIndex >= triggerSummary->sizeFilters())) {
    const trigger::Keys &keys = triggerSummary->filterKeys(filterIndex);
    for (size_t j = 0; j < keys.size(); ++j) {
      trigger::TriggerObject foundObject = triggerObjects[keys[j]];
      //      if(foundObject.id() == 11){ //Electrons check number
      h_triggerElPt->Fill(foundObject.pt());
      h_triggerElEta->Fill(foundObject.eta());
      h_triggerElPhi->Fill(foundObject.phi());
      //      }
    }
  }

  filterIndex = triggerSummary->filterIndex(triggerJetFilter_);
  //  triggerObjects = triggerSummary->getObjects();
  if (!(filterIndex >= triggerSummary->sizeFilters())) {
    const trigger::Keys &keys = triggerSummary->filterKeys(filterIndex);
    for (size_t j = 0; j < keys.size(); ++j) {
      trigger::TriggerObject foundObject = triggerObjects[keys[j]];
      h_triggerJetPt->Fill(foundObject.pt());
      h_triggerJetEta->Fill(foundObject.eta());
      h_triggerJetPhi->Fill(foundObject.phi());
    }
  }
}

void SUSY_HLT_ElecFakes::bookHistos(DQMStore::IBooker &ibooker_) {
  ibooker_.cd();
  ibooker_.setCurrentFolder("HLT/SUSYBSM/" + triggerPath_);

  // online quantities
  h_triggerElPt = ibooker_.book1D("triggerElPt", "Trigger El Pt; GeV", 50, 0.0, 100.0);
  h_triggerElEta = ibooker_.book1D("triggerElEta", "Trigger El Eta", 20, -2.5, 2.5);
  h_triggerElPhi = ibooker_.book1D("triggerElPhi", "Trigger El Phi", 20, -3.5, 3.5);

  h_triggerJetPt = ibooker_.book1D("triggerJetPt", "Trigger Jet Pt; GeV", 20, 0.0, 200.0);
  h_triggerJetEta = ibooker_.book1D("triggerJetEta", "Trigger Jet Eta", 20, -3.0, 3.0);
  h_triggerJetPhi = ibooker_.book1D("triggerJetPhi", "Trigger Jet Phi", 20, -3.5, 3.5);

  //  h_triggerElJetdPhi = ibooker_.book1D("triggerElJetdPhi", "Trigger El,Jet
  //  dPhi", 20, -3.5, 3.5);

  // num and den hists to be divided in harvesting step to make turn on curves
  ibooker_.cd();
}

// define this as a plug-in
DEFINE_FWK_MODULE(SUSY_HLT_ElecFakes);