CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC4_patch1/src/RecoVertex/KinematicFit/plugins/KineExample.cc

Go to the documentation of this file.
00001 #include "RecoVertex/KinematicFit/plugins/KineExample.h"
00002 
00003 #include "DataFormats/VertexReco/interface/Vertex.h"
00004 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00005 #include "DataFormats/TrackReco/interface/Track.h"
00006 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00007 #include "DataFormats/Common/interface/Handle.h"
00008 #include "FWCore/Framework/interface/MakerMacros.h"
00009 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00010 #include "FWCore/Framework/interface/ESHandle.h"
00011 
00012 #include "TrackingTools/TransientTrack/interface/TransientTrackBuilder.h"
00013 #include "TrackingTools/Records/interface/TransientTrackRecord.h"
00014 #include "TrackingTools/TransientTrack/interface/TransientTrack.h"
00015 #include "RecoVertex/VertexPrimitives/interface/TransientVertex.h"
00016 #include "RecoVertex/KalmanVertexFit/interface/KalmanVertexFitter.h"
00017 #include "SimTracker/Records/interface/TrackAssociatorRecord.h"
00018 
00019 #include "RecoVertex/KinematicFitPrimitives/interface/ParticleMass.h"
00020 #include "RecoVertex/KinematicFitPrimitives/interface/MultiTrackKinematicConstraint.h"
00021 #include <RecoVertex/KinematicFitPrimitives/interface/KinematicParticleFactoryFromTransientTrack.h>
00022 // #include "RecoVertex/KinematicFitPrimitives/interface/"
00023 #include "RecoVertex/KinematicFit/interface/KinematicConstrainedVertexFitter.h"
00024 #include "RecoVertex/KinematicFit/interface/TwoTrackMassKinematicConstraint.h"
00025 #include "RecoVertex/KinematicFit/interface/KinematicParticleVertexFitter.h"
00026 #include "RecoVertex/KinematicFit/interface/KinematicParticleFitter.h"
00027 #include "RecoVertex/KinematicFit/interface/MassKinematicConstraint.h"
00028 
00029 #include <iostream>
00030 
00031 using namespace reco;
00032 using namespace edm;
00033 using namespace std;
00034 
00035 KineExample::KineExample(const edm::ParameterSet& iConfig)
00036   : theConfig(iConfig)
00037 {
00038   trackLabel_ = iConfig.getParameter<std::string>("TrackLabel");
00039   outputFile_ = iConfig.getUntrackedParameter<std::string>("outputFile");
00040   kvfPSet = iConfig.getParameter<edm::ParameterSet>("KVFParameters");
00041 //   rootFile_ = TFile::Open(outputFile_.c_str(),"RECREATE");
00042   edm::LogInfo("RecoVertex/KineExample")
00043     << "Initializing KVF TEST analyser  - Output file: " << outputFile_ <<"\n";
00044 }
00045 
00046 
00047 KineExample::~KineExample() {
00048 //   delete rootFile_;
00049 }
00050 
00051 void KineExample::beginRun(Run const& run, EventSetup const& setup){
00052   edm::ESHandle<TrackAssociatorBase> theAssociatorForParamAtPca;
00053   setup.get<TrackAssociatorRecord>().get("TrackAssociatorByChi2",theAssociatorForParamAtPca);
00054   associatorForParamAtPca = (TrackAssociatorByChi2 *) theAssociatorForParamAtPca.product();
00055 
00056 //   tree = new SimpleVertexTree("VertexFitter", associatorForParamAtPca);
00057 }
00058 
00059 
00060 void KineExample::endJob() {
00061 //   delete tree;
00062 }
00063 
00064 //
00065 // member functions
00066 //
00067 
00068 void
00069 KineExample::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup)
00070 {
00071   try {
00072 
00073   cout << "Reconstructing event number: " << iEvent.id() << "\n";
00074 
00075   // get RECO tracks from the event
00076   // `tks` can be used as a ptr to a reco::TrackCollection
00077   edm::Handle<reco::TrackCollection> tks;
00078   iEvent.getByLabel(trackLabel_, tks);
00079   if (!tks.isValid()) {
00080     cout
00081       << "Couln't find track collection: " << iEvent.id()
00082       << "\n";
00083   } else {
00084 
00085     edm::LogInfo("RecoVertex/KineExample")
00086       << "Found: " << (*tks).size() << " reconstructed tracks" << "\n";
00087     cout << "got " << (*tks).size() << " tracks " << endl;
00088 
00089     // Transform Track to TransientTrack
00090     //get the builder:
00091     edm::ESHandle<TransientTrackBuilder> theB;
00092     iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder",theB);
00093     //do the conversion:
00094     vector<TransientTrack> t_tks = (*theB).build(tks);
00095 
00096      cout  << "Found: " << t_tks.size() << " reconstructed tracks" << "\n";
00097 
00098     // Do a KindFit, if >= 4 tracks.
00099     if (t_tks.size() > 3) {
00100 
00101       // For a first test, suppose that the first four tracks are the 2 muons,
00102       // then the 2 kaons. Since this will not be true, the result of the fit
00103       // will not be meaningfull, but at least you will get the idea of how to
00104       // do such a fit.
00105 
00106       //First, to get started, a simple vertex fit:
00107 
00108       vector<TransientTrack> ttv;
00109       ttv.push_back(t_tks[0]); ttv.push_back(t_tks[1]); ttv.push_back(t_tks[2]);ttv.push_back(t_tks[3]);
00110       KalmanVertexFitter kvf(false);
00111       TransientVertex tv = kvf.vertex(ttv);
00112       if (!tv.isValid()) cout << "KVF failed\n";
00113       else std::cout << "KVF fit Position: " << Vertex::Point(tv.position()) << "\n";
00114 
00115 
00116       TransientTrack ttMuPlus = t_tks[0];
00117       TransientTrack ttMuMinus = t_tks[1];
00118       TransientTrack ttKPlus = t_tks[2];
00119       TransientTrack ttKMinus = t_tks[3];
00120 
00121       //the final state muons and kaons from the Bs->J/PsiPhi->mumuKK decay
00122       //Creating a KinematicParticleFactory
00123       KinematicParticleFactoryFromTransientTrack pFactory;
00124 
00125       //The mass of a muon and the insignificant mass sigma to avoid singularities in the covariance matrix.
00126       ParticleMass muon_mass = 0.1056583;
00127       ParticleMass kaon_mass = 0.493677;
00128       float muon_sigma = 0.0000001;
00129       float kaon_sigma = 0.000016;
00130 
00131       //initial chi2 and ndf before kinematic fits. The chi2 of the reconstruction is not considered
00132       float chi = 0.;
00133       float ndf = 0.;
00134 
00135       //making particles
00136       vector<RefCountedKinematicParticle> muonParticles;
00137       vector<RefCountedKinematicParticle> phiParticles;
00138       vector<RefCountedKinematicParticle> allParticles;
00139       muonParticles.push_back(pFactory.particle (ttMuPlus,muon_mass,chi,ndf,muon_sigma));
00140       muonParticles.push_back(pFactory.particle (ttMuMinus,muon_mass,chi,ndf,muon_sigma));
00141       allParticles.push_back(pFactory.particle (ttMuPlus,muon_mass,chi,ndf,muon_sigma));
00142       allParticles.push_back(pFactory.particle (ttMuMinus,muon_mass,chi,ndf,muon_sigma));
00143 
00144       phiParticles.push_back(pFactory.particle (ttKPlus,kaon_mass,chi,ndf,kaon_sigma));
00145       phiParticles.push_back(pFactory.particle (ttKMinus,kaon_mass,chi,ndf,kaon_sigma));
00146       allParticles.push_back(pFactory.particle (ttKPlus,kaon_mass,chi,ndf,kaon_sigma));
00147       allParticles.push_back(pFactory.particle (ttKMinus,kaon_mass,chi,ndf,kaon_sigma));
00148 
00149       /* Example of a simple vertex fit, without other constraints
00150        * The reconstructed decay tree is a result of the kinematic fit
00151        * The KinematicParticleVertexFitter fits the final state particles to their vertex and
00152        * reconstructs the decayed state
00153        */
00154       KinematicParticleVertexFitter fitter;
00155       cout <<"Simple vertex fit with KinematicParticleVertexFitter:\n";
00156       RefCountedKinematicTree vertexFitTree = fitter.fit(allParticles);
00157 
00158       printout(vertexFitTree);
00159 
00161 
00162         //creating the constraint for the J/Psi mass
00163         ParticleMass jpsi = 3.09687;
00164 
00165         //creating the two track mass constraint
00166         MultiTrackKinematicConstraint *  j_psi_c = new  TwoTrackMassKinematicConstraint(jpsi);
00167 
00168         //creating the fitter
00169         KinematicConstrainedVertexFitter kcvFitter;
00170 
00171         //obtaining the resulting tree
00172         RefCountedKinematicTree myTree = kcvFitter.fit(allParticles, j_psi_c);
00173 
00174         cout << "\nGlobal fit done:\n";
00175         printout(myTree);
00176 
00177         //creating the vertex fitter
00178         KinematicParticleVertexFitter kpvFitter;
00179 
00180         //reconstructing a J/Psi decay
00181         RefCountedKinematicTree jpTree = kpvFitter.fit(muonParticles);
00182 
00183         //creating the particle fitter
00184         KinematicParticleFitter csFitter;
00185 
00186         // creating the constraint
00187         float jp_m_sigma = 0.00004;
00188         KinematicConstraint * jpsi_c2 = new MassKinematicConstraint(jpsi,jp_m_sigma);
00189 
00190         //the constrained fit:
00191         jpTree = csFitter.fit(jpsi_c2,jpTree);
00192 
00193         //getting the J/Psi KinematicParticle and putting it together with the kaons.
00194         //The J/Psi KinematicParticle has a pointer to the tree it belongs to
00195         jpTree->movePointerToTheTop();
00196         RefCountedKinematicParticle jpsi_part = jpTree->currentParticle();
00197         phiParticles.push_back(jpsi_part);
00198 
00199         //making a vertex fit and thus reconstructing the Bs parameters
00200         // the resulting tree includes all the final state tracks, the J/Psi meson,
00201         // its decay vertex, the Bs meson and its decay vertex.
00202         RefCountedKinematicTree bsTree = kpvFitter.fit(phiParticles);
00203         cout << "Sequential fit done:\n";
00204         printout(bsTree);
00205 
00206 
00207 
00208 //       // For the analysis: compare to your SimVertex
00209 //       TrackingVertex sv = getSimVertex(iEvent);
00210 //   edm::Handle<TrackingParticleCollection>  TPCollectionH ;
00211 //   iEvent.getByLabel("trackingtruth","TrackTruth",TPCollectionH);
00212 //   const TrackingParticleCollection tPC = *(TPCollectionH.product());
00213 //       reco::RecoToSimCollection recSimColl=associatorForParamAtPca->associateRecoToSim(tks,
00214 //                                                                            TPCollectionH,
00215 //                                                                            &iEvent);
00216 //
00217 //       tree->fill(tv, &sv, &recSimColl);
00218 //     }
00219 
00220     }
00221   }
00222 
00223   }
00224   catch (std::exception & err) {
00225     cout  << "Exception during event number: " << iEvent.id()
00226       << "\n" << err.what() << "\n";
00227   }
00228 
00229 }
00230 
00231 void KineExample::printout(const RefCountedKinematicVertex& myVertex) const
00232 {
00233   if (myVertex->vertexIsValid()) {
00234     cout << "Decay vertex: " << myVertex->position() <<myVertex->chiSquared()<< " "<<myVertex->degreesOfFreedom()<<endl;
00235   } else cout << "Decay vertex Not valid\n";
00236 }
00237 
00238 void KineExample::printout(const RefCountedKinematicParticle& myParticle) const
00239 {
00240   cout << "Particle: \n";
00241 //accessing the reconstructed Bs meson parameters:
00242 //SK: uncomment if needed  AlgebraicVector7 bs_par = myParticle->currentState().kinematicParameters().vector();
00243 
00244 //and their joint covariance matrix:
00245 //SK:uncomment if needed  AlgebraicSymMatrix77 bs_er = myParticle->currentState().kinematicParametersError().matrix();
00246   cout << "Momentum at vertex: " << myParticle->currentState().globalMomentum ()<<endl;
00247   cout << "Parameters at vertex: " << myParticle->currentState().kinematicParameters().vector()<<endl;
00248 }
00249 
00250 void KineExample::printout(const RefCountedKinematicTree& myTree) const
00251 {
00252   if (!myTree->isValid()) {
00253     cout <<"Tree is invalid. Fit failed.\n";
00254     return;
00255   }
00256 
00257 //accessing the tree components, move pointer to top
00258   myTree->movePointerToTheTop();
00259 
00260 //We are now at the top of the decay tree getting the B_s reconstructed KinematicPartlcle
00261   RefCountedKinematicParticle b_s = myTree->currentParticle();
00262   printout(b_s);
00263 
00264 // The B_s decay vertex
00265   RefCountedKinematicVertex b_dec_vertex = myTree->currentDecayVertex();
00266   printout(b_dec_vertex);
00267 
00268   // Get all the children of Bs:
00269   //In this way, the pointer is not moved
00270   vector< RefCountedKinematicParticle > bs_children = myTree->finalStateParticles();
00271 
00272   for (unsigned int i=0;i< bs_children.size();++i) {
00273     printout(bs_children[i]);
00274   }
00275 
00276 //Now navigating down the tree , pointer is moved:
00277   bool child = myTree->movePointerToTheFirstChild();
00278 
00279   if(child) while (myTree->movePointerToTheNextChild()) {
00280     RefCountedKinematicParticle aChild = myTree->currentParticle();
00281     printout(aChild);
00282   }
00283 }
00284 
00285 //Returns the first vertex in the list.
00286 
00287 TrackingVertex KineExample::getSimVertex(const edm::Event& iEvent) const
00288 {
00289    // get the simulated vertices
00290   edm::Handle<TrackingVertexCollection>  TVCollectionH ;
00291   iEvent.getByLabel("trackingtruth","VertexTruth",TVCollectionH);
00292   const TrackingVertexCollection tPC = *(TVCollectionH.product());
00293 
00294 //    Handle<edm::SimVertexContainer> simVtcs;
00295 //    iEvent.getByLabel("g4SimHits", simVtcs);
00296 //    std::cout << "SimVertex " << simVtcs->size() << std::endl;
00297 //    for(edm::SimVertexContainer::const_iterator v=simVtcs->begin();
00298 //        v!=simVtcs->end(); ++v){
00299 //      std::cout << "simvtx "
00300 //             << v->position().x() << " "
00301 //             << v->position().y() << " "
00302 //             << v->position().z() << " "
00303 //             << v->parentIndex() << " "
00304 //             << v->noParent() << " "
00305 //               << std::endl;
00306 //    }
00307    return *(tPC.begin());
00308 }
00309 DEFINE_FWK_MODULE(KineExample);