CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch13/src/PhysicsTools/RecoUtils/src/CandKinematicVertexFitter.cc

Go to the documentation of this file.
00001 #include "PhysicsTools/RecoUtils/interface/CandKinematicVertexFitter.h"
00002 #include "RecoVertex/KinematicFit/interface/KinematicParticleFitter.h"
00003 #include "RecoVertex/KinematicFit/interface/MassKinematicConstraint.h"
00004 #include "DataFormats/GsfTrackReco/interface/GsfTrack.h"
00005 #include "DataFormats/Candidate/interface/VertexCompositeCandidate.h"
00006 #include "FWCore/Utilities/interface/EDMException.h"
00007 #include "SimGeneral/HepPDTRecord/interface/ParticleDataTable.h"
00008 #include "FWCore/Framework/interface/ESHandle.h" 
00009 #include <sstream>
00010 #include <iostream>
00011 using namespace reco;
00012 using namespace std;
00013 
00014 // perform the kinematic fit
00015 bool CandKinematicVertexFitter::fit(const vector<RefCountedKinematicParticle> & particles) const {
00016   try {
00017     tree_ = fitter_.fit(particles);
00018   } catch (std::exception & err) {
00019     std::cerr << ">>> exception thrown by KinematicParticleVertexFitter:\n"
00020               << err.what() << "\n"
00021               << ">>> candidate not fitted to common vertex" << std::endl;
00022     return false;
00023   }
00024   //check tree_ is valid here!
00025   if (tree_->isValid())
00026     return true;
00027   else return false;
00028 }
00029 
00030 // main method called by CandProducer sets the VertexCompositeCandidate
00031 void CandKinematicVertexFitter::set(VertexCompositeCandidate & c) const {
00032   if(bField_ == 0)
00033     throw edm::Exception(edm::errors::InvalidReference)
00034       << "B-Field was not set up CandKinematicVertexFitter.\n"
00035       << "the following method must be called before fitting a candidate:\n"
00036       << " CandKinematicVertexFitter:.set( const MagneticField * )" << endl;
00037   vector<RefCountedKinematicParticle> particles;
00038   vector<Candidate *> daughters;
00039   vector<RecoCandidate::TrackType> trackTypes;
00040   // fill particles with KinematicParticles and daughters with Candidates of the daughters of c
00041   fill(particles, daughters, trackTypes, c);
00042   assert(particles.size() == daughters.size());
00043 
00044   // attempt to fit the KinematicParticles, particles  
00045   if(fit(particles)) {
00046     // after the fit, tree_ contains the KinematicTree from the fit
00047     tree_->movePointerToTheTop();
00048     // set the kinematic properties of the daughters from the fit
00049     RefCountedKinematicVertex vertex = tree_->currentDecayVertex();
00050     if(vertex->vertexIsValid()) {
00051       Candidate::Point vtx(vertex->position());
00052       c.setVertex(vtx);
00053       vector<RefCountedKinematicParticle> treeParticles = tree_->daughterParticles();
00054       vector<RefCountedKinematicParticle>::const_iterator particleIt = treeParticles.begin();
00055       vector<Candidate *>::const_iterator daughterIt = daughters.begin(), daughtersEnd = daughters.end();
00056       vector<RecoCandidate::TrackType>::const_iterator trackTypeIt = trackTypes.begin();
00057       Candidate::LorentzVector mp4(0, 0, 0, 0);
00058       for(; daughterIt != daughtersEnd; ++ particleIt, ++ daughterIt, ++trackTypeIt) {
00059         Candidate & daughter = * * daughterIt;
00060         GlobalVector p3 = (*particleIt)->currentState().globalMomentum();
00061         double px = p3.x(), py = p3.y(), pz = p3.z(), p = p3.mag();
00062         double energy;
00063 
00064         if(!daughter.longLived()) daughter.setVertex(vtx);
00065         double scale;
00066         switch(*trackTypeIt) {
00067         case RecoCandidate::gsfTrackType :
00068           //gsf used for electron tracks
00069           energy = daughter.energy();
00070           scale = energy / p;
00071           px *= scale; py *= scale; pz *= scale; 
00072         default:
00073           double mass = (*particleIt)->currentState().mass();
00074           energy = sqrt(p*p + mass*mass);
00075         };
00076         Candidate::LorentzVector dp4(px, py, pz, energy);
00077         daughter.setP4(dp4);
00078         mp4 += dp4;
00079       }
00080       c.setP4(mp4);
00081       c.setChi2AndNdof(chi2_ = vertex->chiSquared(), ndof_ = vertex->degreesOfFreedom());
00082       GlobalError err = vertex->error();
00083       cov_(0,0) = err.cxx();
00084       cov_(0,1) = err.cyx();
00085       cov_(0,2) = err.czx();
00086       cov_(1,2) = err.czy();
00087       cov_(1,1) = err.cyy();
00088       cov_(2,2) = err.czz();
00089       c.setCovariance(cov_);
00090     }
00091   } else {
00092     c.setChi2AndNdof(chi2_ = -1, ndof_ = 0);
00093     c.setCovariance(cov_ = CovarianceMatrix(ROOT::Math::SMatrixIdentity())); 
00094   }
00095 }
00096 
00097 // methond to fill the properties of a CompositeCandidate's daughters
00098 void CandKinematicVertexFitter::fill(vector<RefCountedKinematicParticle> & particles, 
00099                                      vector<Candidate *> & daughters,
00100                                      vector<RecoCandidate::TrackType> & trackTypes,
00101                                      Candidate & c) const {
00102   size_t nDau = c.numberOfDaughters();
00103   // loop through CompositeCandidate daughters
00104   for(unsigned int j = 0; j < nDau ; ++j) {
00105     Candidate * d = c.daughter(j);
00106     if(d == 0) {
00107       ostringstream message;
00108       message << "Can't access in write mode candidate daughters. "
00109               << "pdgId = " << c.pdgId() << ".\n";
00110       const Candidate * d1 = c.daughter(j);
00111       if(d1 == 0)
00112         message << "Null daughter also found in read-only mode\n";
00113       else
00114         message << "Daughter found in read-only mode with id: " << d1->pdgId() << "\n";
00115       throw edm::Exception(edm::errors::InvalidReference) << message.str();
00116     }
00117     //check for a daughter which itself is a composite
00118     if(d->numberOfDaughters() > 0) {
00119       //try to cast to VertexCompositeCandiate
00120       VertexCompositeCandidate * vtxDau = dynamic_cast<VertexCompositeCandidate*>(d);
00121       if( vtxDau!=0 && vtxDau->vertexChi2()>0 ) {
00122         // if VertexCompositeCandidate refit vtxDau via the set method
00123         (*this).set(*vtxDau);
00124         // if mass constraint is desired, do it here
00125         if ( vtxDau->massConstraint() ) {
00126           KinematicParticleFitter csFitter;
00127           //get particle mass from pdg table via pdgid number
00128           const ParticleData *data = pdt_->particle(vtxDau->pdgId());
00129           ParticleMass mass = data->mass();
00130           float mass_sigma = mass*0.000001; //needs a sigma for the fit
00131           // create a KinematicConstraint and refit the tree with it
00132           //KinematicConstraint * mass_c = new MassKinematicConstraint(mass,mass_sigma);
00133           MassKinematicConstraint mkc(mass,mass_sigma);
00134           KinematicConstraint * mass_c(&mkc);
00135           tree_ = csFitter.fit(mass_c,tree_);
00136           //CHECK THIS! the following works, but might not be safe
00137           //tree_ = csFitter.fit(&(MassKinematicConstraint(mass,mass_sigma)),tree_);
00138         }
00139         // add the kinematic particle from the fit to particles
00140         RefCountedKinematicParticle current = (*this).currentParticle();
00141         particles.push_back(current);
00142         daughters.push_back(d);
00143         trackTypes.push_back(RecoCandidate::noTrackType);       
00144       } else {
00145         fill(particles, daughters, trackTypes, *d);
00146       }
00147     } else {
00148       //get track, make KinematicParticle and add to particles so it can be fit
00149       TrackRef trk = d->get<TrackRef>();
00150       RecoCandidate::TrackType type = d->get<RecoCandidate::TrackType>();
00151       if (!trk.isNull()){
00152         TransientTrack trTrk(trk, bField_);
00153         float chi2 = 0, ndof = 0;
00154         ParticleMass mass = d->mass();
00155         float sigma = mass *1.e-6;
00156         particles.push_back(factory_.particle(trTrk, mass, chi2, ndof, sigma));
00157         daughters.push_back(d);
00158         trackTypes.push_back(type);
00159       } else {
00160         cerr << ">>> warning: candidate of type " << d->pdgId() 
00161              << " has no track reference." << endl;
00162       }
00163     }
00164   }
00165 }