CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PF_PU_AssoMapAlgos.h
Go to the documentation of this file.
1 #ifndef PF_PU_AssoMapAlgos_h
2 #define PF_PU_AssoMapAlgos_h
3 
4 
9 //
10 // Original Author: Matthias Geisler,32 4-B20,+41227676487,
11 // $Id: PF_PU_AssoMapAlgos.h,v 1.7 2012/11/21 09:45:04 mgeisler Exp $
12 //
13 //
14 
15 #include <string>
16 
21 
23 
27 
30 
35 
38 
41 
45 
47 
49 
50 //
51 // constants, enums and typedefs
52 //
53 
54 const double kMass = 0.49765;
55 const double lamMass = 1.11568;
56 
57 /*limits for the quality criteria*/
58 
59 const double tw_90 = 1.e-2;
60 const double tw_70 = 1.e-1;
61 const double tw_50 = 2.e-1;
62 
63 const double sec_70 = 5.;
64 const double sec_50 = 19.;
65 
66 const double fin_70 = 1.e-1;
67 const double fin_50 = 3.e-1;
68 
71 
72 typedef std::pair<reco::TrackRef, int> TrackQualityPair;
73 typedef std::vector<TrackQualityPair > TrackQualityPairVector;
74 
75 typedef std::pair<reco::VertexRef, int> VertexStepPair;
76 
77 typedef std::pair<reco::VertexRef, TrackQualityPair> VertexTrackQuality;
78 
79 typedef std::pair <reco::VertexRef, float> VertexPtsumPair;
80 typedef std::vector< VertexPtsumPair > VertexPtsumVector;
81 
83 
84  public:
85 
86  //dedicated constructor for the algorithms
88  PF_PU_AssoMapAlgos(iConfig, iC) {};
90 
91  //get all needed collections at the beginning
92  virtual void GetInputCollections(edm::Event&, const edm::EventSetup&);
93 
94  //create the track to vertex association map
95  std::auto_ptr<TrackToVertexAssMap> CreateTrackToVertexMap(edm::Handle<reco::TrackCollection>, const edm::EventSetup&);
96 
97  //create the vertex to track association map
98  std::auto_ptr<VertexToTrackAssMap> CreateVertexToTrackMap(edm::Handle<reco::TrackCollection>, const edm::EventSetup&);
99 
100  //function to sort the vertices in the AssociationMap by the sum of (pT - pT_Error)**2
101  std::unique_ptr<TrackToVertexAssMap> SortAssociationMap(TrackToVertexAssMap*);
102 
103  protected:
104  //protected functions
105 
106  //create helping vertex vector to remove associated vertices
107  std::vector<reco::VertexRef>* CreateVertexVector(edm::Handle<reco::VertexCollection>);
108 
109  //erase one vertex from the vertex vector
110  void EraseVertex(std::vector<reco::VertexRef>*, reco::VertexRef);
111 
112  //find an association for a certain track
113  VertexStepPair FindAssociation(const reco::TrackRef&, std::vector<reco::VertexRef>*,
116 
117  //get the quality for a certain association
118  int DefineQuality(int, int, double);
119 
120  private:
121 
122  // private methods for internal usage
123 
124  //function to find the closest vertex in z for a certain track
125  static reco::VertexRef FindClosestZ(const reco::TrackRef, std::vector<reco::VertexRef>*, double tWeight = 0.);
126 
127  //function to find the closest vertex in 3D for a certain track
128  static reco::VertexRef FindClosest3D(reco::TransientTrack, std::vector<reco::VertexRef>*, double tWeight = 0.);
129 
130  //function to calculate the deltaR between a vector and a vector connecting two points
131  static double dR(const math::XYZPoint&, const math::XYZVector&, edm::Handle<reco::BeamSpot>);
132 
133  //function to filter the conversion collection
134  static std::auto_ptr<reco::ConversionCollection> GetCleanedConversions(edm::Handle<reco::ConversionCollection>,
136 
137  //function to find out if the track comes from a gamma conversion
139 
142  edm::Handle<reco::BeamSpot>, std::vector<reco::VertexRef>*, double);
143 
144  //function to filter the Kshort collection
145  static std::auto_ptr<reco::VertexCompositeCandidateCollection> GetCleanedKshort(edm::Handle<reco::VertexCompositeCandidateCollection>, edm::Handle<reco::BeamSpot>, bool);
146 
147  //function to filter the Lambda collection
148  static std::auto_ptr<reco::VertexCompositeCandidateCollection> GetCleanedLambda(edm::Handle<reco::VertexCompositeCandidateCollection>, edm::Handle<reco::BeamSpot>, bool);
149 
150  //function to find out if the track comes from a V0 decay
153 
156  edm::Handle<reco::BeamSpot>, std::vector<reco::VertexRef>*, double);
157 
158  //function to filter the nuclear interaction collection
159  static std::auto_ptr<reco::PFDisplacedVertexCollection> GetCleanedNI(edm::Handle<reco::PFDisplacedVertexCollection>, edm::Handle<reco::BeamSpot>, bool);
160 
161  //function to find out if the track comes from a nuclear interaction
163 
166  edm::Handle<reco::BeamSpot>, std::vector<reco::VertexRef>*, double);
167 
168  //function to find the vertex with the highest TrackWeight for a certain track
169  static reco::VertexRef TrackWeightAssociation(const reco::TrackBaseRef&, std::vector<reco::VertexRef>*);
170 
171 
172  // ----------member data ---------------------------
173 
175 
178 
181 
183 
186 
189  std::auto_ptr<reco::ConversionCollection> cleanedConvCollP;
190 
193  std::auto_ptr<reco::VertexCompositeCandidateCollection> cleanedKshortCollP;
194 
197  std::auto_ptr<reco::VertexCompositeCandidateCollection> cleanedLambdaCollP;
198 
201  std::auto_ptr<reco::PFDisplacedVertexCollection> cleanedNICollP;
202 
204 
206  bool missingColls; // is true if there is a diplaced vertex collection in the event
207 
209 
210  int maxNumWarnings_; // CV: print Warning if TrackExtra objects don't exist in input file,
211  int numWarnings_; // but only a few times
212 
213 };
214 
215 #endif
edm::EDGetTokenT< reco::VertexCompositeCandidateCollection > LambdaCollectionToken_
VertexStepPair FindAssociation(const reco::TrackRef &, std::vector< reco::VertexRef > *, edm::ESHandle< MagneticField >, const edm::EventSetup &, edm::Handle< reco::BeamSpot >, int)
std::vector< PFDisplacedVertex > PFDisplacedVertexCollection
collection of PFDisplacedVertex objects
std::vector< TrackQualityPair > TrackQualityPairVector
edm::ESHandle< MagneticField > bFieldH
edm::EDGetTokenT< reco::BeamSpot > token_BeamSpot_
const double lamMass
std::vector< VertexCompositeCandidate > VertexCompositeCandidateCollection
collection of Candidate objects
static reco::VertexRef FindNIVertex(const reco::TrackRef, const reco::PFDisplacedVertex &, edm::ESHandle< MagneticField >, const edm::EventSetup &, edm::Handle< reco::BeamSpot >, std::vector< reco::VertexRef > *, double)
std::unique_ptr< TrackToVertexAssMap > SortAssociationMap(TrackToVertexAssMap *)
const double kMass
static bool ComesFromNI(const reco::TrackRef, const reco::PFDisplacedVertexCollection &, reco::PFDisplacedVertex *)
static reco::VertexRef FindConversionVertex(const reco::TrackRef, const reco::Conversion &, edm::ESHandle< MagneticField >, const edm::EventSetup &, edm::Handle< reco::BeamSpot >, std::vector< reco::VertexRef > *, double)
int DefineQuality(int, int, double)
static std::auto_ptr< reco::VertexCompositeCandidateCollection > GetCleanedKshort(edm::Handle< reco::VertexCompositeCandidateCollection >, edm::Handle< reco::BeamSpot >, bool)
void EraseVertex(std::vector< reco::VertexRef > *, reco::VertexRef)
static reco::VertexRef FindClosest3D(reco::TransientTrack, std::vector< reco::VertexRef > *, double tWeight=0.)
const double tw_70
std::pair< reco::VertexRef, TrackQualityPair > VertexTrackQuality
const double fin_70
std::auto_ptr< TrackToVertexAssMap > CreateTrackToVertexMap(edm::Handle< reco::TrackCollection >, const edm::EventSetup &)
std::pair< reco::VertexRef, int > VertexStepPair
static std::auto_ptr< reco::VertexCompositeCandidateCollection > GetCleanedLambda(edm::Handle< reco::VertexCompositeCandidateCollection >, edm::Handle< reco::BeamSpot >, bool)
std::vector< Conversion > ConversionCollection
collectin of Conversion objects
Definition: ConversionFwd.h:9
static std::auto_ptr< reco::PFDisplacedVertexCollection > GetCleanedNI(edm::Handle< reco::PFDisplacedVertexCollection >, edm::Handle< reco::BeamSpot >, bool)
static reco::VertexRef FindV0Vertex(const reco::TrackRef, const reco::VertexCompositeCandidate &, edm::ESHandle< MagneticField >, const edm::EventSetup &, edm::Handle< reco::BeamSpot >, std::vector< reco::VertexRef > *, double)
std::vector< VertexPtsumPair > VertexPtsumVector
static reco::VertexRef FindClosestZ(const reco::TrackRef, std::vector< reco::VertexRef > *, double tWeight=0.)
PF_PU_AssoMapAlgos(const edm::ParameterSet &iConfig, edm::ConsumesCollector &&iC)
edm::Handle< reco::PFDisplacedVertexCollection > displVertexCollH
virtual void GetInputCollections(edm::Event &, const edm::EventSetup &)
std::auto_ptr< reco::VertexCompositeCandidateCollection > cleanedKshortCollP
static reco::VertexRef TrackWeightAssociation(const reco::TrackBaseRef &, std::vector< reco::VertexRef > *)
static bool ComesFromV0Decay(const reco::TrackRef, const reco::VertexCompositeCandidateCollection &, const reco::VertexCompositeCandidateCollection &, reco::VertexCompositeCandidate *)
edm::AssociationMap< edm::OneToManyWithQuality< reco::VertexCollection, reco::TrackCollection, int > > TrackToVertexAssMap
std::auto_ptr< reco::VertexCompositeCandidateCollection > cleanedLambdaCollP
static double dR(const math::XYZPoint &, const math::XYZVector &, edm::Handle< reco::BeamSpot >)
edm::EDGetTokenT< reco::PFDisplacedVertexCollection > NIVertexCollectionToken_
edm::Handle< reco::VertexCollection > vtxcollH
std::auto_ptr< reco::PFDisplacedVertexCollection > cleanedNICollP
const double sec_50
Block of elements.
edm::EDGetTokenT< reco::VertexCompositeCandidateCollection > KshortCollectionToken_
std::pair< reco::VertexRef, float > VertexPtsumPair
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
edm::EDGetTokenT< reco::ConversionCollection > ConversionsCollectionToken_
edm::Handle< reco::VertexCompositeCandidateCollection > vertCompCandCollLambdaH
edm::EDGetTokenT< reco::VertexCollection > token_VertexCollection_
std::vector< reco::VertexRef > * CreateVertexVector(edm::Handle< reco::VertexCollection >)
const double tw_50
const double sec_70
static std::auto_ptr< reco::ConversionCollection > GetCleanedConversions(edm::Handle< reco::ConversionCollection >, edm::Handle< reco::BeamSpot >, bool)
static bool ComesFromConversion(const reco::TrackRef, const reco::ConversionCollection &, reco::Conversion *)
std::pair< reco::TrackRef, int > TrackQualityPair
edm::Handle< reco::BeamSpot > beamspotH
edm::Handle< reco::ConversionCollection > convCollH
std::auto_ptr< VertexToTrackAssMap > CreateVertexToTrackMap(edm::Handle< reco::TrackCollection >, const edm::EventSetup &)
edm::AssociationMap< edm::OneToManyWithQuality< reco::TrackCollection, reco::VertexCollection, int > > VertexToTrackAssMap
const double tw_90
edm::Handle< reco::VertexCompositeCandidateCollection > vertCompCandCollKshortH
std::auto_ptr< reco::ConversionCollection > cleanedConvCollP
const double fin_50