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 {
16 
17  input_MaxNumAssociations_ = iConfig.getParameter<int>("MaxNumberOfAssociations");
18 
19  token_VertexCollection_= iC.consumes<VertexCollection>(iConfig.getParameter<InputTag>("VertexCollection"));
20 
21  token_BeamSpot_= iC.consumes<BeamSpot>(iConfig.getParameter<InputTag>("BeamSpot"));
22 
23 }
24 
25 /*************************************************************************************/
26 /* get all needed collections at the beginning */
27 /*************************************************************************************/
28 
29 void
31 {
32 
34 
35  //get the offline beam spot
37 
38  //get the input vertex collection
40 
41  iSetup.get<IdealMagneticFieldRecord>().get(bFieldH);
42 
43 }
44 
45 /*************************************************************************************/
46 /* create the pf candidate to vertex association map */
47 /*************************************************************************************/
48 
49 std::auto_ptr<PFCandToVertexAssMap>
51 {
52 
53  auto_ptr<PFCandToVertexAssMap> pfcand2vertex(new PFCandToVertexAssMap(vtxcollH, pfCandH));
54 
55  int num_vertices = vtxcollH->size();
56  if ( num_vertices < input_MaxNumAssociations_) input_MaxNumAssociations_ = num_vertices;
57 
58  for( unsigned i=0; i<pfCandH->size(); i++ ) {
59 
60  PFCandidateRef candref(pfCandH, i);
61 
62  vector<VertexRef>* vtxColl_help = CreateVertexVector(vtxcollH);
63 
64  VertexPfcQuality VtxPfcQual;
65 
66  TrackRef PFCtrackref = candref->trackRef();
67 
68  if ( PFCtrackref.isNull() ){
69 
70  for ( int assoc_ite = 0; assoc_ite < input_MaxNumAssociations_; ++assoc_ite ) {
71 
72  int quality = -1 - assoc_ite;
73 
74  // Insert the best vertex and the pair of track and the quality of this association in the map
75  pfcand2vertex->insert( vtxColl_help->at(0), make_pair(candref, quality) );
76 
77  PF_PU_AssoMapAlgos::EraseVertex(vtxColl_help, vtxColl_help->at(0));
78 
79  }
80 
81  } else {
82 
83  TransientTrack transtrk(PFCtrackref, &(*bFieldH) );
84  transtrk.setBeamSpot(*beamspotH);
85  transtrk.setES(iSetup);
86 
87  for ( int assoc_ite = 0; assoc_ite < input_MaxNumAssociations_; ++assoc_ite ) {
88 
89  VertexStepPair assocVtx = FindAssociation(PFCtrackref, vtxColl_help, bFieldH, iSetup, beamspotH, assoc_ite);
90  int step = assocVtx.second;
91  double distance = ( IPTools::absoluteImpactParameter3D( transtrk, *(assocVtx.first) ) ).second.value();
92 
93  int quality = DefineQuality(assoc_ite, step, distance);
94 
95  // Insert the best vertex and the pair of track and the quality of this association in the map
96  pfcand2vertex->insert( assocVtx.first, make_pair(candref, quality) );
97 
98  PF_PU_AssoMapAlgos::EraseVertex(vtxColl_help, assocVtx.first);
99 
100  }
101 
102  }
103 
104  delete vtxColl_help;
105  }
106 
107  return pfcand2vertex;
108 
109 }
110 
111 /*************************************************************************************/
112 /* create the vertex to pf candidate association map */
113 /*************************************************************************************/
114 
115 std::auto_ptr<VertexToPFCandAssMap>
117 {
118 
119  auto_ptr<VertexToPFCandAssMap> vertex2pfcand(new VertexToPFCandAssMap(pfCandH, vtxcollH));
120 
121  int num_vertices = vtxcollH->size();
122  if ( num_vertices < input_MaxNumAssociations_) input_MaxNumAssociations_ = num_vertices;
123 
124  for( unsigned i=0; i<pfCandH->size(); i++ ) {
125 
126  PFCandidateRef candref(pfCandH, i);
127 
128  vector<VertexRef>* vtxColl_help = CreateVertexVector(vtxcollH);
129 
130  VertexPfcQuality VtxPfcQual;
131 
132  TrackRef PFCtrackref = candref->trackRef();
133 
134  if ( PFCtrackref.isNull() ){
135 
136  for ( int assoc_ite = 0; assoc_ite < input_MaxNumAssociations_; ++assoc_ite ) {
137 
138  int quality = -1 - assoc_ite;
139 
140  // Insert the best vertex and the pair of track and the quality of this association in the map
141  vertex2pfcand->insert( candref, make_pair(vtxColl_help->at(0), quality) );
142 
143  PF_PU_AssoMapAlgos::EraseVertex(vtxColl_help, vtxColl_help->at(0));
144 
145  }
146 
147  } else {
148 
149  TransientTrack transtrk(PFCtrackref, &(*bFieldH) );
150  transtrk.setBeamSpot(*beamspotH);
151  transtrk.setES(iSetup);
152 
153  for ( int assoc_ite = 0; assoc_ite < input_MaxNumAssociations_; ++assoc_ite ) {
154 
155  VertexStepPair assocVtx = FindAssociation(PFCtrackref, vtxColl_help, bFieldH, iSetup, beamspotH, assoc_ite);
156  int step = assocVtx.second;
157  double distance = ( IPTools::absoluteImpactParameter3D( transtrk, *(assocVtx.first) ) ).second.value();
158 
159  int quality = DefineQuality(assoc_ite, step, distance);
160 
161  // Insert the best vertex and the pair of track and the quality of this association in the map
162  vertex2pfcand->insert( candref, make_pair(assocVtx.first, quality) );
163 
164  PF_PU_AssoMapAlgos::EraseVertex(vtxColl_help, assocVtx.first);
165 
166  }
167 
168  }
169 
170  delete vtxColl_help;
171  }
172 
173  return vertex2pfcand;
174 }
175 
176 /*************************************************************************************/
177 /* create the vertex to pf candidate association map */
178 /*************************************************************************************/
179 
180 std::unique_ptr<PFCandToVertexAssMap>
182  edm::EDProductGetter const* getter)
183 {
184  //create a new PFCandVertexAssMap for the Output which will be sorted
185  unique_ptr<PFCandToVertexAssMap> pfcvertexassOutput(new PFCandToVertexAssMap(getter) );
186 
187  //Create and fill a vector of pairs of vertex and the summed (pT)**2 of the pfcandidates associated to the vertex
188  VertexPtsumVector vertexptsumvector;
189 
190  //loop over all vertices in the association map
191  for(PFCandToVertexAssMap::const_iterator assomap_ite=pfcvertexassInput->begin(); assomap_ite!=pfcvertexassInput->end(); assomap_ite++){
192 
193  const VertexRef assomap_vertexref = assomap_ite->key;
194  const PFCandQualityPairVector pfccoll = assomap_ite->val;
195 
196  float ptsum = 0;
197 
198  PFCandidateRef pfcandref;
199 
200  //get the pfcandidates associated to the vertex and calculate the pT**2
201  for(unsigned int pfccoll_ite=0; pfccoll_ite<pfccoll.size(); pfccoll_ite++){
202 
203  pfcandref = pfccoll[pfccoll_ite].first;
204  int quality = pfccoll[pfccoll_ite].second;
205 
206  if ( (quality<=2) && (quality!=-1) ) continue;
207 
208  double man_pT = pfcandref->pt();
209  if(man_pT>0.) ptsum+=man_pT*man_pT;
210 
211  }
212 
213  vertexptsumvector.push_back(make_pair(assomap_vertexref,ptsum));
214 
215  }
216 
217  while (vertexptsumvector.size()!=0){
218 
219  VertexRef vertexref_highestpT;
220  float highestpT = 0.;
221  int highestpT_index = 0;
222 
223  for(unsigned int vtxptsumvec_ite=0; vtxptsumvec_ite<vertexptsumvector.size(); vtxptsumvec_ite++){
224 
225  if(vertexptsumvector[vtxptsumvec_ite].second>highestpT){
226 
227  vertexref_highestpT = vertexptsumvector[vtxptsumvec_ite].first;
228  highestpT = vertexptsumvector[vtxptsumvec_ite].second;
229  highestpT_index = vtxptsumvec_ite;
230 
231  }
232 
233  }
234 
235  //loop over all vertices in the association map
236  for(PFCandToVertexAssMap::const_iterator assomap_ite=pfcvertexassInput->begin(); assomap_ite!=pfcvertexassInput->end(); assomap_ite++){
237 
238  const VertexRef assomap_vertexref = assomap_ite->key;
239  const PFCandQualityPairVector pfccoll = assomap_ite->val;
240 
241  //if the vertex from the association map the vertex with the highest pT
242  //insert all associated pfcandidates in the output Association Map
243  if(assomap_vertexref==vertexref_highestpT)
244  for(unsigned int pfccoll_ite=0; pfccoll_ite<pfccoll.size(); pfccoll_ite++)
245  pfcvertexassOutput->insert(assomap_vertexref,pfccoll[pfccoll_ite]);
246 
247  }
248 
249  vertexptsumvector.erase(vertexptsumvector.begin()+highestpT_index);
250 
251  }
252 
253  return pfcvertexassOutput;
254 }
VertexStepPair FindAssociation(const reco::TrackRef &, std::vector< reco::VertexRef > *, edm::ESHandle< MagneticField >, const edm::EventSetup &, edm::Handle< reco::BeamSpot >, int)
T getParameter(std::string const &) const
const_iterator end() const
last iterator over the map (read only)
void setBeamSpot(const reco::BeamSpot &beamSpot)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
edm::AssociationMap< edm::OneToManyWithQuality< reco::PFCandidateCollection, reco::VertexCollection, int > > VertexToPFCandAssMap
int DefineQuality(int, int, double)
std::pair< reco::VertexRef, PFCandQualityPair > VertexPfcQuality
void EraseVertex(std::vector< reco::VertexRef > *, reco::VertexRef)
std::pair< bool, Measurement1D > absoluteImpactParameter3D(const reco::TransientTrack &transientTrack, const reco::Vertex &vertex)
Definition: IPTools.cc:37
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
key_type key() const
Accessor for product key.
Definition: Ref.h:264
edm::AssociationMap< edm::OneToManyWithQuality< reco::VertexCollection, reco::PFCandidateCollection, int > > PFCandToVertexAssMap
std::pair< reco::VertexRef, int > VertexStepPair
U second(std::pair< T, U > const &p)
std::vector< VertexPtsumPair > VertexPtsumVector
int iEvent
Definition: GenABIO.cc:230
std::auto_ptr< VertexToPFCandAssMap > CreateVertexToPFCandMap(edm::Handle< reco::PFCandidateCollection >, const edm::EventSetup &)
std::auto_ptr< PFCandToVertexAssMap > CreatePFCandToVertexMap(edm::Handle< reco::PFCandidateCollection >, const edm::EventSetup &)
virtual void GetInputCollections(edm::Event &, const edm::EventSetup &)
std::vector< PFCandQualityPair > PFCandQualityPairVector
edm::Handle< reco::BeamSpot > beamspotH
bool isNull() const
Checks for null.
Definition: Ref.h:249
void setES(const edm::EventSetup &es)
void GetInputCollections(edm::Event &, const edm::EventSetup &)
edm::ESHandle< MagneticField > bFieldH
const T & get() const
Definition: EventSetup.h:56
std::vector< reco::VertexRef > * CreateVertexVector(edm::Handle< reco::VertexCollection >)
fixed size matrix
HLT enums.
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
const_iterator begin() const
first iterator over the map (read only)
edm::EDGetTokenT< reco::BeamSpot > token_BeamSpot_