CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Attributes
PFCand_AssoMapAlgos Class Reference

#include <PFCand_AssoMapAlgos.h>

Inheritance diagram for PFCand_AssoMapAlgos:
PF_PU_AssoMapAlgos PFCand_AssoMap

Public Member Functions

std::auto_ptr
< PFCandToVertexAssMap
CreatePFCandToVertexMap (edm::Handle< reco::PFCandidateCollection >, const edm::EventSetup &)
 
std::auto_ptr
< VertexToPFCandAssMap
CreateVertexToPFCandMap (edm::Handle< reco::PFCandidateCollection >, const edm::EventSetup &)
 
void GetInputCollections (edm::Event &, const edm::EventSetup &)
 
 PFCand_AssoMapAlgos (const edm::ParameterSet &)
 
std::auto_ptr
< PFCandToVertexAssMap
SortPFCandAssociationMap (PFCandToVertexAssMap *)
 
- Public Member Functions inherited from PF_PU_AssoMapAlgos
std::auto_ptr
< TrackToVertexAssMap
CreateTrackToVertexMap (edm::Handle< reco::TrackCollection >, const edm::EventSetup &)
 
std::auto_ptr
< VertexToTrackAssMap
CreateVertexToTrackMap (edm::Handle< reco::TrackCollection >, const edm::EventSetup &)
 
void GetInputCollections (edm::Event &, const edm::EventSetup &)
 
 PF_PU_AssoMapAlgos (const edm::ParameterSet &)
 
std::auto_ptr
< TrackToVertexAssMap
SortAssociationMap (TrackToVertexAssMap *)
 

Private Attributes

edm::Handle< reco::BeamSpotbeamspotH
 
edm::ESHandle< MagneticFieldbFieldH
 
edm::InputTag input_BeamSpot_
 
int input_MaxNumAssociations_
 
edm::InputTag input_VertexCollection_
 
edm::Handle
< reco::VertexCollection
vtxcollH
 

Additional Inherited Members

- Protected Member Functions inherited from PF_PU_AssoMapAlgos
std::vector< reco::VertexRef > * CreateVertexVector (edm::Handle< reco::VertexCollection >)
 
int DefineQuality (int, int, double)
 
void EraseVertex (std::vector< reco::VertexRef > *, reco::VertexRef)
 
VertexStepPair FindAssociation (const reco::TrackRef &, std::vector< reco::VertexRef > *, edm::ESHandle< MagneticField >, const edm::EventSetup &, edm::Handle< reco::BeamSpot >, int)
 

Detailed Description

Definition at line 34 of file PFCand_AssoMapAlgos.h.

Constructor & Destructor Documentation

PFCand_AssoMapAlgos::PFCand_AssoMapAlgos ( const edm::ParameterSet iConfig)

Definition at line 14 of file PFCand_AssoMapAlgos.cc.

References edm::ParameterSet::getParameter(), input_BeamSpot_, input_MaxNumAssociations_, and input_VertexCollection_.

14  :PF_PU_AssoMapAlgos(iConfig)
15 {
16 
17  input_MaxNumAssociations_ = iConfig.getParameter<int>("MaxNumberOfAssociations");
18 
19  input_VertexCollection_= iConfig.getParameter<InputTag>("VertexCollection");
20 
21  input_BeamSpot_= iConfig.getParameter<InputTag>("BeamSpot");
22 
23 }
T getParameter(std::string const &) const
edm::InputTag input_BeamSpot_
edm::InputTag input_VertexCollection_
PF_PU_AssoMapAlgos(const edm::ParameterSet &)

Member Function Documentation

std::auto_ptr< PFCandToVertexAssMap > PFCand_AssoMapAlgos::CreatePFCandToVertexMap ( edm::Handle< reco::PFCandidateCollection pfCandH,
const edm::EventSetup iSetup 
)

Definition at line 50 of file PFCand_AssoMapAlgos.cc.

References IPTools::absoluteImpactParameter3D(), beamspotH, bFieldH, PF_PU_AssoMapAlgos::CreateVertexVector(), PF_PU_AssoMapAlgos::DefineQuality(), PF_PU_AssoMapAlgos::EraseVertex(), PF_PU_AssoMapAlgos::FindAssociation(), i, input_MaxNumAssociations_, edm::Ref< C, T, F >::isNull(), edm::second(), reco::TransientTrack::setBeamSpot(), reco::TransientTrack::setES(), relval_parameters_module::step, and vtxcollH.

