CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
PFNoPileUpPacked.cc
Go to the documentation of this file.
1 // system include files
2 #include <memory>
3 #include <string>
4 
5 // user include files
9 
12 
16 
18 
21 
24 
25 using namespace std;
26 using namespace edm;
27 using namespace reco;
28 
39 public:
42 
43  explicit PFNoPileUpPacked(const edm::ParameterSet&);
44 
45  ~PFNoPileUpPacked() override = default;
46 
47  void produce(edm::Event&, const edm::EventSetup&) override;
48 
49  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
50 
51 private:
57 };
58 
60  tokenCandidatesView_ = consumes<CandidateView>(iConfig.getParameter<InputTag>("candidates"));
61  vertexAssociationQuality_ = iConfig.getParameter<int>("vertexAssociationQuality");
62  tokenVertexAssociation_ = consumes<CandToVertex>(iConfig.getParameter<edm::InputTag>("vertexAssociation"));
63  tokenVertexAssociationQuality_ =
64  consumes<edm::ValueMap<int>>(iConfig.getParameter<edm::InputTag>("vertexAssociation"));
65  produces<edm::PtrVector<reco::Candidate>>();
66 }
67 
69  unique_ptr<edm::PtrVector<reco::Candidate>> pOutput(new edm::PtrVector<reco::Candidate>);
70  Handle<CandidateView> candidateView;
71  iEvent.getByToken(tokenCandidatesView_, candidateView);
72  const edm::Association<reco::VertexCollection>& associatedPV = iEvent.get(tokenVertexAssociation_);
73  const edm::ValueMap<int>& associationQuality = iEvent.get(tokenVertexAssociationQuality_);
74  for (const auto& p : candidateView->ptrs()) {
75  const reco::VertexRef& PVOrig = associatedPV[p];
76  int quality = associationQuality[p];
77  if (!(PVOrig.isNonnull() && (PVOrig.key() > 0) && (quality >= vertexAssociationQuality_)))
78  pOutput->push_back(p);
79  }
80  iEvent.put(std::move(pOutput));
81 }
82 
85  desc.add<edm::InputTag>("candidates", edm::InputTag("packedPFCandidates"));
86  desc.add<int>("vertexAssociationQuality", 7);
87  desc.add<edm::InputTag>("vertexAssociation", edm::InputTag("packedPrimaryVertexAssociationJME", "original"));
88  descriptions.addWithDefaultLabel(desc);
89 }
90 
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
Identifies pile-up candidates from a collection of Candidates, and produces the corresponding collect...
void produce(edm::Event &, const edm::EventSetup &) override
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:238
edm::EDGetTokenT< reco::VertexCollection > tokenVertices_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
uint32_t const *__restrict__ Quality * quality
PFNoPileUpPacked(const edm::ParameterSet &)
edm::EDGetTokenT< CandToVertex > tokenVertexAssociation_
key_type key() const
Accessor for product key.
Definition: Ref.h:250
edm::EDGetTokenT< CandidateView > tokenCandidatesView_
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
int iEvent
Definition: GenABIO.cc:224
def move
Definition: eostools.py:511
bool get(ProductID const &oid, Handle< PROD > &result) const
Definition: Event.h:346
edm::Association< reco::VertexCollection > CandToVertex
ParameterDescriptionBase * add(U const &iLabel, T const &value)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
edm::EDGetTokenT< edm::ValueMap< int > > tokenVertexAssociationQuality_
edm::View< reco::Candidate > CandidateView
T getParameter(std::string const &) const
Definition: ParameterSet.h:303