CMS 3D CMS Logo

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