CMS 3D CMS Logo

PFCand_AssoMapAlgos.cc
Go to the documentation of this file.
1 
3 
5 
6 using namespace edm;
7 using namespace std;
8 using namespace reco;
9 
10 /*************************************************************************************/
11 /* dedicated constructor for the algorithms */
12 /*************************************************************************************/
13 
15  : PF_PU_AssoMapAlgos(iConfig, iC), token_bField_(iC.esConsumes()), token_TrackingGeometry_(iC.esConsumes()) {
16  input_MaxNumAssociations_ = iConfig.getParameter<int>("MaxNumberOfAssociations");
17 
18  token_VertexCollection_ = iC.consumes<VertexCollection>(iConfig.getParameter<InputTag>("VertexCollection"));
19 
20  token_BeamSpot_ = iC.consumes<BeamSpot>(iConfig.getParameter<InputTag>("BeamSpot"));
21 }
22 
23 /*************************************************************************************/
24 /* get all needed collections at the beginning */
25 /*************************************************************************************/
26 
29 
30  //get the offline beam spot
31  iEvent.getByToken(token_BeamSpot_, beamspotH);
32 
33  //get the input vertex collection
35 
38 }
39 
40 /*************************************************************************************/
41 /* create the pf candidate to vertex association and the inverse map */
42 /*************************************************************************************/
43 std::pair<std::unique_ptr<PFCandToVertexAssMap>, std::unique_ptr<VertexToPFCandAssMap>>
45  unique_ptr<PFCandToVertexAssMap> pfcand2vertex(new PFCandToVertexAssMap(vtxcollH, pfCandH));
46  unique_ptr<VertexToPFCandAssMap> vertex2pfcand(new VertexToPFCandAssMap(pfCandH, vtxcollH));
47 
48  int num_vertices = vtxcollH->size();
49  if (num_vertices < input_MaxNumAssociations_)
50  input_MaxNumAssociations_ = num_vertices;
51  vector<VertexRef> vtxColl_help;
53  vtxColl_help = CreateVertexVector(vtxcollH);
54 
55  for (unsigned i = 0; i < pfCandH->size(); i++) {
56  PFCandidateRef candref(pfCandH, i);
57 
59  vtxColl_help = CreateVertexVector(vtxcollH);
60 
61  VertexPfcQuality VtxPfcQual;
62 
63  TrackRef PFCtrackref = candref->trackRef();
64 
65  if (PFCtrackref.isNull()) {
66  for (int assoc_ite = 0; assoc_ite < input_MaxNumAssociations_; ++assoc_ite) {
67  int quality = -1 - assoc_ite;
68 
69  // Insert the best vertex and the pair of track and the quality of this association in the map
70  pfcand2vertex->insert(vtxColl_help.at(0), make_pair(candref, quality));
71  vertex2pfcand->insert(candref, make_pair(vtxColl_help.at(0), quality));
72 
73  //cleanup only if multiple iterations are made
75  PF_PU_AssoMapAlgos::EraseVertex(vtxColl_help, vtxColl_help.at(0));
76  }
77 
78  } else {
79  TransientTrack transtrk(PFCtrackref, &(*bFieldH));
80  transtrk.setBeamSpot(*beamspotH);
82 
83  for (int assoc_ite = 0; assoc_ite < input_MaxNumAssociations_; ++assoc_ite) {
84  VertexStepPair assocVtx =
85  FindAssociation(PFCtrackref, vtxColl_help, bFieldH, trackingGeometryH, beamspotH, assoc_ite);
86  int step = assocVtx.second;
87  double distance = (IPTools::absoluteImpactParameter3D(transtrk, *(assocVtx.first))).second.value();
88 
89  int quality = DefineQuality(assoc_ite, step, distance);
90 
91  // Insert the best vertex and the pair of track and the quality of this association in the map
92  pfcand2vertex->insert(assocVtx.first, make_pair(candref, quality));
93  vertex2pfcand->insert(candref, make_pair(assocVtx.first, quality));
94 
95  //cleanup only if multiple iterations are made
97  PF_PU_AssoMapAlgos::EraseVertex(vtxColl_help, assocVtx.first);
98  }
99 
100  } //check PFCtrackref.isNull
101 
102  } //i on pfCandH
103 
104  return {std::move(pfcand2vertex), std::move(vertex2pfcand)};
105 }
106 
107 /*************************************************************************************/
108 /* create the pf candidate to vertex association map */
109 /*************************************************************************************/
110 
111 std::unique_ptr<PFCandToVertexAssMap> PFCand_AssoMapAlgos::CreatePFCandToVertexMap(
113  return createMappings(pfCandH).first;
114 }
115 
116 /*************************************************************************************/
117 /* create the vertex to pf candidate association map */
118 /*************************************************************************************/
119 
120 std::unique_ptr<VertexToPFCandAssMap> PFCand_AssoMapAlgos::CreateVertexToPFCandMap(
122  return createMappings(pfCandH).second;
123 }
124 
125 /*************************************************************************************/
126 /* create the vertex to pf candidate association map */
127 /*************************************************************************************/
128 
129 std::unique_ptr<PFCandToVertexAssMap> PFCand_AssoMapAlgos::SortPFCandAssociationMap(
130  PFCandToVertexAssMap* pfcvertexassInput, edm::EDProductGetter const* getter) {
131  //create a new PFCandVertexAssMap for the Output which will be sorted
132  unique_ptr<PFCandToVertexAssMap> pfcvertexassOutput(new PFCandToVertexAssMap(getter));
133 
134  //Create and fill a vector of pairs of vertex and the summed (pT)**2 of the pfcandidates associated to the vertex
135  VertexPtsumVector vertexptsumvector;
136 
137  //loop over all vertices in the association map
138  for (PFCandToVertexAssMap::const_iterator assomap_ite = pfcvertexassInput->begin();
139  assomap_ite != pfcvertexassInput->end();
140  assomap_ite++) {
141  const VertexRef assomap_vertexref = assomap_ite->key;
142  const PFCandQualityPairVector pfccoll = assomap_ite->val;
143 
144  float ptsum = 0;
145 
146  PFCandidateRef pfcandref;
147 
148  //get the pfcandidates associated to the vertex and calculate the pT**2
149  for (unsigned int pfccoll_ite = 0; pfccoll_ite < pfccoll.size(); pfccoll_ite++) {
150  pfcandref = pfccoll[pfccoll_ite].first;
151  int quality = pfccoll[pfccoll_ite].second;
152 
153  if ((quality <= 2) && (quality != -1))
154  continue;
155 
156  double man_pT = pfcandref->pt();
157  if (man_pT > 0.)
158  ptsum += man_pT * man_pT;
159  }
160 
161  vertexptsumvector.push_back(make_pair(assomap_vertexref, ptsum));
162  }
163 
164  while (!vertexptsumvector.empty()) {
165  VertexRef vertexref_highestpT;
166  float highestpT = 0.;
167  int highestpT_index = 0;
168 
169  for (unsigned int vtxptsumvec_ite = 0; vtxptsumvec_ite < vertexptsumvector.size(); vtxptsumvec_ite++) {
170  if (vertexptsumvector[vtxptsumvec_ite].second > highestpT) {
171  vertexref_highestpT = vertexptsumvector[vtxptsumvec_ite].first;
172  highestpT = vertexptsumvector[vtxptsumvec_ite].second;
173  highestpT_index = vtxptsumvec_ite;
174  }
175  }
176 
177  //loop over all vertices in the association map
178  for (PFCandToVertexAssMap::const_iterator assomap_ite = pfcvertexassInput->begin();
179  assomap_ite != pfcvertexassInput->end();
180  assomap_ite++) {
181  const VertexRef assomap_vertexref = assomap_ite->key;
182  const PFCandQualityPairVector pfccoll = assomap_ite->val;
183 
184  //if the vertex from the association map the vertex with the highest pT
185  //insert all associated pfcandidates in the output Association Map
186  if (assomap_vertexref == vertexref_highestpT)
187  for (unsigned int pfccoll_ite = 0; pfccoll_ite < pfccoll.size(); pfccoll_ite++)
188  pfcvertexassOutput->insert(assomap_vertexref, pfccoll[pfccoll_ite]);
189  }
190 
191  vertexptsumvector.erase(vertexptsumvector.begin() + highestpT_index);
192  }
193 
194  return pfcvertexassOutput;
195 }
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
void GetInputCollections(edm::Event &, const edm::EventSetup &) override
void setBeamSpot(const reco::BeamSpot &beamSpot)
int DefineQuality(int, int, double)
std::pair< reco::VertexRef, PFCandQualityPair > VertexPfcQuality
VertexStepPair FindAssociation(const reco::TrackRef &, const std::vector< reco::VertexRef > &, edm::ESHandle< MagneticField >, edm::ESHandle< GlobalTrackingGeometry >, edm::Handle< reco::BeamSpot >, int)
std::pair< bool, Measurement1D > absoluteImpactParameter3D(const reco::TransientTrack &transientTrack, const reco::Vertex &vertex)
Definition: IPTools.cc:38
const edm::ESGetToken< GlobalTrackingGeometry, GlobalTrackingGeometryRecord > token_TrackingGeometry_
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
edm::ESHandle< GlobalTrackingGeometry > trackingGeometryH
std::pair< reco::VertexRef, int > VertexStepPair
key_type key() const
Accessor for product key.
Definition: Ref.h:250
std::vector< VertexPtsumPair > VertexPtsumVector
const_iterator end() const
last iterator over the map (read only)
U second(std::pair< T, U > const &p)
string quality
int iEvent
Definition: GenABIO.cc:224
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > token_bField_
virtual void GetInputCollections(edm::Event &, const edm::EventSetup &)
edm::AssociationMap< edm::OneToManyWithQuality< reco::VertexCollection, reco::PFCandidateCollection, int > > PFCandToVertexAssMap
edm::AssociationMap< edm::OneToManyWithQuality< reco::PFCandidateCollection, reco::VertexCollection, int > > VertexToPFCandAssMap
bool isNull() const
Checks for null.
Definition: Ref.h:235
edm::Handle< reco::BeamSpot > beamspotH
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:130
std::vector< reco::VertexRef > CreateVertexVector(edm::Handle< reco::VertexCollection >)
std::vector< PFCandQualityPair > PFCandQualityPairVector
std::unique_ptr< PFCandToVertexAssMap > CreatePFCandToVertexMap(edm::Handle< reco::PFCandidateCollection >)
edm::ESHandle< MagneticField > bFieldH
std::unique_ptr< VertexToPFCandAssMap > CreateVertexToPFCandMap(edm::Handle< reco::PFCandidateCollection >)
void EraseVertex(std::vector< reco::VertexRef > &, reco::VertexRef)
const_iterator begin() const
first iterator over the map (read only)
fixed size matrix
HLT enums.
void setTrackingGeometry(const edm::ESHandle< GlobalTrackingGeometry > &tg)
edm::EDGetTokenT< reco::VertexCollection > token_VertexCollection_
edm::Handle< reco::VertexCollection > vtxcollH
PFCand_AssoMapAlgos(const edm::ParameterSet &, edm::ConsumesCollector &&)
std::unique_ptr< PFCandToVertexAssMap > SortPFCandAssociationMap(PFCandToVertexAssMap *, edm::EDProductGetter const *getter)
step
Definition: StallMonitor.cc:98
edm::EDGetTokenT< reco::BeamSpot > token_BeamSpot_
std::pair< std::unique_ptr< PFCandToVertexAssMap >, std::unique_ptr< VertexToPFCandAssMap > > createMappings(edm::Handle< reco::PFCandidateCollection > pfCandH)
def move(src, dest)
Definition: eostools.py:511