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 
18 
20 
21 #include <algorithm>
22 
24  public:
26 
27  virtual void produce(edm::Event &event, const edm::EventSetup &es);
28 
29  private:
33  // double maxFraction;
34  // double minSignificance;
35 
37 
38  // a vertex proxy for sorting using stl
39  struct vertexProxy{
41  double invm;
42  size_t ntracks;
43  };
44 
45  // comparison operator for vertexProxy, used in sorting
46  friend bool operator<(vertexProxy v1, vertexProxy v2){
47  if(v1.ntracks>2 && v2.ntracks<3) return true;
48  if(v1.ntracks<3 && v2.ntracks>2) return false;
49  return (v1.invm>v2.invm);
50  }
51 
52  void resolveBtoDchain(std::vector<vertexProxy>& coll, unsigned int i, unsigned int k);
54 };
55 
56 
57 
59  token_primaryVertex(consumes<reco::VertexCollection>(params.getParameter<edm::InputTag>("primaryVertices"))),
60  token_secondaryVertex(consumes<reco::VertexCollection>(params.getParameter<edm::InputTag>("secondaryVertices"))),
61  minDRForUnique(params.getUntrackedParameter<double>("minDRUnique",0.4)),
62  vecSumIMCUTForUnique(params.getUntrackedParameter<double>("minvecSumIMifsmallDRUnique",5.5)),
63  minCosPAtomerge(params.getUntrackedParameter<double>("minCosPAtomerge",0.99)),
64  maxPtreltomerge(params.getUntrackedParameter<double>("maxPtreltomerge",7777.0))
65  // maxFraction(params.getParameter<double>("maxFraction")),
66  // minSignificance(params.getParameter<double>("minSignificance"))
67 {
68  produces<reco::VertexCollection>();
69 }
70 //-----------------------
72 
73  using namespace reco;
74  // PV
76  iEvent.getByToken(token_primaryVertex, PVcoll);
77 
78  if(PVcoll->size()!=0) {
79 
80  const reco::VertexCollection pvc = *( PVcoll.product());
81  pv = pvc[0];
82 
83  // get the IVF collection
84  edm::Handle<reco::VertexCollection> secondaryVertices;
85  iEvent.getByToken(token_secondaryVertex, secondaryVertices);
86 
87 
88 
89  //loop over vertices, fill into own collection for sorting
90  std::vector<vertexProxy> vertexProxyColl;
91  for(std::vector<reco::Vertex>::const_iterator sv = secondaryVertices->begin();
92  sv != secondaryVertices->end(); ++sv) {
93  vertexProxy vtx = {*sv,(*sv).p4().M(),(*sv).tracksSize()};
94  vertexProxyColl.push_back( vtx );
95  }
96 
97  // sort the vertices by mass and track multiplicity
98  sort( vertexProxyColl.begin(), vertexProxyColl.end());
99 
100 
101  // loop forward over all vertices
102  for(unsigned int iVtx=0; iVtx < vertexProxyColl.size(); iVtx++){
103 
104  // nested loop backwards (in order to start from light masses)
105  // check all vertices against each other for B->D chain
106  for(unsigned int kVtx=vertexProxyColl.size()-1; kVtx>iVtx; kVtx--){
107  // remove D vertices from the collection and add the tracks to the original one
108  resolveBtoDchain(vertexProxyColl, iVtx, kVtx);
109  }
110  }
111 
112  // now create new vertex collection and add to event
113  VertexCollection *bvertices = new VertexCollection();
114  for(std::vector<vertexProxy>::iterator it=vertexProxyColl.begin(); it!=vertexProxyColl.end(); it++) bvertices->push_back((*it).vert);
115  std::auto_ptr<VertexCollection> bvertColl(bvertices);
116  iEvent.put(bvertColl);
117  }
118  else{
119  std::auto_ptr<VertexCollection> bvertCollEmpty(new VertexCollection);
120  iEvent.put(bvertCollEmpty);
121  }
122 }
123 
124 
125 //------------------------------------
126 void BtoCharmDecayVertexMerger::resolveBtoDchain(std::vector<vertexProxy>& coll, unsigned int i, unsigned int k){
127  using namespace reco;
128 
129  GlobalVector momentum1 = GlobalVector(coll[i].vert.p4().X(), coll[i].vert.p4().Y(), coll[i].vert.p4().Z());
130  GlobalVector momentum2 = GlobalVector(coll[k].vert.p4().X(), coll[k].vert.p4().Y(), coll[k].vert.p4().Z());
131 
132  reco::SecondaryVertex sv1(pv, coll[i].vert, momentum1 , true);
133  reco::SecondaryVertex sv2(pv, coll[k].vert, momentum2 , true);
134 
135 
136  // find out which one is near and far
137  reco::SecondaryVertex svNear = sv1;
138  reco::SecondaryVertex svFar = sv2;
139  GlobalVector momentumNear = momentum1;
140  GlobalVector momentumFar = momentum2;
141 
142  // swap if it is the other way around
143  if(sv1.dist3d().value() >= sv2.dist3d().value()){
144  svNear = sv2;
145  svFar = sv1;
146  momentumNear = momentum2;
147  momentumFar = momentum1;
148  }
149  GlobalVector nearToFar = flightDirection( svNear, svFar);
150  GlobalVector pvToNear = flightDirection( pv, svNear);
151 
152  double cosPA = nearToFar.dot(momentumFar) / momentumFar.mag()/ nearToFar.mag();
153  double cosa = pvToNear. dot(momentumFar) / pvToNear.mag() / momentumFar.mag();
154  double ptrel = sqrt(1.0 - cosa*cosa)* momentumFar.mag();
155 
156  double vertexDeltaR = Geom::deltaR(flightDirection(pv, sv1), flightDirection(pv, sv2) );
157 
158  // create a set of all tracks from both vertices, avoid double counting by using a std::set<>
159  std::set<reco::TrackRef> trackrefs;
160  // first vertex
161  for(reco::Vertex::trackRef_iterator ti = sv1.tracks_begin(); ti!=sv1.tracks_end(); ti++){
162  if(sv1.trackWeight(*ti)>0.5){
163  reco::TrackRef t = ti->castTo<reco::TrackRef>();
164  trackrefs.insert(t);
165  }
166  }
167  // second vertex
168  for(reco::Vertex::trackRef_iterator ti = sv2.tracks_begin(); ti!=sv2.tracks_end(); ti++){
169  if(sv2.trackWeight(*ti)>0.5){
170  reco::TrackRef t = ti->castTo<reco::TrackRef>();
171  trackrefs.insert(t);
172  }
173  }
174 
175  // now calculate one LorentzVector from the track momenta
176  ROOT::Math::LorentzVector<ROOT::Math::PxPyPzM4D<double> > mother;
177  for(std::set<reco::TrackRef>::const_iterator it = trackrefs.begin(); it!= trackrefs.end(); it++){
178  ROOT::Math::LorentzVector<ROOT::Math::PxPyPzM4D<double> > temp ( (*it)->px(),(*it)->py(),(*it)->pz(), 0.13957 );
179  mother += temp;
180  }
181 
182 
183  // check if criteria for merging are fulfilled
184  if(vertexDeltaR<minDRForUnique && mother.M()<vecSumIMCUTForUnique && cosPA>minCosPAtomerge && fabs(ptrel)<maxPtreltomerge ) {
185 
186 
187  // add tracks of second vertex which are missing to first vertex
188  // loop over the second
189  bool bFoundDuplicate=false;
190  for(reco::Vertex::trackRef_iterator ti = sv2.tracks_begin(); ti!=sv2.tracks_end(); ti++){
191  reco::Vertex::trackRef_iterator it = find(sv1.tracks_begin(), sv1.tracks_end(), *ti);
192  if (it==sv1.tracks_end()) coll[i].vert.add( *ti, sv2.refittedTrack(*ti), sv2.trackWeight(*ti) );
193  else bFoundDuplicate=true;
194  }
195  // 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.
196  // the weight must be the larger one, otherwise we may have outliers which are not real outliers
197 
198  if(bFoundDuplicate){
199  // create backup track containers from main vertex
200  std::vector<TrackBaseRef > tracks_;
201  std::vector<Track> refittedTracks_;
202  std::vector<float> weights_;
203  for(reco::Vertex::trackRef_iterator it = coll[i].vert.tracks_begin(); it!=coll[i].vert.tracks_end(); it++) {
204  tracks_.push_back( *it);
205  refittedTracks_.push_back( coll[i].vert.refittedTrack(*it));
206  weights_.push_back( coll[i].vert.trackWeight(*it) );
207  }
208  // delete tracks and add all tracks back, and check in which vertex the weight is larger
209  coll[i].vert.removeTracks();
210  std::vector<Track>::iterator it2 = refittedTracks_.begin();
211  std::vector<float>::iterator it3 = weights_.begin();
212  for(reco::Vertex::trackRef_iterator it = tracks_.begin(); it!=tracks_.end(); it++, it2++, it3++){
213  float weight = *it3;
214  float weight2= sv2.trackWeight(*it);
215  Track refittedTrackWithLargerWeight = *it2;
216  if( weight2 >weight) {
217  weight = weight2;
218  refittedTrackWithLargerWeight = sv2.refittedTrack(*it);
219  }
220  coll[i].vert.add(*it , refittedTrackWithLargerWeight , weight);
221  }
222  }
223 
224  // remove the second vertex from the collection
225  coll.erase( coll.begin() + k );
226  }
227 
228 }
229 //-------------
230 
231 
234  GlobalVector res(sv.position().X() - pv.position().X(),
235  sv.position().Y() - pv.position().Y(),
236  sv.position().Z() - pv.position().Z());
237  return res;
238 }
239 
240 
241 
int i
Definition: DBlmapReader.cc:9
const std::vector< reco::PFCandidatePtr > & tracks_
edm::EDGetTokenT< reco::VertexCollection > token_primaryVertex
trackRef_iterator tracks_end() const
last iterator over tracks
Definition: Vertex.cc:44
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
Measurement1D dist3d() const
Track refittedTrack(const TrackBaseRef &track) const
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
void resolveBtoDchain(std::vector< vertexProxy > &coll, unsigned int i, unsigned int k)
virtual void produce(edm::Event &event, const edm::EventSetup &es)
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
const Point & position() const
position
Definition: Vertex.h:92
int iEvent
Definition: GenABIO.cc:243
GlobalVector flightDirection(reco::Vertex &pv, reco::Vertex &sv)
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
float trackWeight(const TrackBaseRef &r) const
returns the weight with which a Track has contributed to the vertex-fit.
T sqrt(T t)
Definition: SSEVec.h:48
friend bool operator<(vertexProxy v1, vertexProxy v2)
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
int k[5][pyjets_maxn]
BtoCharmDecayVertexMerger(const edm::ParameterSet &params)
double value() const
Definition: Measurement1D.h:28
double deltaR(const Vector1 &v1, const Vector2 &v2)
Definition: VectorUtil.h:84
T dot(const Basic3DVector &v) const
Scalar product, or &quot;dot&quot; product, with a vector of same type.
std::vector< TrackBaseRef >::const_iterator trackRef_iterator
The iteratator for the vector&lt;TrackRef&gt;
Definition: Vertex.h:37
trackRef_iterator tracks_begin() const
first iterator over tracks
Definition: Vertex.cc:39
int weight
Definition: histoStyle.py:50
edm::EDGetTokenT< reco::VertexCollection > token_secondaryVertex
Global3DVector GlobalVector
Definition: GlobalVector.h:10