Skip to content

Commit

Permalink
Merge pull request #45373 from stahlleiton/HLT_UPC_CMSSW_14_1_X
Browse files Browse the repository at this point in the history
Update HLTPixelTrackFilter to store candidates
  • Loading branch information
cmsbuild authored Jul 9, 2024
2 parents db3008e + a07f4db commit 70af9f3
Showing 1 changed file with 24 additions and 11 deletions.
35 changes: 24 additions & 11 deletions HLTrigger/special/plugins/HLTPixelTrackFilter.cc
Original file line number Diff line number Diff line change
@@ -1,5 +1,5 @@
#include "FWCore/Framework/interface/Frameworkfwd.h"
#include "FWCore/Framework/interface/stream/EDFilter.h"
#include "HLTrigger/HLTcore/interface/HLTFilter.h"

#include "FWCore/Framework/interface/Event.h"
#include "FWCore/Framework/interface/EventSetup.h"
Expand All @@ -11,35 +11,37 @@
#include "FWCore/ParameterSet/interface/ParameterSet.h"
#include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
#include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
#include "DataFormats/TrackReco/interface/Track.h"
#include "DataFormats/RecoCandidate/interface/RecoChargedCandidate.h"
//
// class declaration
//

class HLTPixelTrackFilter : public edm::stream::EDFilter<> {
class HLTPixelTrackFilter : public HLTFilter {
public:
explicit HLTPixelTrackFilter(const edm::ParameterSet&);
~HLTPixelTrackFilter() override;
bool hltFilter(edm::Event&,
const edm::EventSetup&,
trigger::TriggerFilterObjectWithRefs& filterproduct) const override;
static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);

private:
bool filter(edm::Event&, const edm::EventSetup&) override;

edm::InputTag inputTag_; // input tag identifying product containing pixel clusters
unsigned int min_pixelTracks_; // minimum number of clusters
unsigned int max_pixelTracks_; // maximum number of clusters
edm::EDGetTokenT<reco::TrackCollection> inputToken_;
edm::EDGetTokenT<reco::RecoChargedCandidateCollection> inputToken_;
};

//
// constructors and destructor
//

HLTPixelTrackFilter::HLTPixelTrackFilter(const edm::ParameterSet& config)
: inputTag_(config.getParameter<edm::InputTag>("pixelTracks")),
: HLTFilter(config),
inputTag_(config.getParameter<edm::InputTag>("pixelTracks")),
min_pixelTracks_(config.getParameter<unsigned int>("minPixelTracks")),
max_pixelTracks_(config.getParameter<unsigned int>("maxPixelTracks")) {
inputToken_ = consumes<reco::TrackCollection>(inputTag_);
inputToken_ = consumes<reco::RecoChargedCandidateCollection>(inputTag_);
LogDebug("") << "Using the " << inputTag_ << " input collection";
LogDebug("") << "Requesting at least " << min_pixelTracks_ << " PixelTracks";
if (max_pixelTracks_ > 0)
Expand All @@ -50,6 +52,7 @@ HLTPixelTrackFilter::~HLTPixelTrackFilter() = default;

void HLTPixelTrackFilter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
edm::ParameterSetDescription desc;
makeHLTFilterDescription(desc);
desc.add<edm::InputTag>("pixelTracks", edm::InputTag("hltPixelTracks"));
desc.add<unsigned int>("minPixelTracks", 0);
desc.add<unsigned int>("maxPixelTracks", 0);
Expand All @@ -60,12 +63,22 @@ void HLTPixelTrackFilter::fillDescriptions(edm::ConfigurationDescriptions& descr
// member functions
//
// ------------ method called to produce the data ------------
bool HLTPixelTrackFilter::filter(edm::Event& event, const edm::EventSetup& iSetup) {
bool HLTPixelTrackFilter::hltFilter(edm::Event& event,
const edm::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.
if (saveTags())
filterproduct.addCollectionTag(inputTag_);

// get hold of products from Event
edm::Handle<reco::TrackCollection> trackColl;
event.getByToken(inputToken_, trackColl);
const auto& trackColl = event.getHandle(inputToken_);

unsigned int numTracks = trackColl->size();
for (size_t i = 0; i < numTracks; i++)
filterproduct.addObject(trigger::TriggerTrack, reco::RecoChargedCandidateRef(trackColl, i));

LogDebug("") << "Number of tracks accepted: " << numTracks;
bool accept = (numTracks >= min_pixelTracks_);
if (max_pixelTracks_ > 0)
Expand Down

0 comments on commit 70af9f3

Please sign in to comment.