File indexing completed on 2024-04-06 12:01:11
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017 #include "CommonTools/RecoUtils/interface/PFCand_NoPU_WithAM.h"
0018
0019
0020 #include <memory>
0021 #include <vector>
0022
0023
0024
0025 #include "FWCore/Framework/interface/Event.h"
0026 #include "FWCore/Framework/interface/Run.h"
0027 #include "FWCore/Framework/interface/MakerMacros.h"
0028
0029 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0030
0031 #include "DataFormats/Common/interface/View.h"
0032
0033 #include "DataFormats/TrackReco/interface/Track.h"
0034 #include "DataFormats/TrackReco/interface/TrackFwd.h"
0035 #include "DataFormats/TrackReco/interface/TrackBase.h"
0036
0037 using namespace edm;
0038 using namespace std;
0039 using namespace reco;
0040
0041
0042
0043
0044 PFCand_NoPU_WithAM::PFCand_NoPU_WithAM(const edm::ParameterSet& iConfig) {
0045
0046
0047 input_AssociationType_ = iConfig.getParameter<InputTag>("AssociationType");
0048
0049 token_PFCandToVertexAssMap_ =
0050 mayConsume<PFCandToVertexAssMap>(iConfig.getParameter<InputTag>("VertexPFCandAssociationMap"));
0051 token_VertexToPFCandAssMap_ =
0052 mayConsume<VertexToPFCandAssMap>(iConfig.getParameter<InputTag>("VertexPFCandAssociationMap"));
0053
0054 token_VertexCollection_ = mayConsume<VertexCollection>(iConfig.getParameter<InputTag>("VertexCollection"));
0055
0056 input_MinQuality_ = iConfig.getParameter<int>("MinQuality");
0057
0058
0059
0060 if (input_AssociationType_.label() == "PFCandsToVertex") {
0061 produces<PFCandidateCollection>("P2V");
0062 } else {
0063 if (input_AssociationType_.label() == "VertexToPFCands") {
0064 produces<PFCandidateCollection>("V2P");
0065 } else {
0066 if (input_AssociationType_.label() == "Both") {
0067 produces<PFCandidateCollection>("P2V");
0068 produces<PFCandidateCollection>("V2P");
0069 } else {
0070 cout << "No correct InputTag for AssociationType!" << endl;
0071 cout << "Won't produce any PFCandiateCollection!" << endl;
0072 }
0073 }
0074 }
0075 }
0076
0077
0078
0079
0080
0081
0082 void PFCand_NoPU_WithAM::produce(edm::StreamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const {
0083 unique_ptr<PFCandidateCollection> p2v_firstvertex(new PFCandidateCollection());
0084 unique_ptr<PFCandidateCollection> v2p_firstvertex(new PFCandidateCollection());
0085
0086 bool p2vassmap = false;
0087 bool v2passmap = false;
0088
0089
0090 Handle<PFCandToVertexAssMap> p2vAM;
0091 Handle<VertexToPFCandAssMap> v2pAM;
0092
0093 string asstype = input_AssociationType_.label();
0094
0095 if ((asstype == "PFCandsToVertex") || (asstype == "Both")) {
0096 if (iEvent.getByToken(token_PFCandToVertexAssMap_, p2vAM)) {
0097 p2vassmap = true;
0098 }
0099 }
0100
0101 if ((asstype == "VertexToPFCands") || (asstype == "Both")) {
0102 if (iEvent.getByToken(token_VertexToPFCandAssMap_, v2pAM)) {
0103 v2passmap = true;
0104 }
0105 }
0106
0107 if (!p2vassmap && !v2passmap) {
0108 cout << "No input collection could be found" << endl;
0109 return;
0110 }
0111
0112 int negativeQuality = 0;
0113 if (input_MinQuality_ >= 2) {
0114 negativeQuality = -1;
0115 } else {
0116 if (input_MinQuality_ == 1) {
0117 negativeQuality = -2;
0118 } else {
0119 negativeQuality = -3;
0120 }
0121 }
0122
0123 if (p2vassmap) {
0124 const PFCandQualityPairVector pfccoll = p2vAM->begin()->val;
0125
0126
0127 for (unsigned int pfccoll_ite = 0; pfccoll_ite < pfccoll.size(); pfccoll_ite++) {
0128 PFCandidateRef pfcand = pfccoll[pfccoll_ite].first;
0129 int quality = pfccoll[pfccoll_ite].second;
0130
0131 if ((quality >= input_MinQuality_) || ((quality < 0) && (quality >= negativeQuality))) {
0132 p2v_firstvertex->push_back(*pfcand);
0133 }
0134 }
0135
0136 iEvent.put(std::move(p2v_firstvertex), "P2V");
0137 }
0138
0139 if (v2passmap) {
0140
0141 Handle<VertexCollection> input_vtxcollH;
0142 iEvent.getByToken(token_VertexCollection_, input_vtxcollH);
0143
0144 VertexRef firstVertexRef(input_vtxcollH, 0);
0145
0146 VertexToPFCandAssMap::const_iterator v2p_ite;
0147
0148 for (v2p_ite = v2pAM->begin(); v2p_ite != v2pAM->end(); v2p_ite++) {
0149 PFCandidateRef pfcand = v2p_ite->key;
0150
0151 for (unsigned v_ite = 0; v_ite < (v2p_ite->val).size(); v_ite++) {
0152 VertexRef vtxref = (v2p_ite->val)[v_ite].first;
0153 int quality = (v2p_ite->val)[v_ite].second;
0154
0155 if ((vtxref == firstVertexRef) &&
0156 ((quality >= input_MinQuality_) || ((quality < 0) && (quality >= negativeQuality)))) {
0157 v2p_firstvertex->push_back(*pfcand);
0158 }
0159 }
0160 }
0161
0162 iEvent.put(std::move(v2p_firstvertex), "V2P");
0163 }
0164 }
0165
0166
0167 void PFCand_NoPU_WithAM::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0168
0169
0170 edm::ParameterSetDescription desc;
0171
0172 desc.add<InputTag>("AssociationType");
0173 desc.add<InputTag>("VertexPFCandAssociationMap");
0174 desc.add<InputTag>("VertexCollection");
0175 desc.add<int>("MinQuality");
0176
0177 descriptions.addDefault(desc);
0178 }
0179
0180
0181 DEFINE_FWK_MODULE(PFCand_NoPU_WithAM);