CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_5/src/RecoVertex/KinematicFit/src/FinalTreeBuilder.cc

Go to the documentation of this file.
00001 #include "RecoVertex/KinematicFit/interface/FinalTreeBuilder.h"
00002 #include "RecoVertex/KinematicFitPrimitives/interface/KinematicRefittedTrackState.h"
00003 #include "RecoVertex/KinematicFitPrimitives/interface/KinematicStatePropagator.h"
00004 #include "RecoVertex/KinematicFitPrimitives/interface/KinematicState.h"
00005 //#include "Vertex/KinematicFitPrimitives/interface/KinematicLinearizedTrackState.h"
00006 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00007 
00008 FinalTreeBuilder::FinalTreeBuilder()
00009 {
00010  kvFactory = new KinematicVertexFactory();
00011  KinematicStatePropagator * ksp = 0;
00012  pFactory = new VirtualKinematicParticleFactory(ksp);
00013 }
00014 
00015 FinalTreeBuilder::~FinalTreeBuilder()
00016 {
00017  delete kvFactory;
00018  delete pFactory;
00019 }
00020 
00021 RefCountedKinematicTree FinalTreeBuilder::buildTree(const CachingVertex<6>& vtx,
00022                              const std::vector<RefCountedKinematicParticle> &input) const
00023 {
00024 //creating resulting kinematic particle
00025  AlgebraicVector7 par;
00026  AlgebraicMatrix cov(7,7,0);
00027  par(0) = vtx.position().x();
00028  par(1) = vtx.position().y();
00029  par(2) = vtx.position().z();
00030  double en = 0.;
00031  int ch = 0;
00032 
00033 //new particle momentum calculation and refitted kinematic states
00034  std::vector<KinematicRefittedTrackState *> rStates;
00035  std::vector<RefCountedVertexTrack> refTracks = vtx.tracks();
00036  for(std::vector<RefCountedVertexTrack>::const_iterator i = refTracks.begin();i !=refTracks.end();++i)
00037  {
00038   KinematicRefittedTrackState * rs = dynamic_cast<KinematicRefittedTrackState *>(&(*((*i)->refittedState())));
00039   AlgebraicVector4 f_mom = rs->kinematicMomentumVector();
00040   par(3) += f_mom(0);
00041   par(4) += f_mom(1);
00042   par(5) += f_mom(2);
00043   en += sqrt(f_mom(1)*f_mom(1)+f_mom(2)*f_mom(2)+f_mom(3)*f_mom(3) + f_mom(0)*f_mom(0));
00044   ch += (*i)->linearizedTrack()->charge();
00045   rStates.push_back(rs);
00046  }
00047 
00048 //math precision check (numerical stability)
00049  double differ = en*en - (par(3)*par(3)+par(4)*par(4)+par(5)*par(5));
00050  if(differ>0.)
00051  {
00052   par(6) = sqrt(differ);
00053  }else{
00054    LogDebug("FinalTreeBuilder")
00055         << "Fit failed: Current precision does not allow to calculate the mass\n";
00056    return ReferenceCountingPointer<KinematicTree>(new KinematicTree());
00057  }
00058 
00059 // covariance matrix calculation: momentum-momentum components part (4x4)
00060 // and position-momentum components part:
00061  AlgebraicMatrix m_all = momentumPart(vtx,par);
00062 
00063 //position-position components part (3x3)
00064  // AlgebraicMatrix x_X = vtx.error().matrix();
00065 
00066 //making new matrix itself
00067  cov.sub(1,1,m_all);
00068 
00069 //covariance sym matrix
00070  AlgebraicSymMatrix77 sCov;
00071  for(int i = 1; i<8; i++)
00072  {
00073   for(int j = 1; j<8; j++)
00074   {
00075    if(i<=j) sCov(i-1,j-1) = cov(i,j);
00076   }
00077  }
00078 
00079 //valid decay vertex for our resulting particle
00080  RefCountedKinematicVertex dVrt = kvFactory->vertex(vtx);
00081 
00082 //invalid production vertex for our resulting particle
00083  RefCountedKinematicVertex pVrt = kvFactory->vertex();
00084 
00085 //new born kinematic particle
00086  KinematicParameters kPar(par);
00087  KinematicParametersError kEr(sCov);
00088  const MagneticField* field=input.front()->magneticField();
00089  KinematicState nState(kPar, kEr, ch, field);
00090 
00091 //invalid previous particle and empty constraint:
00092  KinematicParticle * zp = 0;
00093  RefCountedKinematicParticle pPart = ReferenceCountingPointer<KinematicParticle>(zp);
00094 
00095  float vChi = vtx.totalChiSquared();
00096  float vNdf = vtx.degreesOfFreedom();
00097  RefCountedKinematicParticle nPart = pFactory->particle(nState, vChi, vNdf, pPart);
00098 
00099 //adding top particle to the tree
00100  RefCountedKinematicTree resTree = ReferenceCountingPointer<KinematicTree>(new KinematicTree());
00101  resTree->addParticle(pVrt,dVrt,nPart);
00102 
00103 //making refitted kinematic particles and putting them to the tree
00104  std::vector<RefCountedKinematicParticle> rrP;
00105 
00106  std::vector<RefCountedKinematicParticle>::const_iterator j;
00107  std::vector<RefCountedVertexTrack>::const_iterator i;
00108  for(j=input.begin(), i=refTracks.begin(); j !=input.end(), i !=refTracks.end();++j, ++i)
00109  {
00110   RefCountedLinearizedTrackState lT = (*i)->linearizedTrack();
00111   KinematicRefittedTrackState * rS= dynamic_cast<KinematicRefittedTrackState *>(&(*((*i)->refittedState())));
00112 
00113 //   RefCountedRefittedTrackState rS = (*i)->refittedState();
00114   AlgebraicVector7 lPar = rS->kinematicParameters();
00115   KinematicParameters lkPar(lPar);
00116   AlgebraicSymMatrix77 lCov = rS->kinematicParametersCovariance();
00117   KinematicParametersError lkCov(lCov);
00118   TrackCharge lch = lT->charge();
00119   KinematicState nState(lkPar,lkCov,lch, field);
00120   RefCountedKinematicParticle nPart = (*j)->refittedParticle(nState,vChi,vNdf);
00121   rrP.push_back(nPart);
00122   if((*j)->correspondingTree() != 0)
00123   {
00124 
00125 //here are the particles having trees after them
00126    KinematicTree * tree = (*j)->correspondingTree();
00127    tree->movePointerToTheTop();
00128    tree->replaceCurrentParticle(nPart);
00129    RefCountedKinematicVertex cdVertex = resTree->currentDecayVertex();
00130    resTree->addTree(cdVertex, tree);
00131   }else{
00132 
00133 //here are just particles fitted to this tree
00134    RefCountedKinematicVertex nV = kvFactory->vertex();
00135    resTree->addParticle(dVrt,nV,nPart);
00136   }
00137  }
00138  return resTree;
00139 }
00140 
00141 //method returning the full covariance matrix
00142 //of new born kinematic particle
00143 AlgebraicMatrix FinalTreeBuilder::momentumPart(const CachingVertex<6>& vtx,
00144         const AlgebraicVector7& par) const
00145 {
00146  std::vector<RefCountedVertexTrack> refTracks  = vtx.tracks();
00147  int size = refTracks.size();
00148  AlgebraicMatrix cov(7+4*(size-1),7+4*(size-1));
00149  AlgebraicMatrix jac(7,7+4*(size-1));
00150  jac(1,1) = 1.;
00151  jac(2,2) = 1.;
00152  jac(3,3) = 1.;
00153 
00154  double energy_total = sqrt(par(3)*par(3)+par(6)*par(6) + par(5)*par(5)+par(4)*par(4));
00155 
00156  std::vector<KinematicRefittedTrackState *>::const_iterator rs;
00157  std::vector<RefCountedVertexTrack>::const_iterator rt_i;
00158  int i_int = 0;
00159  for(rt_i = refTracks.begin(); rt_i != refTracks.end(); rt_i++)
00160  {
00161   double a;
00162   AlgebraicVector6 param = (**rt_i).refittedState()->parameters(); // rho, theta, phi,tr_im, z_im, mass
00163   double rho = param[0];
00164   double theta = param[1];
00165   double phi = param[2];
00166   double mass = param[5];
00167 
00168   if ((**rt_i).linearizedTrack()->charge()!=0) {
00169       a = -(**rt_i).refittedState()->freeTrajectoryState().parameters().magneticFieldInInverseGeV(vtx.position()).z()
00170                 * (**rt_i).refittedState()->freeTrajectoryState().parameters ().charge();
00171     if (a==0.) throw cms::Exception("FinalTreeBuilder", "Field is 0");
00172   } else {
00173     a = 1;
00174   }
00175 
00176   AlgebraicMatrix jc_el(4,4,0);
00177   jc_el(1,1) = -a*cos(phi)/(rho*rho); //dpx/d rho
00178   jc_el(2,1) = -a*sin(phi)/(rho*rho); //dpy/d rho
00179   jc_el(3,1) = -a/(rho*rho*tan(theta)); //dpz/d rho
00180 
00181   jc_el(3,2) = -a/(rho*sin(theta)*sin(theta)); //dpz/d theta
00182 
00183   jc_el(1,3) = -a*sin(phi)/rho; //dpx/d phi
00184   jc_el(2,3) = a*cos(phi)/rho; //dpy/d
00185 
00186 //non-trival elements: mass correlations:
00187   double energy_local  = sqrt(a*a/(rho*rho)*(1+1/(tan(theta)*tan(theta)))  + mass*mass);
00188 
00189   jc_el(4,4) = energy_total*mass/(par(6)*energy_local); // dm/dm
00190 
00191   jc_el(4,1) = (-(energy_total/energy_local*a*a/(rho*rho*rho*sin(theta)*sin(theta)) )
00192                 + par(3)*a/(rho*rho)*cos(phi) + par(4)*a/(rho*rho)*sin(phi)
00193                 + par(5)*a/(rho*rho*tan(theta)) )/par(6);       //dm / drho
00194 
00195   jc_el(4,2) = (-(energy_total/energy_local*a*a/(rho*rho*sin(theta)*sin(theta)*tan(theta)) )
00196                 + par(5)*a/(rho*sin(theta)*sin(theta)) )/par(6);//dm d theta
00197 
00198   jc_el(4,3) = ( par(3)*sin(phi) - par(4)*cos(phi) )*a/(rho*par(6));    //dm/dphi
00199 
00200   jac.sub(4,i_int*4+4,jc_el);
00201 
00202 //top left corner elements
00203   if(i_int == 0) {
00204    cov.sub(1,1,asHepMatrix<7>((**rt_i).fullCovariance()));
00205   } else {
00206 //4-momentum corellatons: diagonal elements of the matrix
00207    AlgebraicMatrix fullCovMatrix(asHepMatrix<7>((**rt_i).fullCovariance()));
00208    AlgebraicMatrix m_m_cov = fullCovMatrix.sub(4,7,4,7);
00209    AlgebraicMatrix x_p_cov = fullCovMatrix.sub(1,3,4,7);
00210    AlgebraicMatrix p_x_cov = fullCovMatrix.sub(4,7,1,3);
00211 
00212 // cout << "Full covariance: \n"<< (**rt_i).fullCovariance()<<endl;
00213 // cout << "Full m_m_cov: "<< m_m_cov<<endl;
00214 //  cout<<"p_x_cov"<< p_x_cov<<endl;
00215 //  cout<<"x_p_cov"<< x_p_cov<<endl;
00216 
00217 //putting everything to the joint covariance matrix:
00218 //diagonal momentum-momentum elements:
00219    cov.sub(i_int*4 + 4, i_int*4 + 4,m_m_cov);
00220 
00221 //position momentum elements:
00222    cov.sub(1,i_int*4 + 4,x_p_cov);
00223    cov.sub(i_int*4 + 4,1,p_x_cov);
00224 
00225 //off diagonal elements: corellations
00226 // track momentum - track momentum
00227   }
00228    int j_int = 0;
00229    for(std::vector<RefCountedVertexTrack>::const_iterator rt_j = refTracks.begin(); rt_j != refTracks.end(); rt_j++)
00230    {
00231     if(i_int < j_int)
00232     {
00233      AlgebraicMatrix i_k_cov_m = asHepMatrix<4,4>(vtx.tkToTkCovariance((*rt_i),(*rt_j)));
00234 //     cout<<"i_k_cov_m"<<i_k_cov_m <<endl;
00235      cov.sub(i_int*4 + 4, j_int*4 + 4,i_k_cov_m);
00236      cov.sub(j_int*4 + 4, i_int*4 + 4,i_k_cov_m.T());
00237     }
00238     j_int++;
00239    }
00240   i_int++;
00241  }
00242 // cout<<"jac"<<jac<<endl;
00243 // cout<<"cov"<<cov<<endl;
00244 //  cout << "final result new"<<jac*cov*jac.T()<<endl;
00245 
00246  return jac*cov*jac.T();
00247 }
00248