CMS 3D CMS Logo

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  // virtual destructor needed when virtual functions declared
92 
93  //get all needed collections at the beginning
94  virtual void GetInputCollections(edm::Event&, const edm::EventSetup&);
95 
96  //create the track-to-vertex and vertex-to-track maps in one go
97  std::pair<std::unique_ptr<TrackToVertexAssMap>, std::unique_ptr<VertexToTrackAssMap>>
99 
100  //create the track to vertex association map
101  std::unique_ptr<TrackToVertexAssMap> CreateTrackToVertexMap(edm::Handle<reco::TrackCollection>, const edm::EventSetup&);
102 
103  //create the vertex to track association map
104  std::unique_ptr<VertexToTrackAssMap> CreateVertexToTrackMap(edm::Handle<reco::TrackCollection>, const edm::EventSetup&);
105 
106  //function to sort the vertices in the AssociationMap by the sum of (pT - pT_Error)**2
107  std::unique_ptr<TrackToVertexAssMap> SortAssociationMap(TrackToVertexAssMap*,edm::Handle<reco::TrackCollection>);
108 
109  protected:
110  //protected functions
111 
112  //create helping vertex vector to remove associated vertices
113  std::vector<reco::VertexRef> CreateVertexVector(edm::Handle<reco::VertexCollection>);
114 
115  //erase one vertex from the vertex vector
116  void EraseVertex(std::vector<reco::VertexRef>&, reco::VertexRef);
117 
118  //find an association for a certain track
119  VertexStepPair FindAssociation(const reco::TrackRef&, const std::vector<reco::VertexRef>&,
122 
123  //get the quality for a certain association
124  int DefineQuality(int, int, double);
125 
126  private:
127 
128  // private methods for internal usage
129 
130  //function to find the closest vertex in z for a certain track
131  static reco::VertexRef FindClosestZ(const reco::TrackRef, const std::vector<reco::VertexRef>&, double tWeight = 0.);
132 
133  //function to find the closest vertex in 3D for a certain track
134  static reco::VertexRef FindClosest3D(reco::TransientTrack, const std::vector<reco::VertexRef>&, double tWeight = 0.);
135 
136  //function to calculate the deltaR between a vector and a vector connecting two points
137  static double dR(const math::XYZPoint&, const math::XYZVector&, edm::Handle<reco::BeamSpot>);
138 
139  //function to filter the conversion collection
140  static std::unique_ptr<reco::ConversionCollection> GetCleanedConversions(edm::Handle<reco::ConversionCollection>,
142 
143  //function to find out if the track comes from a gamma conversion
145 
148  edm::Handle<reco::BeamSpot>, const std::vector<reco::VertexRef>&, double);
149 
150  //function to filter the Kshort collection
151  static std::unique_ptr<reco::VertexCompositeCandidateCollection> GetCleanedKshort(edm::Handle<reco::VertexCompositeCandidateCollection>, edm::Handle<reco::BeamSpot>, bool);
152 
153  //function to filter the Lambda collection
154  static std::unique_ptr<reco::VertexCompositeCandidateCollection> GetCleanedLambda(edm::Handle<reco::VertexCompositeCandidateCollection>, edm::Handle<reco::BeamSpot>, bool);
155 
156  //function to find out if the track comes from a V0 decay
159 
162  edm::Handle<reco::BeamSpot>, const std::vector<reco::VertexRef>&, double);
163 
164  //function to filter the nuclear interaction collection
165  static std::unique_ptr<reco::PFDisplacedVertexCollection> GetCleanedNI(edm::Handle<reco::PFDisplacedVertexCollection>, edm::Handle<reco::BeamSpot>, bool);
166 
167  //function to find out if the track comes from a nuclear interaction
169 
172  edm::Handle<reco::BeamSpot>, const std::vector<reco::VertexRef>&, double);
173 
174  //function to find the vertex with the highest TrackWeight for a certain track
175  template<typename TREF> static reco::VertexRef TrackWeightAssociation(const TREF&, const std::vector<reco::VertexRef>&);
176 
177 
178  // ----------member data ---------------------------
179 
181 
184 
187 
189 
192 
195  std::unique_ptr<reco::ConversionCollection> cleanedConvCollP;
196 
199  std::unique_ptr<reco::VertexCompositeCandidateCollection> cleanedKshortCollP;
200 
203  std::unique_ptr<reco::VertexCompositeCandidateCollection> cleanedLambdaCollP;
204 
207  std::unique_ptr<reco::PFDisplacedVertexCollection> cleanedNICollP;
208 
210 
212  bool missingColls; // is true if there is a diplaced vertex collection in the event
213 
215 
216  int maxNumWarnings_; // CV: print Warning if TrackExtra objects don't exist in input file,
217  int numWarnings_; // but only a few times
218 
219 };
220 
221 #endif
edm::EDGetTokenT< reco::VertexCompositeCandidateCollection > LambdaCollectionToken_
static std::unique_ptr< reco::VertexCompositeCandidateCollection > GetCleanedLambda(edm::Handle< reco::VertexCompositeCandidateCollection >, edm::Handle< reco::BeamSpot >, bool)
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 std::unique_ptr< reco::PFDisplacedVertexCollection > GetCleanedNI(edm::Handle< reco::PFDisplacedVertexCollection >, edm::Handle< reco::BeamSpot >, bool)
std::unique_ptr< TrackToVertexAssMap > SortAssociationMap(TrackToVertexAssMap *, edm::Handle< reco::TrackCollection >)
const double kMass
static bool ComesFromNI(const reco::TrackRef, const reco::PFDisplacedVertexCollection &, reco::PFDisplacedVertex *)
static reco::VertexRef FindV0Vertex(const reco::TrackRef, const reco::VertexCompositeCandidate &, edm::ESHandle< MagneticField >, const edm::EventSetup &, edm::Handle< reco::BeamSpot >, const std::vector< reco::VertexRef > &, double)
static std::unique_ptr< reco::VertexCompositeCandidateCollection > GetCleanedKshort(edm::Handle< reco::VertexCompositeCandidateCollection >, edm::Handle< reco::BeamSpot >, bool)
int DefineQuality(int, int, double)
static reco::VertexRef TrackWeightAssociation(const TREF &, const std::vector< reco::VertexRef > &)
const double tw_70
std::pair< reco::VertexRef, TrackQualityPair > VertexTrackQuality
const double fin_70
std::unique_ptr< TrackToVertexAssMap > CreateTrackToVertexMap(edm::Handle< reco::TrackCollection >, const edm::EventSetup &)
std::unique_ptr< reco::PFDisplacedVertexCollection > cleanedNICollP
std::unique_ptr< reco::VertexCompositeCandidateCollection > cleanedLambdaCollP
std::pair< reco::VertexRef, int > VertexStepPair
std::vector< Conversion > ConversionCollection
collectin of Conversion objects
Definition: ConversionFwd.h:9
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
static std::unique_ptr< reco::ConversionCollection > GetCleanedConversions(edm::Handle< reco::ConversionCollection >, edm::Handle< reco::BeamSpot >, bool)
PF_PU_AssoMapAlgos(const edm::ParameterSet &iConfig, edm::ConsumesCollector &&iC)
edm::Handle< reco::PFDisplacedVertexCollection > displVertexCollH
virtual void GetInputCollections(edm::Event &, const edm::EventSetup &)
std::pair< std::unique_ptr< TrackToVertexAssMap >, std::unique_ptr< VertexToTrackAssMap > > createMappings(edm::Handle< reco::TrackCollection > trkcollH, const edm::EventSetup &iSetup)
static bool ComesFromV0Decay(const reco::TrackRef, const reco::VertexCompositeCandidateCollection &, const reco::VertexCompositeCandidateCollection &, reco::VertexCompositeCandidate *)
edm::AssociationMap< edm::OneToManyWithQuality< reco::VertexCollection, reco::TrackCollection, int > > TrackToVertexAssMap
static reco::VertexRef FindConversionVertex(const reco::TrackRef, const reco::Conversion &, edm::ESHandle< MagneticField >, const edm::EventSetup &, edm::Handle< reco::BeamSpot >, const std::vector< reco::VertexRef > &, double)
static double dR(const math::XYZPoint &, const math::XYZVector &, edm::Handle< reco::BeamSpot >)
edm::EDGetTokenT< reco::PFDisplacedVertexCollection > NIVertexCollectionToken_
edm::Handle< reco::VertexCollection > vtxcollH
const double sec_50
#define noexcept
static reco::VertexRef FindClosest3D(reco::TransientTrack, const std::vector< reco::VertexRef > &, double tWeight=0.)
std::vector< reco::VertexRef > CreateVertexVector(edm::Handle< reco::VertexCollection >)
Block of elements.
std::unique_ptr< VertexToTrackAssMap > CreateVertexToTrackMap(edm::Handle< reco::TrackCollection >, const edm::EventSetup &)
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::unique_ptr< reco::VertexCompositeCandidateCollection > cleanedKshortCollP
void EraseVertex(std::vector< reco::VertexRef > &, reco::VertexRef)
const double tw_50
const double sec_70
static bool ComesFromConversion(const reco::TrackRef, const reco::ConversionCollection &, reco::Conversion *)
virtual ~PF_PU_AssoMapAlgos()(false)
std::pair< reco::TrackRef, int > TrackQualityPair
edm::Handle< reco::BeamSpot > beamspotH
edm::Handle< reco::ConversionCollection > convCollH
static reco::VertexRef FindClosestZ(const reco::TrackRef, const std::vector< reco::VertexRef > &, double tWeight=0.)
edm::AssociationMap< edm::OneToManyWithQuality< reco::TrackCollection, reco::VertexCollection, int > > VertexToTrackAssMap
const double tw_90
edm::Handle< reco::VertexCompositeCandidateCollection > vertCompCandCollKshortH
static reco::VertexRef FindNIVertex(const reco::TrackRef, const reco::PFDisplacedVertex &, edm::ESHandle< MagneticField >, const edm::EventSetup &, edm::Handle< reco::BeamSpot >, const std::vector< reco::VertexRef > &, double)
const double fin_50
std::unique_ptr< reco::ConversionCollection > cleanedConvCollP