CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
BtoCharmDecayVertexMerger.cc
Go to the documentation of this file.
1 // this code is partially taken from BCandidateProducer of L. Wehrli
2 
3 // first revised prototype version for CMSSW 29.Aug.2012
4 
10 
20 
22 
23 #include <algorithm>
24 
27  const std::vector<reco::CandidatePtr> & tracks = vtx.daughterPtrVector();
28 
29  for(std::vector<reco::CandidatePtr>::const_iterator track = tracks.begin(); track != tracks.end(); ++track) {
30  ROOT::Math::LorentzVector<ROOT::Math::PxPyPzM4D<double> > vec;
31  vec.SetPx((*track)->px());
32  vec.SetPy((*track)->py());
33  vec.SetPz((*track)->pz());
35  sum += vec;
36  }
37  return sum;
38 }
39 
40 template<typename VTX>
42  public:
44 
45  virtual void produce(edm::Event &event, const edm::EventSetup &es);
46 
48 
49  private:
53  // double maxFraction;
54  // double minSignificance;
55 
57 
58  // a vertex proxy for sorting using stl
59  struct VertexProxy{
60  VTX vert;
61  double invm;
62  size_t ntracks;
63  };
64 
65  // comparison operator for VertexProxy, used in sorting
66  friend bool operator<(VertexProxy v1, VertexProxy v2){
67  if(v1.ntracks>2 && v2.ntracks<3) return true;
68  if(v1.ntracks<3 && v2.ntracks>2) return false;
69  return (v1.invm>v2.invm);
70  }
71 
72  VertexProxy buildVertexProxy(const VTX & vtx);
73 
74  void resolveBtoDchain(std::vector<VertexProxy>& coll, unsigned int i, unsigned int k);
78 };
79 
80 
81 template<typename VTX>
83  token_primaryVertex(consumes<reco::VertexCollection>(params.getParameter<edm::InputTag>("primaryVertices"))),
84  token_secondaryVertex(consumes<edm::View<VTX> >(params.getParameter<edm::InputTag>("secondaryVertices"))),
85  maxDRForUnique(params.getParameter<double>("maxDRUnique")),
86  maxvecSumIMCUTForUnique(params.getParameter<double>("maxvecSumIMifsmallDRUnique")),
87  minCosPAtomerge(params.getParameter<double>("minCosPAtomerge")),
88  maxPtreltomerge(params.getParameter<double>("maxPtreltomerge"))
89  // maxFraction(params.getParameter<double>("maxFraction")),
90  // minSignificance(params.getParameter<double>("minSignificance"))
91 {
92  produces<std::vector<VTX> >();
93 }
94 //------------------------------------
95 template<typename VTX>
97 
98  using namespace reco;
99  // PV
101  iEvent.getByToken(token_primaryVertex, PVcoll);
102 
103  if(PVcoll->size()!=0) {
104 
105  const reco::VertexCollection pvc = *( PVcoll.product());
106  pv = pvc[0];
107 
108  // get the IVF collection
110  iEvent.getByToken(token_secondaryVertex, secondaryVertices);
111 
112 
113 
114  //loop over vertices, fill into own collection for sorting
115  std::vector<VertexProxy> vertexProxyColl;
116  for(typename edm::View<VTX>::const_iterator sv = secondaryVertices->begin();
117  sv != secondaryVertices->end(); ++sv) {
118  vertexProxyColl.push_back( buildVertexProxy(*sv) );
119  }
120 
121  // sort the vertices by mass and track multiplicity
122  sort( vertexProxyColl.begin(), vertexProxyColl.end());
123 
124 
125  // loop forward over all vertices
126  for(unsigned int iVtx=0; iVtx < vertexProxyColl.size(); ++iVtx){
127 
128  // nested loop backwards (in order to start from light masses)
129  // check all vertices against each other for B->D chain
130  for(unsigned int kVtx=vertexProxyColl.size()-1; kVtx>iVtx; --kVtx){
131  // remove D vertices from the collection and add the tracks to the original one
132  resolveBtoDchain(vertexProxyColl, iVtx, kVtx);
133  }
134  }
135 
136  // now create new vertex collection and add to event
137  std::auto_ptr<std::vector<VTX> > bvertColl(new std::vector<VTX>);
138  for(typename std::vector<VertexProxy>::const_iterator it=vertexProxyColl.begin(); it!=vertexProxyColl.end(); ++it) bvertColl->push_back((*it).vert);
139  iEvent.put(bvertColl);
140  }
141  else{
142  std::auto_ptr<std::vector<VTX> > bvertCollEmpty(new std::vector<VTX>);
143  iEvent.put(bvertCollEmpty);
144  }
145 }
146 //------------------------------------
147 template<typename VTX>
150  return GlobalVector(v2.x() - v1.x(),
151  v2.y() - v1.y(),
152  v2.z() - v1.z());
153 }
154 
155 template<typename VTX>
158  return GlobalVector(v2.vertex().x() - v1.x(),
159  v2.vertex().y() - v1.y(),
160  v2.vertex().z() - v1.z());
161 }
162 
163 template<typename VTX>
166  return GlobalVector(v2.vertex().x() - v1.vertex().x(),
167  v2.vertex().y() - v1.vertex().y(),
168  v2.vertex().z() - v1.vertex().z());
169 }
170 //-------------- template specializations --------------------
171 
172 // -------------- buildVertexProxy ----------------
173 template<>
176  VertexProxy vtxProxy = {vtx,vtx.p4().M(),vtx.tracksSize()};
177  return vtxProxy;
178 }
179 
180 template<>
183  VertexProxy vtxProxy = {vtx,vtxP4(vtx).M(),vtx.numberOfSourceCandidatePtrs()};
184  return vtxProxy;
185 }
186 
187 // -------------- resolveBtoDchain ----------------
188 template<>
189 void BtoCharmDecayVertexMergerT<reco::Vertex>::resolveBtoDchain(std::vector<VertexProxy>& coll, unsigned int i, unsigned int k){
190  using namespace reco;
191 
192  GlobalVector momentum1 = GlobalVector(coll[i].vert.p4().X(), coll[i].vert.p4().Y(), coll[i].vert.p4().Z());
193  GlobalVector momentum2 = GlobalVector(coll[k].vert.p4().X(), coll[k].vert.p4().Y(), coll[k].vert.p4().Z());
194 
195  SecondaryVertex sv1(pv, coll[i].vert, momentum1 , true);
196  SecondaryVertex sv2(pv, coll[k].vert, momentum2 , true);
197 
198 
199  // find out which one is near and far
200  SecondaryVertex svNear = sv1;
201  SecondaryVertex svFar = sv2;
202  GlobalVector momentumNear = momentum1;
203  GlobalVector momentumFar = momentum2;
204 
205  // swap if it is the other way around
206  if(sv1.dist3d().value() >= sv2.dist3d().value()){
207  svNear = sv2;
208  svFar = sv1;
209  momentumNear = momentum2;
210  momentumFar = momentum1;
211  }
212  GlobalVector nearToFar = flightDirection( svNear, svFar);
213  GlobalVector pvToNear = flightDirection( pv, svNear);
214  GlobalVector pvToFar = flightDirection( pv, svFar);
215 
216  double cosPA = nearToFar.dot(momentumFar) / momentumFar.mag()/ nearToFar.mag();
217  double cosa = pvToNear. dot(momentumFar) / pvToNear.mag() / momentumFar.mag();
218  double ptrel = sqrt(1.0 - cosa*cosa)* momentumFar.mag();
219 
220  double vertexDeltaR = Geom::deltaR(pvToNear, pvToFar);
221 
222  // create a set of all tracks from both vertices, avoid double counting by using a std::set<>
223  std::set<reco::TrackRef> trackrefs;
224  // first vertex
225  for(reco::Vertex::trackRef_iterator ti = sv1.tracks_begin(); ti!=sv1.tracks_end(); ++ti){
226  if(sv1.trackWeight(*ti)>0.5){
227  reco::TrackRef t = ti->castTo<reco::TrackRef>();
228  trackrefs.insert(t);
229  }
230  }
231  // second vertex
232  for(reco::Vertex::trackRef_iterator ti = sv2.tracks_begin(); ti!=sv2.tracks_end(); ++ti){
233  if(sv2.trackWeight(*ti)>0.5){
234  reco::TrackRef t = ti->castTo<reco::TrackRef>();
235  trackrefs.insert(t);
236  }
237  }
238 
239  // now calculate one LorentzVector from the track momenta
240  ROOT::Math::LorentzVector<ROOT::Math::PxPyPzM4D<double> > mother;
241  for(std::set<reco::TrackRef>::const_iterator it = trackrefs.begin(); it!= trackrefs.end(); ++it){
242  ROOT::Math::LorentzVector<ROOT::Math::PxPyPzM4D<double> > temp ( (*it)->px(),(*it)->py(),(*it)->pz(), ParticleMasses::piPlus );
243  mother += temp;
244  }
245 
246 
247  // check if criteria for merging are fulfilled
248  if(vertexDeltaR<maxDRForUnique && mother.M()<maxvecSumIMCUTForUnique && cosPA>minCosPAtomerge && std::abs(ptrel)<maxPtreltomerge ) {
249 
250 
251  // add tracks of second vertex which are missing to first vertex
252  // loop over the second
253  bool bFoundDuplicate=false;
254  for(reco::Vertex::trackRef_iterator ti = sv2.tracks_begin(); ti!=sv2.tracks_end(); ++ti){
255  reco::Vertex::trackRef_iterator it = find(sv1.tracks_begin(), sv1.tracks_end(), *ti);
256  if (it==sv1.tracks_end()) coll[i].vert.add( *ti, sv2.refittedTrack(*ti), sv2.trackWeight(*ti) );
257  else bFoundDuplicate=true;
258  }
259  // in case a duplicate track is found, need to create the full track list in the vertex from scratch because we need to modify the weight.
260  // the weight must be the larger one, otherwise we may have outliers which are not real outliers
261 
262  if(bFoundDuplicate){
263  // create backup track containers from main vertex
264  std::vector<TrackBaseRef > tracks_;
265  std::vector<Track> refittedTracks_;
266  std::vector<float> weights_;
267  for(reco::Vertex::trackRef_iterator it = coll[i].vert.tracks_begin(); it!=coll[i].vert.tracks_end(); ++it) {
268  tracks_.push_back( *it);
269  refittedTracks_.push_back( coll[i].vert.refittedTrack(*it));
270  weights_.push_back( coll[i].vert.trackWeight(*it) );
271  }
272  // delete tracks and add all tracks back, and check in which vertex the weight is larger
273  coll[i].vert.removeTracks();
274  std::vector<Track>::iterator it2 = refittedTracks_.begin();
275  std::vector<float>::iterator it3 = weights_.begin();
276  for(reco::Vertex::trackRef_iterator it = tracks_.begin(); it!=tracks_.end(); ++it, ++it2, ++it3){
277  float weight = *it3;
278  float weight2= sv2.trackWeight(*it);
279  Track refittedTrackWithLargerWeight = *it2;
280  if( weight2 >weight) {
281  weight = weight2;
282  refittedTrackWithLargerWeight = sv2.refittedTrack(*it);
283  }
284  coll[i].vert.add(*it , refittedTrackWithLargerWeight , weight);
285  }
286  }
287 
288  // remove the second vertex from the collection
289  coll.erase( coll.begin() + k );
290  }
291 
292 }
293 
294 template<>
295 void BtoCharmDecayVertexMergerT<reco::VertexCompositePtrCandidate>::resolveBtoDchain(std::vector<VertexProxy>& coll, unsigned int i, unsigned int k){
296  using namespace reco;
297 
298  Candidate::LorentzVector vtx1P4 = vtxP4(coll[i].vert);
299  Candidate::LorentzVector vtx2P4 = vtxP4(coll[k].vert);
300  GlobalVector momentum1 = GlobalVector(vtx1P4.X(), vtx1P4.Y(), vtx1P4.Z());
301  GlobalVector momentum2 = GlobalVector(vtx2P4.X(), vtx2P4.Y(), vtx2P4.Z());
302 
303  SecondaryVertex sv1(pv, coll[i].vert, momentum1 , true);
304  SecondaryVertex sv2(pv, coll[k].vert, momentum2 , true);
305 
306 
307  // find out which one is near and far
308  SecondaryVertex svNear = sv1;
309  SecondaryVertex svFar = sv2;
310  GlobalVector momentumNear = momentum1;
311  GlobalVector momentumFar = momentum2;
312 
313  // swap if it is the other way around
314  if(sv1.dist3d().value() >= sv2.dist3d().value()){
315  svNear = sv2;
316  svFar = sv1;
317  momentumNear = momentum2;
318  momentumFar = momentum1;
319  }
320  GlobalVector nearToFar = flightDirection( svNear, svFar);
321  GlobalVector pvToNear = flightDirection( pv, svNear);
322  GlobalVector pvToFar = flightDirection( pv, svFar);
323 
324  double cosPA = nearToFar.dot(momentumFar) / momentumFar.mag()/ nearToFar.mag();
325  double cosa = pvToNear. dot(momentumFar) / pvToNear.mag() / momentumFar.mag();
326  double ptrel = sqrt(1.0 - cosa*cosa)* momentumFar.mag();
327 
328  double vertexDeltaR = Geom::deltaR(pvToNear, pvToFar);
329 
330  // create a set of all tracks from both vertices, avoid double counting by using a std::set<>
331  std::set<reco::CandidatePtr> trackrefs;
332  // first vertex
333  for(size_t i=0; i < sv1.numberOfSourceCandidatePtrs(); ++i)
334  trackrefs.insert(sv1.daughterPtr(i));
335  // second vertex
336  for(size_t i=0; i < sv2.numberOfSourceCandidatePtrs(); ++i)
337  trackrefs.insert(sv2.daughterPtr(i));
338 
339  // now calculate one LorentzVector from the track momenta
341  for(std::set<reco::CandidatePtr>::const_iterator it = trackrefs.begin(); it!= trackrefs.end(); ++it){
342  Candidate::LorentzVector temp ( (*it)->px(),(*it)->py(),(*it)->pz(), ParticleMasses::piPlus );
343  mother += temp;
344  }
345 
346  // check if criteria for merging are fulfilled
347  if(vertexDeltaR<maxDRForUnique && mother.M()<maxvecSumIMCUTForUnique && cosPA>minCosPAtomerge && std::abs(ptrel)<maxPtreltomerge ) {
348 
349  // add tracks of second vertex which are missing to first vertex
350  // loop over the second
351  const std::vector<reco::CandidatePtr> & tracks1 = sv1.daughterPtrVector();
352  const std::vector<reco::CandidatePtr> & tracks2 = sv2.daughterPtrVector();
353  for(std::vector<reco::CandidatePtr>::const_iterator ti = tracks2.begin(); ti!=tracks2.end(); ++ti){
354  std::vector<reco::CandidatePtr>::const_iterator it = find(tracks1.begin(), tracks1.end(), *ti);
355  if (it==tracks1.end()) {
356  coll[i].vert.addDaughter( *ti );
357  coll[i].vert.setP4( (*ti)->p4() + coll[i].vert.p4() );
358  }
359  }
360 
361  // remove the second vertex from the collection
362  coll.erase( coll.begin() + k );
363  }
364 }
365 
366 
367 // define specific instances of the templated BtoCharmDecayVertexMergerT class
370 
371 // define plugins
const double piPlus
Definition: ParticleMasses.h:9
tuple t
Definition: tree.py:139
int i
Definition: DBlmapReader.cc:9
reco::TemplatedSecondaryVertex< VTX > SecondaryVertex
const std::vector< reco::PFCandidatePtr > & tracks_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
virtual const Point & vertex() const
vertex position (overwritten by PF...)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
double y() const
y coordinate
Definition: Vertex.h:110
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
BtoCharmDecayVertexMergerT< reco::Vertex > BtoCharmDecayVertexMerger
VertexProxy buildVertexProxy(const VTX &vtx)
int iEvent
Definition: GenABIO.cc:230
GlobalVector flightDirection(const reco::Vertex &pv, const reco::Vertex &sv)
friend bool operator<(VertexProxy v1, VertexProxy v2)
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:115
T sqrt(T t)
Definition: SSEVec.h:48
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double z() const
y coordinate
Definition: Vertex.h:112
virtual size_type numberOfSourceCandidatePtrs() const
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
double x() const
x coordinate
Definition: Vertex.h:108
JetCorrectorParametersCollection coll
Definition: classes.h:10
tuple tracks
Definition: testEve_cfg.py:39
BtoCharmDecayVertexMergerT(const edm::ParameterSet &params)
math::XYZTLorentzVectorD p4(float mass=0.13957018, float minWeight=0.5) const
Returns the four momentum of the sum of the tracks, assuming the given mass for the decay products...
Definition: Vertex.cc:113
void resolveBtoDchain(std::vector< VertexProxy > &coll, unsigned int i, unsigned int k)
GlobalVector flightDirection(const reco::Vertex &v1, const reco::Vertex &v2)
double value() const
Definition: Measurement1D.h:28
virtual void produce(edm::Event &event, const edm::EventSetup &es)
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:37
BtoCharmDecayVertexMergerT< reco::VertexCompositePtrCandidate > CandidateBtoCharmDecayVertexMerger
reco::Candidate::LorentzVector vtxP4(const reco::VertexCompositePtrCandidate &vtx)
T dot(const Basic3DVector &v) const
Scalar product, or &quot;dot&quot; product, with a vector of same type.
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:85
edm::EDGetTokenT< edm::View< VTX > > token_secondaryVertex
const daughters & daughterPtrVector() const
references to daughtes
std::vector< TrackBaseRef >::const_iterator trackRef_iterator
The iteratator for the vector&lt;TrackRef&gt;
Definition: Vertex.h:37
int weight
Definition: histoStyle.py:50
edm::EDGetTokenT< reco::VertexCollection > token_primaryVertex
size_t tracksSize() const
number of tracks
Definition: Vertex.cc:34
Global3DVector GlobalVector
Definition: GlobalVector.h:10