Referenced by PFCand_AssoMap::produce().

51 {
52 
53  auto_ptr<PFCandToVertexAssMap> pfcand2vertex(new PFCandToVertexAssMap());
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 }
VertexStepPair FindAssociation(const reco::TrackRef &, std::vector< reco::VertexRef > *, edm::ESHandle< MagneticField >, const edm::EventSetup &, edm::Handle< reco::BeamSpot >, int)
int i
Definition: DBlmapReader.cc:9
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
edm::AssociationMap< edm::OneToManyWithQuality< reco::VertexCollection, reco::PFCandidateCollection, int > > PFCandToVertexAssMap
std::pair< reco::VertexRef, int > VertexStepPair
U second(std::pair< T, U > const &p)
bool isNull() const
Checks for null.
Definition: Ref.h:247
edm::Handle< reco::BeamSpot > beamspotH
edm::ESHandle< MagneticField > bFieldH
std::vector< reco::VertexRef > * CreateVertexVector(edm::Handle< reco::VertexCollection >)
edm::Handle< reco::VertexCollection > vtxcollH
std::auto_ptr< VertexToPFCandAssMap > PFCand_AssoMapAlgos::CreateVertexToPFCandMap ( edm::Handle< reco::PFCandidateCollection pfCandH,
const edm::EventSetup iSetup 
)

Definition at line 116 of file PFCand_AssoMapAlgos.cc.

References IPTools::absoluteImpactParameter3D(), beamspotH, bFieldH, PF_PU_AssoMapAlgos::CreateVertexVector(), PF_PU_AssoMapAlgos::DefineQuality(), PF_PU_AssoMapAlgos::EraseVertex(), PF_PU_AssoMapAlgos::FindAssociation(), i, input_MaxNumAssociations_, edm::Ref< C, T, F >::isNull(), edm::second(), reco::TransientTrack::setBeamSpot(), reco::TransientTrack::setES(), relval_parameters_module::step, and vtxcollH.

Referenced by PFCand_AssoMap::produce().

117 {
118 
119  auto_ptr<VertexToPFCandAssMap> vertex2pfcand(new VertexToPFCandAssMap());
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 }
VertexStepPair FindAssociation(const reco::TrackRef &, std::vector< reco::VertexRef > *, edm::ESHandle< MagneticField >, const edm::EventSetup &, edm::Handle< reco::BeamSpot >, int)
int i
Definition: DBlmapReader.cc:9
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::pair< reco::VertexRef, int > VertexStepPair
U second(std::pair< T, U > const &p)
bool isNull() const
Checks for null.
Definition: Ref.h:247
edm::Handle< reco::BeamSpot > beamspotH
edm::ESHandle< MagneticField > bFieldH
std::vector< reco::VertexRef > * CreateVertexVector(edm::Handle< reco::VertexCollection >)
edm::Handle< reco::VertexCollection > vtxcollH
void PFCand_AssoMapAlgos::GetInputCollections ( edm::Event iEvent,
const edm::EventSetup iSetup 
)

Definition at line 30 of file PFCand_AssoMapAlgos.cc.

References beamspotH, bFieldH, edm::EventSetup::get(), edm::Event::getByLabel(), PF_PU_AssoMapAlgos::GetInputCollections(), input_BeamSpot_, input_VertexCollection_, and vtxcollH.

Referenced by PFCand_AssoMap::produce().

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 }
edm::InputTag input_BeamSpot_
edm::InputTag input_VertexCollection_
void GetInputCollections(edm::Event &, const edm::EventSetup &)
edm::Handle< reco::BeamSpot > beamspotH
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:390
edm::ESHandle< MagneticField > bFieldH
const T & get() const
Definition: EventSetup.h:55
edm::Handle< reco::VertexCollection > vtxcollH
std::auto_ptr< PFCandToVertexAssMap > PFCand_AssoMapAlgos::SortPFCandAssociationMap ( PFCandToVertexAssMap pfcvertexassInput)

Definition at line 181 of file PFCand_AssoMapAlgos.cc.

References edm::AssociationMap< Tag >::begin(), edm::AssociationMap< Tag >::end(), edm::Ref< C, T, F >::key(), and edm::second().

Referenced by PFCand_AssoMap::produce().

