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 and the inverse map */
47 /*************************************************************************************/
48 std::pair<std::unique_ptr<PFCandToVertexAssMap>, std::unique_ptr<VertexToPFCandAssMap>>
50  unique_ptr<PFCandToVertexAssMap> pfcand2vertex(new PFCandToVertexAssMap(vtxcollH, pfCandH));
51  unique_ptr<VertexToPFCandAssMap> vertex2pfcand(new VertexToPFCandAssMap(pfCandH, vtxcollH));
52 
53  int num_vertices = vtxcollH->size();
54  if ( num_vertices < input_MaxNumAssociations_) input_MaxNumAssociations_ = num_vertices;
55  vector<VertexRef> vtxColl_help;
56  if (input_MaxNumAssociations_ == 1) vtxColl_help = CreateVertexVector(vtxcollH);
57 
58  for( unsigned i=0; i<pfCandH->size(); i++ ) {
59 
60  PFCandidateRef candref(pfCandH, i);
61 
62  if (input_MaxNumAssociations_ > 1) 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  vertex2pfcand->insert( candref, make_pair(vtxColl_help.at(0), quality) );
77 
78  //cleanup only if multiple iterations are made
79  if (input_MaxNumAssociations_ > 1) PF_PU_AssoMapAlgos::EraseVertex(vtxColl_help, vtxColl_help.at(0));
80 
81  }
82 
83  } else {
84 
85  TransientTrack transtrk(PFCtrackref, &(*bFieldH) );
86  transtrk.setBeamSpot(*beamspotH);
87  transtrk.setES(iSetup);
88 
89  for ( int assoc_ite = 0; assoc_ite < input_MaxNumAssociations_; ++assoc_ite ) {
90 
91  VertexStepPair assocVtx = FindAssociation(PFCtrackref, vtxColl_help, bFieldH, iSetup, beamspotH, assoc_ite);
92  int step = assocVtx.second;
93  double distance = ( IPTools::absoluteImpactParameter3D( transtrk, *(assocVtx.first) ) ).second.value();
94 
95  int quality = DefineQuality(assoc_ite, step, distance);
96 
97  // Insert the best vertex and the pair of track and the quality of this association in the map
98  pfcand2vertex->insert( assocVtx.first, make_pair(candref, quality) );
99  vertex2pfcand->insert( candref, make_pair(assocVtx.first, quality) );
100 
101  //cleanup only if multiple iterations are made
102  if (input_MaxNumAssociations_ > 2) PF_PU_AssoMapAlgos::EraseVertex(vtxColl_help, assocVtx.first);
103 
104  }
105 
106  }//check PFCtrackref.isNull
107 
108  }//i on pfCandH
109 
110  return {std::move(pfcand2vertex), std::move(vertex2pfcand)};
111 }
112 
113 /*************************************************************************************/
114 /* create the pf candidate to vertex association map */
115 /*************************************************************************************/
116 
117 std::unique_ptr<PFCandToVertexAssMap>
119 {
120  return createMappings(pfCandH, iSetup).first;
121 }
122 
123 /*************************************************************************************/
124 /* create the vertex to pf candidate association map */
125 /*************************************************************************************/
126 
127 std::unique_ptr<VertexToPFCandAssMap>
129 {
130  return createMappings(pfCandH, iSetup).second;
131 }
132 
133 /*************************************************************************************/
134 /* create the vertex to pf candidate association map */
135 /*************************************************************************************/
136 
137 std::unique_ptr<PFCandToVertexAssMap>
139  edm::EDProductGetter const* getter)
140 {
141  //create a new PFCandVertexAssMap for the Output which will be sorted
142  unique_ptr<PFCandToVertexAssMap> pfcvertexassOutput(new PFCandToVertexAssMap(getter) );
143 
144  //Create and fill a vector of pairs of vertex and the summed (pT)**2 of the pfcandidates associated to the vertex
145  VertexPtsumVector vertexptsumvector;
146 
147  //loop over all vertices in the association map
148  for(PFCandToVertexAssMap::const_iterator assomap_ite=pfcvertexassInput->begin(); assomap_ite!=pfcvertexassInput->end(); assomap_ite++){
149 
150  const VertexRef assomap_vertexref = assomap_ite->key;
151  const PFCandQualityPairVector pfccoll = assomap_ite->val;
152 
153  float ptsum = 0;
154 
155  PFCandidateRef pfcandref;
156 
157  //get the pfcandidates associated to the vertex and calculate the pT**2
158  for(unsigned int pfccoll_ite=0; pfccoll_ite<pfccoll.size(); pfccoll_ite++){
159 
160  pfcandref = pfccoll[pfccoll_ite].first;
161  int quality = pfccoll[pfccoll_ite].second;
162 
163  if ( (quality<=2) && (quality!=-1) ) continue;
164 
165  double man_pT = pfcandref->pt();
166  if(man_pT>0.) ptsum+=man_pT*man_pT;
167 
168  }
169 
170  vertexptsumvector.push_back(make_pair(assomap_vertexref,ptsum));
171 
172  }
173 
174  while (!vertexptsumvector.empty()){
175 
176  VertexRef vertexref_highestpT;
177  float highestpT = 0.;
178  int highestpT_index = 0;
179 
180  for(unsigned int vtxptsumvec_ite=0; vtxptsumvec_ite<vertexptsumvector.size(); vtxptsumvec_ite++){
181 
182  if(vertexptsumvector[vtxptsumvec_ite].second>highestpT){
183 
184  vertexref_highestpT = vertexptsumvector[vtxptsumvec_ite].first;
185  highestpT = vertexptsumvector[vtxptsumvec_ite].second;
186  highestpT_index = vtxptsumvec_ite;
187 
188  }
189 
190  }
191 
192  //loop over all vertices in the association map
193  for(PFCandToVertexAssMap::const_iterator assomap_ite=pfcvertexassInput->begin(); assomap_ite!=pfcvertexassInput->end(); assomap_ite++){
194 
195  const VertexRef assomap_vertexref = assomap_ite->key;
196  const PFCandQualityPairVector pfccoll = assomap_ite->val;
197 
198  //if the vertex from the association map the vertex with the highest pT
199  //insert all associated pfcandidates in the output Association Map
200  if(assomap_vertexref==vertexref_highestpT)
201  for(unsigned int pfccoll_ite=0; pfccoll_ite<pfccoll.size(); pfccoll_ite++)
202  pfcvertexassOutput->insert(assomap_vertexref,pfccoll[pfccoll_ite]);
203 
204  }
205 
206  vertexptsumvector.erase(vertexptsumvector.begin()+highestpT_index);
207 
208  }
209 
210  return pfcvertexassOutput;
211 }
T getParameter(std::string const &) const
void GetInputCollections(edm::Event &, const edm::EventSetup &) override
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:517
edm::AssociationMap< edm::OneToManyWithQuality< reco::PFCandidateCollection, reco::VertexCollection, int > > VertexToPFCandAssMap
int DefineQuality(int, int, double)
std::unique_ptr< PFCandToVertexAssMap > CreatePFCandToVertexMap(edm::Handle< reco::PFCandidateCollection >, const edm::EventSetup &)
std::pair< reco::VertexRef, PFCandQualityPair > VertexPfcQuality
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:263
edm::AssociationMap< edm::OneToManyWithQuality< reco::VertexCollection, reco::PFCandidateCollection, int > > PFCandToVertexAssMap
std::pair< reco::VertexRef, int > VertexStepPair
std::pair< std::unique_ptr< PFCandToVertexAssMap >, std::unique_ptr< VertexToPFCandAssMap > > createMappings(edm::Handle< reco::PFCandidateCollection > pfCandH, const edm::EventSetup &iSetup)
U second(std::pair< T, U > const &p)
VertexStepPair FindAssociation(const reco::TrackRef &, const std::vector< reco::VertexRef > &, edm::ESHandle< MagneticField >, const edm::EventSetup &, edm::Handle< reco::BeamSpot >, int)
std::vector< VertexPtsumPair > VertexPtsumVector
int iEvent
Definition: GenABIO.cc:224
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:248
void setES(const edm::EventSetup &es)
std::vector< reco::VertexRef > CreateVertexVector(edm::Handle< reco::VertexCollection >)
edm::ESHandle< MagneticField > bFieldH
void EraseVertex(std::vector< reco::VertexRef > &, reco::VertexRef)
fixed size matrix
HLT enums.
T get() const
Definition: EventSetup.h:71
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:94
const_iterator begin() const
first iterator over the map (read only)
edm::EDGetTokenT< reco::BeamSpot > token_BeamSpot_
std::unique_ptr< VertexToPFCandAssMap > CreateVertexToPFCandMap(edm::Handle< reco::PFCandidateCollection >, const edm::EventSetup &)
def move(src, dest)
Definition: eostools.py:511