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 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199
/** \class HLTMuonL3PreFilter
 *
 * See header file for documentation
 *
 *  \author J-R Vlimant
 *
 */

#include "HLTMuonL1toL3TkPreFilter.h"

#include "DataFormats/Common/interface/Handle.h"
#include "DataFormats/Common/interface/RefToBase.h"

#include "DataFormats/HLTReco/interface/TriggerFilterObjectWithRefs.h"
#include "DataFormats/HLTReco/interface/TriggerRefsCollections.h"

#include "FWCore/MessageLogger/interface/MessageLogger.h"
#include "DataFormats/TrackReco/interface/Track.h"
#include "DataFormats/RecoCandidate/interface/RecoChargedCandidate.h"
#include "DataFormats/RecoCandidate/interface/RecoChargedCandidateFwd.h"
#include "DataFormats/MuonReco/interface/MuonTrackLinks.h"
#include "DataFormats/BeamSpot/interface/BeamSpot.h"
#include "DataFormats/MuonSeed/interface/L3MuonTrajectorySeed.h"
#include "DataFormats/MuonSeed/interface/L3MuonTrajectorySeedCollection.h"

#include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"

//
// constructors and destructor
//
using namespace std;
using namespace edm;
using namespace reco;
using namespace trigger;

HLTMuonL1toL3TkPreFilter::HLTMuonL1toL3TkPreFilter(const ParameterSet& iConfig)
    : HLTFilter(iConfig),
      beamspotTag_(iConfig.getParameter<edm::InputTag>("BeamSpotTag")),
      beamspotToken_(consumes<reco::BeamSpot>(beamspotTag_)),
      candTag_(iConfig.getParameter<InputTag>("CandTag")),
      candToken_(consumes<reco::RecoChargedCandidateCollection>(candTag_)),
      previousCandTag_(iConfig.getParameter<InputTag>("PreviousCandTag")),
      previousCandToken_(consumes<trigger::TriggerFilterObjectWithRefs>(previousCandTag_)),
      min_N_(iConfig.getParameter<int>("MinN")),
      max_Eta_(iConfig.getParameter<double>("MaxEta")),
      min_Nhits_(iConfig.getParameter<int>("MinNhits")),
      max_Dr_(iConfig.getParameter<double>("MaxDr")),
      max_Dz_(iConfig.getParameter<double>("MaxDz")),
      min_Pt_(iConfig.getParameter<double>("MinPt")),
      nsigma_Pt_(iConfig.getParameter<double>("NSigmaPt")) {
  LogDebug("HLTMuonL1toL3TkPreFilter") << " CandTag/MinN/MaxEta/MinNhits/MaxDr/MaxDz/MinPt/NSigmaPt : "
                                       << candTag_.encode() << " " << min_N_ << " " << max_Eta_ << " " << min_Nhits_
                                       << " " << max_Dr_ << " " << max_Dz_ << " " << min_Pt_ << " " << nsigma_Pt_;

  //register your products
  produces<TriggerFilterObjectWithRefs>();
}

HLTMuonL1toL3TkPreFilter::~HLTMuonL1toL3TkPreFilter() = default;

//
// member functions
//
void HLTMuonL1toL3TkPreFilter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
  edm::ParameterSetDescription desc;
  makeHLTFilterDescription(desc);
  desc.add<edm::InputTag>("BeamSpotTag", edm::InputTag("hltBeamSpotTag"));
  desc.add<edm::InputTag>("CandTag", edm::InputTag("hltCandTag"));
  desc.add<edm::InputTag>("PreviousCandTag", edm::InputTag("hltPreviousCandTag"));
  desc.add<int>("MinN", 0);
  desc.add<double>("MaxEta", 9999.0);
  desc.add<int>("MinNhits", 0);
  desc.add<double>("MaxDr", 9999.0);
  desc.add<double>("MaxDz", 9999.0);
  desc.add<double>("MinPt", 0.0);
  desc.add<double>("NSigmaPt", 9999.0);
  descriptions.add("hltMuonL1toL3TkPreFilter", desc);
}