182 {
183  //create a new PFCandVertexAssMap for the Output which will be sorted
184  auto_ptr<PFCandToVertexAssMap> pfcvertexassOutput(new PFCandToVertexAssMap() );
185 
186  //Create and fill a vector of pairs of vertex and the summed (pT)**2 of the pfcandidates associated to the vertex
187  VertexPtsumVector vertexptsumvector;
188 
189  //loop over all vertices in the association map
190  for(PFCandToVertexAssMap::const_iterator assomap_ite=pfcvertexassInput->begin(); assomap_ite!=pfcvertexassInput->end(); assomap_ite++){
191 
192  const VertexRef assomap_vertexref = assomap_ite->key;
193  const PFCandQualityPairVector pfccoll = assomap_ite->val;
194 
195  float ptsum = 0;
196 
197  PFCandidateRef pfcandref;
198 
199  //get the pfcandidates associated to the vertex and calculate the pT**2
200  for(unsigned int pfccoll_ite=0; pfccoll_ite<pfccoll.size(); pfccoll_ite++){
201 
202  pfcandref = pfccoll[pfccoll_ite].first;
203  int quality = pfccoll[pfccoll_ite].second;
204 
205  if ( (quality<=2) && (quality!=-1) ) continue;
206 
207  double man_pT = pfcandref->pt();
208  if(man_pT>0.) ptsum+=man_pT*man_pT;
209 
210  }
211 
212  vertexptsumvector.push_back(make_pair(assomap_vertexref,ptsum));
213 
214  }
215 
216  while (vertexptsumvector.size()!=0){
217 
218  VertexRef vertexref_highestpT;
219  float highestpT = 0.;
220  int highestpT_index = 0;
221 
222  for(unsigned int vtxptsumvec_ite=0; vtxptsumvec_ite<vertexptsumvector.size(); vtxptsumvec_ite++){
223 
224  if(vertexptsumvector[vtxptsumvec_ite].second>highestpT){
225 
226  vertexref_highestpT = vertexptsumvector[vtxptsumvec_ite].first;
227  highestpT = vertexptsumvector[vtxptsumvec_ite].second;
228  highestpT_index = vtxptsumvec_ite;
229 
230  }
231 
232  }
233 
234  //loop over all vertices in the association map
235  for(PFCandToVertexAssMap::const_iterator assomap_ite=pfcvertexassInput->begin(); assomap_ite!=pfcvertexassInput->end(); assomap_ite++){
236 
237  const VertexRef assomap_vertexref = assomap_ite->key;
238  const PFCandQualityPairVector pfccoll = assomap_ite->val;
239 
240  //if the vertex from the association map the vertex with the highest pT
241  //insert all associated pfcandidates in the output Association Map
242  if(assomap_vertexref==vertexref_highestpT)
243  for(unsigned int pfccoll_ite=0; pfccoll_ite<pfccoll.size(); pfccoll_ite++)
244  pfcvertexassOutput->insert(assomap_vertexref,pfccoll[pfccoll_ite]);
245 
246  }
247 
248  vertexptsumvector.erase(vertexptsumvector.begin()+highestpT_index);
249 
250  }
251 
252  return pfcvertexassOutput;
253 }
const_iterator end() const
last iterator over the map (read only)
edm::AssociationMap< edm::OneToManyWithQuality< reco::VertexCollection, reco::PFCandidateCollection, int > > PFCandToVertexAssMap
U second(std::pair< T, U > const &p)
std::vector< VertexPtsumPair > VertexPtsumVector
std::vector< PFCandQualityPair > PFCandQualityPairVector
key_type key() const
Accessor for product key.
Definition: Ref.h:266
const_iterator begin() const
first iterator over the map (read only)

Member Data Documentation

edm::Handle<reco::BeamSpot> PFCand_AssoMapAlgos::beamspotH
private
edm::ESHandle<MagneticField> PFCand_AssoMapAlgos::bFieldH
private
edm::InputTag PFCand_AssoMapAlgos::input_BeamSpot_
private

Definition at line 64 of file PFCand_AssoMapAlgos.h.

Referenced by GetInputCollections(), and PFCand_AssoMapAlgos().

int PFCand_AssoMapAlgos::input_MaxNumAssociations_
private
edm::InputTag PFCand_AssoMapAlgos::input_VertexCollection_
private

Definition at line 61 of file PFCand_AssoMapAlgos.h.

Referenced by GetInputCollections(), and PFCand_AssoMapAlgos().

edm::Handle<reco::VertexCollection> PFCand_AssoMapAlgos::vtxcollH
private