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 
8 //
9 // Original Author: Matthias Geisler,32 4-B20,+41227676487,
10 // $Id: PF_PU_AssoMapAlgos.h,v 1.7 2012/11/21 09:45:04 mgeisler Exp $
11 //
12 //
13 
14 #include <string>
15 
20 
22 
26 
29 
34 
37 
40 
44 
46 
48 
49 //
50 // constants, enums and typedefs
51 //
52 
53 const double kMass = 0.49765;
54 const double lamMass = 1.11568;
55 
56 /*limits for the quality criteria*/
57 
58 const double tw_90 = 1.e-2;
59 const double tw_70 = 1.e-1;
60 const double tw_50 = 2.e-1;
61 
62 const double sec_70 = 5.;
63 const double sec_50 = 19.;
64 
65 const double fin_70 = 1.e-1;
66 const double fin_50 = 3.e-1;
67 
72 
73 typedef std::pair<reco::TrackRef, int> TrackQualityPair;
74 typedef std::vector<TrackQualityPair> TrackQualityPairVector;
75 
76 typedef std::pair<reco::VertexRef, int> VertexStepPair;
77 
78 typedef std::pair<reco::VertexRef, TrackQualityPair> VertexTrackQuality;
79 
80 typedef std::pair<reco::VertexRef, float> VertexPtsumPair;
81 typedef std::vector<VertexPtsumPair> VertexPtsumVector;
82 
84 public:
85  //dedicated constructor for the algorithms
88  // virtual destructor needed when virtual functions declared
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 and vertex-to-track maps in one go
95  std::pair<std::unique_ptr<TrackToVertexAssMap>, std::unique_ptr<VertexToTrackAssMap>> createMappings(
96  edm::Handle<reco::TrackCollection> trkcollH, const edm::EventSetup& iSetup);
97 
98  //create the track to vertex association map
99  std::unique_ptr<TrackToVertexAssMap> CreateTrackToVertexMap(edm::Handle<reco::TrackCollection>,
100  const edm::EventSetup&);
101 
102  //create the vertex to track association map
103  std::unique_ptr<VertexToTrackAssMap> CreateVertexToTrackMap(edm::Handle<reco::TrackCollection>,
104  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
120  const std::vector<reco::VertexRef>&,
122  const edm::EventSetup&,
124  int);
125 
126  //get the quality for a certain association
127  int DefineQuality(int, int, double);
128 
129 private:
130  // private methods for internal usage
131 
132  //function to find the closest vertex in z for a certain track
133  static reco::VertexRef FindClosestZ(const reco::TrackRef, const std::vector<reco::VertexRef>&, double tWeight = 0.);
134 
135  //function to find the closest vertex in 3D for a certain track
136  static reco::VertexRef FindClosest3D(reco::TransientTrack, const std::vector<reco::VertexRef>&, double tWeight = 0.);
137 
138  //function to calculate the deltaR between a vector and a vector connecting two points
139  static double dR(const math::XYZPoint&, const math::XYZVector&, edm::Handle<reco::BeamSpot>);
140 
141  //function to filter the conversion collection
142  static std::unique_ptr<reco::ConversionCollection> GetCleanedConversions(edm::Handle<reco::ConversionCollection>,
144  bool);
145 
146  //function to find out if the track comes from a gamma conversion
148 
150  const reco::Conversion&,
152  const edm::EventSetup&,
154  const std::vector<reco::VertexRef>&,
155  double);
156 
157  //function to filter the Kshort collection
158  static std::unique_ptr<reco::VertexCompositeCandidateCollection> GetCleanedKshort(
160 
161  //function to filter the Lambda collection
162  static std::unique_ptr<reco::VertexCompositeCandidateCollection> GetCleanedLambda(
164 
165  //function to find out if the track comes from a V0 decay
166  static bool ComesFromV0Decay(const reco::TrackRef,
170 
174  const edm::EventSetup&,
176  const std::vector<reco::VertexRef>&,
177  double);
178 
179  //function to filter the nuclear interaction collection
180  static std::unique_ptr<reco::PFDisplacedVertexCollection> GetCleanedNI(edm::Handle<reco::PFDisplacedVertexCollection>,
182  bool);
183 
184  //function to find out if the track comes from a nuclear interaction
186 
190  const edm::EventSetup&,
192  const std::vector<reco::VertexRef>&,
193  double);
194 
195  //function to find the vertex with the highest TrackWeight for a certain track
196  template <typename TREF>
197  static reco::VertexRef TrackWeightAssociation(const TREF&, const std::vector<reco::VertexRef>&);
198 
199  // ----------member data ---------------------------
200 
202 
205 
208 
210 
213 
216  std::unique_ptr<reco::ConversionCollection> cleanedConvCollP;
217 
220  std::unique_ptr<reco::VertexCompositeCandidateCollection> cleanedKshortCollP;
221 
224  std::unique_ptr<reco::VertexCompositeCandidateCollection> cleanedLambdaCollP;
225 
228  std::unique_ptr<reco::PFDisplacedVertexCollection> cleanedNICollP;
229 
231 
233  bool missingColls; // is true if there is a diplaced vertex collection in the event
234 
236 
237  int maxNumWarnings_; // CV: print Warning if TrackExtra objects don't exist in input file,
238  int numWarnings_; // but only a few times
239 };
240 
241 #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
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
std::vector< VertexPtsumPair > VertexPtsumVector
VertexStepPair FindAssociation(const reco::TrackRef &, const std::vector< reco::VertexRef > &, edm::ESHandle< MagneticField >, const edm::EventSetup &, edm::Handle< reco::BeamSpot >, int)
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
std::vector< TrackQualityPair > TrackQualityPairVector
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
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