CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
PFPileUp.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 
22 
24 
27 
28 using namespace std;
29 using namespace edm;
30 using namespace reco;
31 
43 public:
44  typedef std::vector<edm::FwdPtr<reco::PFCandidate>> PFCollection;
46  typedef std::vector<reco::PFCandidate> PFCollectionByValue;
48 
49  explicit PFPileUp(const edm::ParameterSet&);
50 
51  ~PFPileUp() override;
52 
53  void produce(edm::Event&, const edm::EventSetup&) override;
54 
55 private:
57 
62 
65 
67  bool enable_;
68 
70  bool verbose_;
71 
74 
79 };
80 
82  tokenPFCandidates_ = consumes<PFCollection>(iConfig.getParameter<InputTag>("PFCandidates"));
83  tokenPFCandidatesView_ = mayConsume<PFView>(iConfig.getParameter<InputTag>("PFCandidates"));
84 
85  tokenVertices_ = consumes<VertexCollection>(iConfig.getParameter<InputTag>("Vertices"));
86 
87  fUseVertexAssociation = iConfig.getParameter<bool>("useVertexAssociation");
88  vertexAssociationQuality_ = iConfig.getParameter<int>("vertexAssociationQuality");
89  if (fUseVertexAssociation) {
90  tokenVertexAssociation_ = consumes<CandToVertex>(iConfig.getParameter<edm::InputTag>("vertexAssociation"));
91  tokenVertexAssociationQuality_ =
92  consumes<edm::ValueMap<int>>(iConfig.getParameter<edm::InputTag>("vertexAssociation"));
93  }
94 
95  enable_ = iConfig.getParameter<bool>("Enable");
96 
97  verbose_ = iConfig.getUntrackedParameter<bool>("verbose", false);
98 
99  if (iConfig.exists("checkClosestZVertex")) {
100  checkClosestZVertex_ = iConfig.getParameter<bool>("checkClosestZVertex");
101  } else {
102  checkClosestZVertex_ = false;
103  }
104 
105  // Configure the algo
106  pileUpAlgo_.setVerbose(verbose_);
107  pileUpAlgo_.setCheckClosestZVertex(checkClosestZVertex_);
108 
109  //produces<reco::PFCandidateCollection>();
110  produces<PFCollection>();
111  // produces< PFCollectionByValue > ();
112 }
113 
115 
116 void PFPileUp::produce(Event& iEvent, const EventSetup& iSetup) {
117  // LogDebug("PFPileUp")<<"START event: "<<iEvent.id().event()
118  // <<" in run "<<iEvent.id().run()<<endl;
119 
120  // get PFCandidates
121 
122  unique_ptr<PFCollection> pOutput(new PFCollection);
123 
124  unique_ptr<PFCollectionByValue> pOutputByValue(new PFCollectionByValue);
125 
126  if (enable_) {
127  // get vertices
129  iEvent.getByToken(tokenVertices_, vertices);
130 
131  // get PF Candidates
133  PFCollection const* pfCandidatesRef = nullptr;
134  PFCollection usedIfNoFwdPtrs;
135  bool getFromFwdPtr = iEvent.getByToken(tokenPFCandidates_, pfCandidates);
136  if (getFromFwdPtr) {
137  pfCandidatesRef = pfCandidates.product();
138  }
139  // Maintain backwards-compatibility.
140  // If there is no vector of FwdPtr<PFCandidate> found, then
141  // make a dummy vector<FwdPtr<PFCandidate> > for the PFPileupAlgo,
142  // set the pointer "pfCandidatesRef" to point to it, and
143  // then we can pass it to the PFPileupAlgo.
144  else {
145  Handle<PFView> pfView;
146  bool getFromView = iEvent.getByToken(tokenPFCandidatesView_, pfView);
147  if (!getFromView) {
148  throw cms::Exception(
149  "PFPileUp is misconfigured. This needs to be either vector<FwdPtr<PFCandidate> >, or View<PFCandidate>");
150  }
151  for (edm::View<reco::PFCandidate>::const_iterator viewBegin = pfView->begin(),
152  viewEnd = pfView->end(),
153  iview = viewBegin;
154  iview != viewEnd;
155  ++iview) {
156  usedIfNoFwdPtrs.push_back(
157  edm::FwdPtr<reco::PFCandidate>(pfView->ptrAt(iview - viewBegin), pfView->ptrAt(iview - viewBegin)));
158  }
159  pfCandidatesRef = &usedIfNoFwdPtrs;
160  }
161 
162  if (pfCandidatesRef == nullptr) {
163  throw cms::Exception(
164  "Something went dreadfully wrong with PFPileUp. pfCandidatesRef should never be zero, so this is a logic "
165  "error.");
166  }
167 
168  if (fUseVertexAssociation) {
169  const edm::Association<reco::VertexCollection>& associatedPV = iEvent.get(tokenVertexAssociation_);
170  const edm::ValueMap<int>& associationQuality = iEvent.get(tokenVertexAssociationQuality_);
171  PFCollection pfCandidatesFromPU;
172  for (auto& p : (*pfCandidatesRef)) {
173  const reco::VertexRef& PVOrig = associatedPV[p];
174  int quality = associationQuality[p];
175  if (PVOrig.isNonnull() && (PVOrig.key() > 0) && (quality >= vertexAssociationQuality_))
176  pfCandidatesFromPU.push_back(p);
177  }
178  pOutput->insert(pOutput->end(), pfCandidatesFromPU.begin(), pfCandidatesFromPU.end());
179  } else {
180  pileUpAlgo_.process(*pfCandidatesRef, *vertices);
181  pOutput->insert(
182  pOutput->end(), pileUpAlgo_.getPFCandidatesFromPU().begin(), pileUpAlgo_.getPFCandidatesFromPU().end());
183  }
184  // for ( PFCollection::const_iterator byValueBegin = pileUpAlgo_.getPFCandidatesFromPU().begin(),
185  // byValueEnd = pileUpAlgo_.getPFCandidatesFromPU().end(), ibyValue = byValueBegin;
186  // ibyValue != byValueEnd; ++ibyValue ) {
187  // pOutputByValue->push_back( **ibyValue );
188  // }
189 
190  } // end if enabled
191  // outsize of the loop to fill the collection anyway even when disabled
192  iEvent.put(std::move(pOutput));
193  // iEvent.put(std::move(pOutputByValue));
194 }
T getUntrackedParameter(std::string const &, T const &) const
int vertexAssociationQuality_
Definition: PFPileUp.cc:78
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:238
Identifies pile-up candidates from a collection of PFCandidates, and produces the corresponding colle...
Definition: PFPileUp.cc:42
edm::Association< reco::VertexCollection > CandToVertex
Definition: PFPileUp.cc:47
edm::EDGetTokenT< reco::VertexCollection > tokenVertices_
vertices
Definition: PFPileUp.cc:64
edm::EDGetTokenT< PFView > tokenPFCandidatesView_
fall-back token
Definition: PFPileUp.cc:61
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
bool exists(std::string const &parameterName) const
checks if a parameter exists
key_type key() const
Accessor for product key.
Definition: Ref.h:250
void produce(edm::Event &, const edm::EventSetup &) override
Definition: PFPileUp.cc:116
int iEvent
Definition: GenABIO.cc:224
bool verbose_
verbose ?
Definition: PFPileUp.cc:70
std::vector< reco::PFCandidate > PFCollectionByValue
Definition: PFPileUp.cc:46
def move
Definition: eostools.py:511
bool get(ProductID const &oid, Handle< PROD > &result) const
Definition: Event.h:346
edm::EDGetTokenT< PFCollection > tokenPFCandidates_
PFCandidates to be analyzed.
Definition: PFPileUp.cc:59
std::vector< edm::FwdPtr< reco::PFCandidate > > PFCollection
Definition: PFPileUp.cc:44
edm::EDGetTokenT< CandToVertex > tokenVertexAssociation_
Definition: PFPileUp.cc:75
bool enable_
enable PFPileUp selection
Definition: PFPileUp.cc:67
bool fUseVertexAssociation
Definition: PFPileUp.cc:77
edm::EDGetTokenT< edm::ValueMap< int > > tokenVertexAssociationQuality_
Definition: PFPileUp.cc:76
T const * product() const
Definition: Handle.h:70
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
~PFPileUp() override
Definition: PFPileUp.cc:114
PFPileUpAlgo pileUpAlgo_
Definition: PFPileUp.cc:56
bool checkClosestZVertex_
use the closest z vertex if a track is not in a vertex
Definition: PFPileUp.cc:73
PFPileUp(const edm::ParameterSet &)
Definition: PFPileUp.cc:81
edm::View< reco::PFCandidate > PFView
Definition: PFPileUp.cc:45