CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_7/src/RecoBTag/SecondaryVertex/plugins/BtoCharmDecayVertexMerger.cc

Go to the documentation of this file.
00001 // this code is partially taken from BCandidateProducer of L. Wehrli
00002 
00003 // first revised prototype version for CMSSW 29.Aug.2012
00004 
00005 #include "FWCore/Framework/interface/EDProducer.h"
00006 #include "FWCore/Framework/interface/Event.h"
00007 #include "FWCore/Framework/interface/MakerMacros.h"
00008 #include "FWCore/Utilities/interface/InputTag.h"
00009 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00010 
00011 #include "DataFormats/Common/interface/Handle.h"
00012 #include "DataFormats/TrackReco/interface/Track.h"
00013 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00014 #include "DataFormats/VertexReco/interface/Vertex.h"
00015 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00016 #include "RecoVertex/VertexTools/interface/VertexDistance3D.h"
00017 #include "RecoBTag/SecondaryVertex/interface/SecondaryVertex.h"
00018 
00019 #include "DataFormats/GeometryVector/interface/VectorUtil.h"
00020 
00021 #include <algorithm>
00022 
00023 class BtoCharmDecayVertexMerger : public edm::EDProducer {
00024     public:
00025         BtoCharmDecayVertexMerger(const edm::ParameterSet &params);
00026 
00027         virtual void produce(edm::Event &event, const edm::EventSetup &es);
00028 
00029     private:
00030 
00031  
00032   edm::InputTag                         primaryVertexCollection;
00033   edm::InputTag                         secondaryVertexCollection;
00034   reco::Vertex pv;
00035   //    double                                  maxFraction;
00036   //    double                                  minSignificance;
00037 
00038   double minDRForUnique, vecSumIMCUTForUnique, minCosPAtomerge, maxPtreltomerge; 
00039 
00040   // a vertex proxy for sorting using stl
00041   struct vertexProxy{
00042     reco::Vertex vert; 
00043     double invm;
00044     size_t ntracks;
00045   };
00046 
00047   // comparison operator for vertexProxy, used in sorting
00048   friend bool operator<(vertexProxy v1, vertexProxy v2){
00049     if(v1.ntracks>2 && v2.ntracks<3) return true;
00050     if(v1.ntracks<3 && v2.ntracks>2) return false;
00051     return (v1.invm>v2.invm);
00052   }
00053 
00054   void resolveBtoDchain(std::vector<vertexProxy>& coll, unsigned int i, unsigned int k);
00055   GlobalVector flightDirection(reco::Vertex &pv, reco::Vertex &sv);
00056 };
00057 
00058 
00059 
00060 BtoCharmDecayVertexMerger::BtoCharmDecayVertexMerger(const edm::ParameterSet &params) :
00061         primaryVertexCollection(params.getParameter<edm::InputTag>("primaryVertices")),
00062         secondaryVertexCollection(params.getParameter<edm::InputTag>("secondaryVertices")),
00063         minDRForUnique(params.getUntrackedParameter<double>("minDRUnique",0.4)),
00064         vecSumIMCUTForUnique(params.getUntrackedParameter<double>("minvecSumIMifsmallDRUnique",5.5)), 
00065         minCosPAtomerge(params.getUntrackedParameter<double>("minCosPAtomerge",0.99)),
00066         maxPtreltomerge(params.getUntrackedParameter<double>("maxPtreltomerge",7777.0))
00067                                                     //  maxFraction(params.getParameter<double>("maxFraction")),
00068                                                     //  minSignificance(params.getParameter<double>("minSignificance"))
00069 {
00070         produces<reco::VertexCollection>();
00071 }
00072 //-----------------------
00073 void BtoCharmDecayVertexMerger::produce(edm::Event &iEvent, const edm::EventSetup &iSetup){
00074 
00075   using namespace reco;
00076   // PV
00077   edm::Handle<reco::VertexCollection> PVcoll;
00078   iEvent.getByLabel(primaryVertexCollection, PVcoll);
00079 
00080   if(PVcoll->size()!=0) {
00081 
00082   const reco::VertexCollection pvc = *( PVcoll.product());
00083   pv = pvc[0];
00084 
00085   // get the IVF collection
00086   edm::Handle<VertexCollection> secondaryVertices;
00087   iEvent.getByLabel(secondaryVertexCollection, secondaryVertices);
00088   
00089  
00090 
00091   //loop over vertices,  fill into own collection for sorting 
00092   std::vector<vertexProxy> vertexProxyColl; 
00093   for(std::vector<reco::Vertex>::const_iterator sv = secondaryVertices->begin();
00094       sv != secondaryVertices->end(); ++sv) {
00095     vertexProxy vtx = {*sv,(*sv).p4().M(),(*sv).tracksSize()};
00096     vertexProxyColl.push_back( vtx );
00097   }
00098   
00099   // sort the vertices by mass and track multiplicity
00100   sort( vertexProxyColl.begin(), vertexProxyColl.end()); 
00101   
00102   
00103   // loop forward over all vertices
00104   for(unsigned int iVtx=0; iVtx < vertexProxyColl.size(); iVtx++){
00105 
00106     // nested loop backwards (in order to start from light masses)
00107     // check all vertices against each other for B->D chain
00108     for(unsigned int kVtx=vertexProxyColl.size()-1; kVtx>iVtx; kVtx--){
00109       // remove D vertices from the collection and add the tracks to the original one
00110       resolveBtoDchain(vertexProxyColl, iVtx, kVtx); 
00111     }
00112   }
00113 
00114   // now create new vertex collection and add to event
00115   VertexCollection *bvertices = new VertexCollection();
00116   for(std::vector<vertexProxy>::iterator it=vertexProxyColl.begin(); it!=vertexProxyColl.end(); it++) bvertices->push_back((*it).vert); 
00117   std::auto_ptr<VertexCollection> bvertColl(bvertices);
00118   iEvent.put(bvertColl);
00119   }
00120   else{
00121     std::auto_ptr<VertexCollection> bvertCollEmpty(new VertexCollection);
00122     iEvent.put(bvertCollEmpty);
00123   }
00124 }
00125 
00126 
00127 //------------------------------------
00128 void BtoCharmDecayVertexMerger::resolveBtoDchain(std::vector<vertexProxy>& coll, unsigned int i, unsigned int k){
00129   using namespace reco;
00130 
00131   GlobalVector momentum1 = GlobalVector(coll[i].vert.p4().X(), coll[i].vert.p4().Y(), coll[i].vert.p4().Z());
00132   GlobalVector momentum2 = GlobalVector(coll[k].vert.p4().X(), coll[k].vert.p4().Y(), coll[k].vert.p4().Z());
00133 
00134   reco::SecondaryVertex sv1(pv, coll[i].vert, momentum1 , true);
00135   reco::SecondaryVertex sv2(pv, coll[k].vert, momentum2 , true);
00136 
00137  
00138   // find out which one is near and far
00139   reco::SecondaryVertex svNear = sv1;
00140   reco::SecondaryVertex svFar  = sv2;
00141   GlobalVector momentumNear  = momentum1;
00142   GlobalVector momentumFar   = momentum2;
00143   
00144   // swap if it is the other way around
00145   if(sv1.dist3d().value() >= sv2.dist3d().value()){
00146     svNear = sv2;
00147     svFar = sv1;
00148     momentumNear = momentum2;
00149     momentumFar = momentum1;
00150   }
00151   GlobalVector nearToFar = flightDirection( svNear, svFar);
00152   GlobalVector pvToNear  = flightDirection( pv, svNear);
00153 
00154   double cosPA =  nearToFar.dot(momentumFar) / momentumFar.mag()/ nearToFar.mag();  
00155   double cosa  =  pvToNear. dot(momentumFar) / pvToNear.mag()   / momentumFar.mag(); 
00156   double ptrel = sqrt(1.0 - cosa*cosa)* momentumFar.mag();  
00157   
00158   double vertexDeltaR = Geom::deltaR(flightDirection(pv, sv1), flightDirection(pv, sv2) );
00159 
00160   // create a set of all tracks from both vertices, avoid double counting by using a std::set<>
00161   std::set<reco::TrackRef> trackrefs;
00162   // first vertex
00163   for(reco::Vertex::trackRef_iterator ti = sv1.tracks_begin(); ti!=sv1.tracks_end(); ti++){
00164     if(sv1.trackWeight(*ti)>0.5){
00165       reco::TrackRef t = ti->castTo<reco::TrackRef>();
00166       trackrefs.insert(t);
00167     }
00168   }
00169   // second vertex
00170   for(reco::Vertex::trackRef_iterator ti = sv2.tracks_begin(); ti!=sv2.tracks_end(); ti++){
00171     if(sv2.trackWeight(*ti)>0.5){
00172       reco::TrackRef t = ti->castTo<reco::TrackRef>();
00173       trackrefs.insert(t);
00174     }
00175   }
00176   
00177   // now calculate one LorentzVector from the track momenta
00178   ROOT::Math::LorentzVector<ROOT::Math::PxPyPzM4D<double> > mother;
00179   for(std::set<reco::TrackRef>::const_iterator it = trackrefs.begin(); it!= trackrefs.end(); it++){
00180     ROOT::Math::LorentzVector<ROOT::Math::PxPyPzM4D<double> >  temp ( (*it)->px(),(*it)->py(),(*it)->pz(), 0.13957 );
00181     mother += temp;
00182   }
00183  
00184 
00185   //   check if criteria for merging are fulfilled
00186   if(vertexDeltaR<minDRForUnique && mother.M()<vecSumIMCUTForUnique && cosPA>minCosPAtomerge && fabs(ptrel)<maxPtreltomerge ) {
00187     
00188     
00189     // add tracks of second vertex which are missing to first vertex
00190     // loop over the second
00191     bool bFoundDuplicate=false;
00192     for(reco::Vertex::trackRef_iterator ti = sv2.tracks_begin(); ti!=sv2.tracks_end(); ti++){
00193       reco::Vertex::trackRef_iterator it = find(sv1.tracks_begin(), sv1.tracks_end(), *ti);
00194       if (it==sv1.tracks_end()) coll[i].vert.add( *ti, sv2.refittedTrack(*ti), sv2.trackWeight(*ti) );     
00195       else bFoundDuplicate=true;
00196     }
00197     // 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.
00198     // the weight must be the larger one, otherwise we may have outliers which are not real outliers
00199     
00200     if(bFoundDuplicate){
00201       // create backup track containers from main vertex
00202       std::vector<TrackBaseRef > tracks_;
00203       std::vector<Track> refittedTracks_;
00204       std::vector<float> weights_;
00205       for(reco::Vertex::trackRef_iterator it = coll[i].vert.tracks_begin(); it!=coll[i].vert.tracks_end(); it++) {
00206         tracks_.push_back( *it);
00207         refittedTracks_.push_back( coll[i].vert.refittedTrack(*it));
00208         weights_.push_back( coll[i].vert.trackWeight(*it) );
00209       }
00210       // delete tracks and add all tracks back, and check in which vertex the weight is larger
00211       coll[i].vert.removeTracks();
00212       std::vector<Track>::iterator it2 = refittedTracks_.begin();
00213       std::vector<float>::iterator it3 = weights_.begin();
00214       for(reco::Vertex::trackRef_iterator it = tracks_.begin(); it!=tracks_.end(); it++, it2++, it3++){
00215         float weight = *it3;
00216         float weight2= sv2.trackWeight(*it);
00217         Track refittedTrackWithLargerWeight = *it2;
00218         if( weight2 >weight) { 
00219           weight = weight2;
00220           refittedTrackWithLargerWeight = sv2.refittedTrack(*it);
00221         }
00222         coll[i].vert.add(*it , refittedTrackWithLargerWeight  , weight);
00223       }
00224     }
00225     
00226     // remove the second vertex from the collection
00227     coll.erase( coll.begin() + k  );
00228   }
00229   
00230 }
00231 //-------------
00232 
00233 
00234 GlobalVector
00235 BtoCharmDecayVertexMerger::flightDirection(reco::Vertex &pv, reco::Vertex &sv){
00236   GlobalVector res(sv.position().X() - pv.position().X(),
00237                     sv.position().Y() - pv.position().Y(),
00238                     sv.position().Z() - pv.position().Z());
00239   return res;
00240 }
00241 
00242 
00243 
00244 DEFINE_FWK_MODULE(BtoCharmDecayVertexMerger);