// ------------ method called to produce the data  ------------
bool HLTMuonL1toL3TkPreFilter::hltFilter(Event& iEvent,
                                         const EventSetup& iSetup,
                                         trigger::TriggerFilterObjectWithRefs& filterproduct) const {
  // All HLT filters must create and fill an HLT filter object,
  // recording any reconstructed physics objects satisfying (or not)
  // this HLT filter, and place it in the Event.

  // get hold of trks
  Handle<RecoChargedCandidateCollection> mucands;
  iEvent.getByToken(candToken_, mucands);
  if (saveTags())
    filterproduct.addCollectionTag(candTag_);
  // sort them by L2Track
  std::map<l1extra::L1MuonParticleRef, std::vector<RecoChargedCandidateRef> > L1toL3s;
  unsigned int n = 0;
  unsigned int maxN = mucands->size();
  for (; n != maxN; n++) {
    TrackRef tk = (*mucands)[n].track();
    edm::Ref<L3MuonTrajectorySeedCollection> l3seedRef =
        tk->seedRef().castTo<edm::Ref<L3MuonTrajectorySeedCollection> >();
    l1extra::L1MuonParticleRef l1mu = l3seedRef->l1Particle();
    L1toL3s[l1mu].push_back(RecoChargedCandidateRef(mucands, n));
  }

  // additionnal objects needed
  Handle<TriggerFilterObjectWithRefs> previousLevelCands;
  iEvent.getByToken(previousCandToken_, previousLevelCands);
  BeamSpot beamSpot;
  Handle<BeamSpot> recoBeamSpotHandle;
  iEvent.getByToken(beamspotToken_, recoBeamSpotHandle);
  beamSpot = *recoBeamSpotHandle;

  //needed to compare to L1
  vector<l1extra::L1MuonParticleRef> vl1cands;
  previousLevelCands->getObjects(TriggerL1Mu, vl1cands);

  auto L1toL3s_it = L1toL3s.begin();
  auto L1toL3s_end = L1toL3s.end();
  for (; L1toL3s_it != L1toL3s_end; ++L1toL3s_it) {
    if (!triggeredAtL1(L1toL3s_it->first, vl1cands))
      continue;

    //loop over the L3Tk reconstructed for this L1.
    unsigned int iTk = 0;
    unsigned int maxItk = L1toL3s_it->second.size();
    for (; iTk != maxItk; iTk++) {
      RecoChargedCandidateRef& cand = L1toL3s_it->second[iTk];
      TrackRef tk = cand->track();

      if (fabs(tk->eta()) > max_Eta_)
        continue;

      // cut on number of hits
      if (tk->numberOfValidHits() < min_Nhits_)
        continue;

      //dr cut
      //if (fabs(tk->d0())>max_Dr_) continue;
      if (fabs(tk->dxy(beamSpot.position())) > max_Dr_)
        continue;

      //dz cut
      if (fabs(tk->dz()) > max_Dz_)
        continue;

      // Pt threshold cut
      double pt = tk->pt();
      double err0 = tk->error(0);
      double abspar0 = fabs(tk->parameter(0));
      double ptLx = pt;
      // convert 50% efficiency threshold to 90% efficiency threshold
      if (abspar0 > 0)
        ptLx += nsigma_Pt_ * err0 / abspar0 * pt;
      LogTrace("HLTMuonL1toL3TkPreFilter") << " ...Muon in loop, pt= " << pt << ", ptLx= " << ptLx;
      if (ptLx < min_Pt_)
        continue;

      //one good L3Tk
      filterproduct.addObject(TriggerMuon, cand);
      break;  // and go on with the next L1 association
    }

  }  //loop over L1s from L3 grouping

  vector<RecoChargedCandidateRef> vref;
  filterproduct.getObjects(TriggerMuon, vref);
  for (auto& i : vref) {
    TrackRef tk = i->track();
    LogDebug("HLTMuonL1toL3TkPreFilter") << " Track passing filter: pt= " << tk->pt() << ", eta: " << tk->eta();
  }

  // filter decision
  const bool accept((int)n >= min_N_);

  LogDebug("HLTMuonL1toL3TkPreFilter") << " >>>>> Result of HLTMuonL1toL3TkPreFilter is " << accept
                                       << ", number of muons passing thresholds= " << n;

  return accept;
}

bool HLTMuonL1toL3TkPreFilter::triggeredAtL1(const l1extra::L1MuonParticleRef& l1mu,
                                             std::vector<l1extra::L1MuonParticleRef>& vcands) const {
  bool ok = false;

  // compare to previously triggered L1
  for (auto& vcand : vcands) {
    //    l1extra::L1MuonParticleRef candref =  L1MuonParticleRef(vcands[i]);
    if (vcand == l1mu) {
      ok = true;
      LogDebug("HLTMuonL1toL3TkPreFilter") << "The L1 mu triggered";
      break;
    }
  }
  return ok;
}

// declare this class as a framework plugin
#include "FWCore/Framework/interface/MakerMacros.h"
DEFINE_FWK_MODULE(HLTMuonL1toL3TkPreFilter);