CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_5_3_13_patch3/src/RecoVertex/KinematicFit/src/LagrangeParentParticleFitter.cc

Go to the documentation of this file.
00001 #include "RecoVertex/KinematicFit/interface/LagrangeParentParticleFitter.h"
00002 #include "RecoVertex/KinematicFitPrimitives/interface/KinematicVertexFactory.h"
00003 #include "RecoVertex/KinematicFit/interface/InputSort.h"
00004 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00005 
00006 LagrangeParentParticleFitter::LagrangeParentParticleFitter()
00007 {defaultParameters();}
00008 
00009 /*
00010 RefCountedKinematicTree LagrangeParentParticleFitter::fit(RefCountedKinematicTree tree, KinematicConstraint * cs) const
00011 {
00012 //fitting top paticle only
00013  tree->movePointerToTheTop();
00014  RefCountedKinematicParticle particle = tree->currentParticle(); 
00015  RefCountedKinematicVertex inVertex = tree->currentDecayVertex();
00016 
00017 //parameters and covariance matrix of the particle to fit 
00018  AlgebraicVector par = particle->currentState().kinematicParameters().vector();
00019  AlgebraicSymMatrix cov = particle->currentState().kinematicParametersError().matrix();
00020 
00021 //chi2 and ndf
00022  float chi = 0.;
00023  float ndf = 0.;
00024  
00025 //constraint
00026 //here it is defined at 7-point
00027 //taken just at current state parameters
00028  AlgebraicVector vl;
00029  AlgebraicMatrix dr;
00030  AlgebraicVector dev;
00031  
00032 //initial expansion point 
00033  AlgebraicVector exPoint = particle->currentState().kinematicParameters().vector();
00034  
00035 //refitted parameters and covariance matrix:
00036 //simple and symmetric 
00037  AlgebraicVector refPar;
00038  AlgebraicSymMatrix refCovS; 
00039  
00040  int nstep = 0;
00041  double df = 0.;
00042  do{ 
00043    chi = particle->chiSquared();
00044    ndf = particle->degreesOfFreedom();
00045    df =  0.;
00046    
00047 //derivative and value at expansion point  
00048    vl = cs->value(exPoint).first;
00049    dr = cs->derivative(exPoint).first;
00050    dev = cs->deviations();  
00051    
00052 //residual between expansion and current parameters
00053 //0 at the first step   
00054    AlgebraicVector delta_alpha = par - exPoint;  
00055   
00056 //parameters needed for refit
00057 // v_d = (D * V_alpha & * D.T)^(-1)   
00058    AlgebraicMatrix drt = dr.T();
00059    AlgebraicMatrix v_d = dr * cov * drt;
00060    int ifail = 0;
00061    v_d.invert(ifail);
00062    if(ifail != 0) throw VertexException("ParentParticleFitter::error inverting covariance matrix");
00063   
00064 //lagrangian multipliers  
00065 //lambda = V_d * (D * delta_alpha + d)
00066    AlgebraicVector lambda = v_d * (dr*delta_alpha + vl);
00067   
00068 //refitted parameters
00069    refPar = par - (cov * drt * lambda);  
00070   
00071 //refitted covariance: simple and SymMatrix  
00072    refCovS = cov;
00073    AlgebraicMatrix sPart = drt * v_d * dr;
00074    AlgebraicMatrix covF = cov * sPart * cov; 
00075    AlgebraicSymMatrix sCovF(7,0);
00076   
00077    for(int i = 1; i<8; i++)
00078    {
00079      for(int j = 1; j<8; j++)
00080      {if(i<=j) sCovF(i,j) = covF(i,j);}
00081    }
00082    refCovS -= sCovF; 
00083   
00084    for(int i = 1; i < 8; i++)
00085    {refCovS(i,i) += dev(i);}
00086   
00087 //chiSquared  
00088    chi +=  (lambda.T() * (dr*delta_alpha + vl))(1);
00089    ndf +=  cs->numberOfEquations();
00090    
00091 //new expansionPoint
00092    exPoint = refPar;
00093    AlgebraicVector vlp = cs->value(exPoint).first;
00094    for(int i = 1; i< vl.num_row();i++)
00095    {df += std::abs(vlp(i));}
00096    nstep++; 
00097   }while(df>theMaxDiff && nstep<theMaxStep);
00098   
00100 //creating an output KinematicParticle 
00101 //and putting it on its place in the tree
00102   KinematicParameters param(refPar);
00103   KinematicParametersError er(refCovS);
00104   KinematicState kState(param,er,particle->initialState().particleCharge());
00105   RefCountedKinematicParticle refParticle  = particle->refittedParticle(kState,chi,ndf,cs->clone());                                                     
00106   tree->replaceCurrentParticle(refParticle);
00107   
00108 //replacing the  vertex with its refitted version
00109   GlobalPoint nvPos(param.vector()(1), param.vector()(2), param.vector()(3)); 
00110   AlgebraicSymMatrix nvMatrix = er.matrix().sub(1,3);
00111   GlobalError nvError(asSMatrix<3>(nvMatrix));
00112   VertexState vState(nvPos, nvError, 1.0);
00113   KinematicVertexFactory vFactory;
00114   RefCountedKinematicVertex nVertex = vFactory.vertex(vState, inVertex,chi,ndf);
00115   tree->replaceCurrentVertex(nVertex);
00116    
00117   return tree;
00118 }
00119 */
00120 
00121 std::vector<RefCountedKinematicTree>  LagrangeParentParticleFitter::fit(std::vector<RefCountedKinematicTree> trees, 
00122                                                                            KinematicConstraint * cs)const                                         
00123 {
00124  
00125  InputSort iSort;
00126  std::vector<RefCountedKinematicParticle> prt = iSort.sort(trees);
00127  int nStates = prt.size(); 
00128 
00129 //making full initial parameters and covariance
00130  AlgebraicVector part(7*nStates,0);
00131  AlgebraicSymMatrix cov(7*nStates,0);
00132  
00133  AlgebraicVector chi_in(nStates,0);
00134  AlgebraicVector ndf_in(nStates,0);
00135  int l_c=0;
00136  for(std::vector<RefCountedKinematicParticle>::const_iterator i = prt.begin(); i != prt.end(); i++)
00137  {
00138   AlgebraicVector7 lp = (*i)->currentState().kinematicParameters().vector();
00139   for(int j = 1; j != 8; j++){part(7*l_c + j) = lp(j-1);}
00140   AlgebraicSymMatrix lc= asHepMatrix<7>((*i)->currentState().kinematicParametersError().matrix());
00141   cov.sub(7*l_c+1,lc);
00142   chi_in(l_c+1) = (*i)->chiSquared();
00143   ndf_in(l_c+1) = (*i)->degreesOfFreedom();
00144   l_c++;
00145  }
00146 //refitted parameters and covariance matrix:
00147 //simple and symmetric 
00148  AlgebraicVector refPar;
00149  AlgebraicSymMatrix refCovS; 
00150  
00151 //constraint values, derivatives and deviations:  
00152  AlgebraicVector vl;
00153  AlgebraicMatrix dr;
00154  AlgebraicVector dev;
00155  int nstep = 0;
00156  double df = 0.; 
00157  AlgebraicVector exPoint = part;
00158  
00159 //this piece of code should later be refactored:
00160 // The algorithm is the same as above, but the
00161 // number of refitted particles is >1. Smart way of
00162 // refactoring should be chosen for it.
00163  AlgebraicVector chi;
00164  AlgebraicVector ndf;
00165 // cout<<"Starting the main loop"<<endl;
00166  do{
00167   df = 0.;
00168   chi = chi_in;
00169   ndf = ndf_in;
00170 //  cout<<"Iterational linearization point: "<<exPoint<<endl;
00171   vl = cs->value(exPoint).first;
00172   dr = cs->derivative(exPoint).first;
00173   dev = cs->deviations(nStates);
00174 
00175 //  cout<<"The value : "<<vl<<endl;
00176 //  cout<<"The derivative: "<<dr<<endl;
00177 //  cout<<"deviations: "<<dev<<endl;
00178 //  cout<<"covariance "<<cov<<endl;
00179 
00180 //residual between expansion and current parameters
00181 //0 at the first step   
00182   AlgebraicVector delta_alpha = part - exPoint;  
00183 
00184 //parameters needed for refit
00185 // v_d = (D * V_alpha & * D.T)^(-1)   
00186    AlgebraicMatrix drt = dr.T();
00187    AlgebraicMatrix v_d = dr * cov * drt;
00188    int ifail = 0;
00189    v_d.invert(ifail);
00190    if(ifail != 0) {
00191      LogDebug("KinematicConstrainedVertexFitter")
00192         << "Fit failed: unable to invert covariance matrix\n";
00193      return std::vector<RefCountedKinematicTree>();
00194    }
00195   
00196 //lagrangian multipliers  
00197 //lambda = V_d * (D * delta_alpha + d)
00198    AlgebraicVector lambda = v_d * (dr*delta_alpha + vl);
00199   
00200 //refitted parameters
00201    refPar = part - (cov * drt * lambda);  
00202   
00203 //refitted covariance: simple and SymMatrix  
00204    refCovS = cov;
00205    AlgebraicMatrix sPart = drt * v_d * dr;
00206    AlgebraicMatrix covF = cov * sPart * cov; 
00207 
00208 //total covariance symmatrix    
00209   AlgebraicSymMatrix sCovF(nStates*7,0);
00210   for(int i = 1; i< nStates*7 +1; ++i)
00211   {
00212    for(int j = 1; j< nStates*7 +1; j++)
00213    {if(i<=j) sCovF(i,j) = covF(i,j);}
00214   }
00215   
00216   refCovS -= sCovF; 
00217   
00218 //  cout<<"Fitter: refitted covariance "<<refCovS<<endl;
00219   for(int i = 1; i < nStates+1; i++)
00220   {for(int j = 1; j<8; j++){refCovS((i-1)+j,(i-1)+j) += dev(j);}} 
00221   
00222 //chiSquared  
00223   for(int k =1; k<nStates+1; k++)
00224   {
00225    chi(k) +=  (lambda.T() * (dr*delta_alpha + vl))(1);
00226    ndf(k) +=  cs->numberOfEquations();
00227   }  
00228 //new expansionPoint
00229   exPoint = refPar;
00230   AlgebraicVector vlp = cs->value(exPoint).first;
00231   for(int i = 1; i< vl.num_row();i++)
00232   {df += std::abs(vlp(i));}
00233   nstep++; 
00234  }while(df>theMaxDiff && nstep<theMaxStep);
00235 //here math and iterative part is finished, starting an output production
00236 //creating an output KinematicParticle and putting it on its place in the tree
00237  
00238 //vector of refitted particles and trees 
00239  std::vector<RefCountedKinematicParticle> refPart;
00240  std::vector<RefCountedKinematicTree> refTrees = trees;
00241  
00242  int j=1;
00243  std::vector<RefCountedKinematicTree>::const_iterator tr = refTrees.begin();
00244  for(std::vector<RefCountedKinematicParticle>::const_iterator i = prt.begin(); i!= prt.end(); i++)
00245  {
00246   AlgebraicVector7 lRefPar;
00247   for(int k = 1; k<8 ; k++)
00248   {lRefPar(k-1) = refPar((j-1)*7+k);}
00249   AlgebraicSymMatrix77 lRefCovS = asSMatrix<7>(refCovS.sub((j-1)*7 +1,(j-1)*7+7));
00250   
00251 //new refitted parameters and covariance  
00252   KinematicParameters param(lRefPar);
00253   KinematicParametersError er(lRefCovS); 
00254   KinematicState kState(param,er,(*i)->initialState().particleCharge(), (**i).magneticField());
00255   RefCountedKinematicParticle refParticle  = (*i)->refittedParticle(kState,chi(j),ndf(j),cs->clone());
00256   
00257 //replacing particle itself  
00258   (*tr)->findParticle(*i);
00259   RefCountedKinematicVertex inVertex =  (*tr)->currentDecayVertex();
00260   (*tr)->replaceCurrentParticle(refParticle);
00261   
00262 //replacing the  vertex with its refitted version
00263   GlobalPoint nvPos(param.position()); 
00264   AlgebraicSymMatrix nvMatrix = asHepMatrix<7>(er.matrix()).sub(1,3);
00265   GlobalError nvError(asSMatrix<3>(nvMatrix));
00266   VertexState vState(nvPos, nvError, 1.0);
00267   KinematicVertexFactory vFactory;
00268   RefCountedKinematicVertex nVertex = vFactory.vertex(vState,inVertex,chi(j),ndf(j));
00269   (*tr)->replaceCurrentVertex(nVertex);  
00270   tr++;
00271   j++;
00272  }
00273  return refTrees; 
00274 }     
00275 
00276 void LagrangeParentParticleFitter::setParameters(const edm::ParameterSet& pSet)
00277 {
00278   theMaxDiff = pSet.getParameter<double>("maxDistance");
00279   theMaxStep = pSet.getParameter<int>("maxNbrOfIterations");;
00280 }
00281 
00282 void LagrangeParentParticleFitter::defaultParameters()
00283 {
00284   theMaxDiff = 0.001;
00285   theMaxStep = 100;
00286 }