CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC2_patch1/src/PhysicsTools/KinFitter/src/TFitParticleMomDev.cc

Go to the documentation of this file.
00001 // Classname: TFitParticleMomDev
00002 // Author: Jan E. Sundermann, Verena Klose (TU Dresden)      
00003 
00004 
00005 //________________________________________________________________
00006 // 
00007 // TFitParticleMomDev
00008 // --------------------
00009 //
00010 // Particle with special parametrization of the momentum 4vector and
00011 // free mass (4 free parameters). The parametrization is chosen as
00012 // follows:
00013 //
00014 // p = r*|p|*u_r + theta*u_theta + phi*u_phi
00015 // E(fit) =  Sqrt( |p|^2 + d^2*m^2 )
00016 //
00017 // with u_r = p/|p|
00018 //      u_phi = (u_z x u_r)/|u_z x u_r|
00019 //      u_theta = (u_r x u_phi)/|u_r x u_phi|
00020 //
00021 // The initial parameters values are chosen like (r, theta, phi, d) = (1., 0., 0., 1.)
00022 // corresponding to the measured momentum and mass.
00023 //
00024 
00025 #include <iostream>
00026 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00027 #include "PhysicsTools/KinFitter/interface/TFitParticleMomDev.h"
00028 #include "TMath.h"
00029 
00030 
00031 //----------------
00032 // Constructor --
00033 //----------------
00034 TFitParticleMomDev::TFitParticleMomDev()
00035   :TAbsFitParticle()  
00036 {
00037   init(0, 0);
00038 }
00039 
00040 TFitParticleMomDev::TFitParticleMomDev( const TFitParticleMomDev& fitParticle )
00041   :TAbsFitParticle( fitParticle.GetName(), fitParticle.GetTitle() )
00042 {
00043 
00044   _nPar = fitParticle._nPar;
00045   _u1 = fitParticle._u1;
00046   _u2 = fitParticle._u2;
00047   _u3 = fitParticle._u3;
00048   _covMatrix.ResizeTo(  fitParticle._covMatrix );
00049   _covMatrix = fitParticle._covMatrix;
00050   _iniparameters.ResizeTo( fitParticle._iniparameters );
00051   _iniparameters = fitParticle._iniparameters;
00052   _parameters.ResizeTo( fitParticle._parameters );
00053   _parameters = fitParticle._parameters;
00054   _pini = fitParticle._pini;
00055   _pcurr = fitParticle._pcurr;
00056 
00057 }
00058 
00059 TFitParticleMomDev::TFitParticleMomDev(TLorentzVector* pini, const TMatrixD* theCovMatrix)
00060   :TAbsFitParticle()  
00061 {
00062   init(pini, theCovMatrix);
00063 }
00064 
00065 TFitParticleMomDev::TFitParticleMomDev(const TString &name, const TString &title, 
00066                            TLorentzVector* pini, const TMatrixD* theCovMatrix)
00067   :TAbsFitParticle(name, title)  
00068 {
00069   init(pini, theCovMatrix);
00070 }
00071 
00072 TAbsFitParticle* TFitParticleMomDev::clone( TString newname ) const {
00073   // Returns a copy of itself
00074   
00075   TAbsFitParticle* myclone = new TFitParticleMomDev( *this );
00076   if ( newname.Length() > 0 ) myclone->SetName(newname);
00077   return myclone;
00078 
00079 }
00080 
00081 //--------------
00082 // Destructor --
00083 //--------------
00084 TFitParticleMomDev::~TFitParticleMomDev() {
00085 
00086 }
00087 
00088 //--------------
00089 // Operations --
00090 //--------------
00091 void TFitParticleMomDev::init(TLorentzVector* pini, const TMatrixD* theCovMatrix ) {
00092 
00093   _nPar = 4;
00094   _iniparameters.ResizeTo(_nPar,1);
00095   _iniparameters(0,0)=1.;
00096   _iniparameters(1,0)=0.;
00097   _iniparameters(2,0)=0.;
00098   _iniparameters(3,0)=1.;
00099   _parameters.ResizeTo(_nPar,1);
00100   _parameters = _iniparameters;
00101   setIni4Vec(pini);
00102   setCovMatrix(theCovMatrix);
00103 
00104 
00105 }
00106 
00107 void TFitParticleMomDev::setIni4Vec(const TLorentzVector* pini) {
00108   // Set the initial 4vector. Will also set the 
00109   // inital parameter values
00110 
00111   if (pini == 0) {
00112 
00113     _u1.SetXYZ(0., 0., 0.);
00114     _u3.SetXYZ(0., 0., 0.);
00115     _u2.SetXYZ(0., 0., 0.);
00116     _pini.SetXYZT(0., 0., 0., 0.);
00117     _pcurr = _pini;
00118     
00119   } else {
00120     
00121     _pini = (*pini);
00122     _pcurr = _pini;
00123 
00124     _u1 = pini->Vect();
00125     _u1 *= 1./_u1.Mag();
00126     
00127     TVector3 uz(0., 0., 1.);
00128     _u3 = uz.Cross(_u1);
00129     _u3 *= 1./_u3.Mag();
00130     
00131     _u2 = _u3.Cross(_u1);
00132     _u2 *= 1./_u2.Mag();
00133 
00134   }
00135 
00136   // reset parameters
00137   _parameters = _iniparameters;
00138 
00139 }
00140 
00141 TLorentzVector* TFitParticleMomDev::calc4Vec( const TMatrixD* params ) {
00142   // Calculates a 4vector corresponding to the given
00143   // parameter values
00144 
00145   if (params == 0) {
00146     return 0;
00147   }
00148 
00149   if ( params->GetNcols() != 1 || params->GetNrows() !=_nPar ) {
00150     edm::LogError ("WrongMatrixSize")
00151       << GetName() << "::calc4Vec - Parameter matrix has wrong size.";
00152     return 0;
00153   }
00154 
00155   Double_t X = _pini.P() * (*params)(0,0) *_u1.X() +
00156     (*params)(1,0) * _u2.X()+
00157     (*params)(2,0) * _u3.X();
00158   Double_t Y = _pini.P() * (*params)(0,0) *_u1.Y() +
00159     (*params)(1,0) * _u2.Y()+
00160     (*params)(2,0) * _u3.Y();
00161   Double_t Z = _pini.P() * (*params)(0,0) *_u1.Z() +
00162     (*params)(1,0) * _u2.Z()+
00163     (*params)(2,0) * _u3.Z();
00164   Double_t E =  TMath::Sqrt( X*X + Y*Y + Z*Z + (*params)(3,0)*(*params)(3,0)*_pini.M2() );
00165 
00166   TLorentzVector* vec = new TLorentzVector( X, Y, Z, E );
00167   return vec;
00168 
00169 }
00170 
00171 TMatrixD* TFitParticleMomDev::getDerivative() {
00172   // returns derivative dP/dy with P=(p,E) and y=(r, theta, phi, d) 
00173   // the free parameters of the fit. The columns of the matrix contain 
00174   // (dP/dr, dP/dtheta, dP/dphi, dP/dd).
00175 
00176   TMatrixD* DerivativeMatrix = new TMatrixD(4,4);
00177   (*DerivativeMatrix) *= 0.;
00178 
00179   //1st column: dP/dr
00180   (*DerivativeMatrix)(0,0)=_pini.P()*_u1.X();
00181   (*DerivativeMatrix)(1,0)=_pini.P()*_u1.Y();
00182   (*DerivativeMatrix)(2,0)=_pini.P()*_u1.Z();
00183   (*DerivativeMatrix)(3,0)=0.;
00184 
00185   //2nd column: dP/dtheta
00186   (*DerivativeMatrix)(0,1)=_u2.X();
00187   (*DerivativeMatrix)(1,1)=_u2.Y();
00188   (*DerivativeMatrix)(2,1)=_u2.Z();
00189   (*DerivativeMatrix)(3,1)=0.;
00190 
00191    //3rd column: dP/dphi
00192   (*DerivativeMatrix)(0,2)=_u3.X();
00193   (*DerivativeMatrix)(1,2)=_u3.Y();
00194   (*DerivativeMatrix)(2,2)=_u3.Z();
00195   (*DerivativeMatrix)(3,2)=0.;
00196 
00197    //4th column: dP/dm
00198   (*DerivativeMatrix)(0,3)=0.;
00199   (*DerivativeMatrix)(1,3)=0.;
00200   (*DerivativeMatrix)(2,3)=0.;
00201   (*DerivativeMatrix)(3,3)=_pini.M()*_pini.M()*_parameters(3,0)/_pcurr.E();
00202 
00203   return DerivativeMatrix;  
00204 }
00205 
00206 TMatrixD* TFitParticleMomDev::transform(const TLorentzVector& vec) {
00207   // Returns the parameters corresponding to the given 
00208   // 4vector wrt. to the current base vectors u_r, u_theta, and u_phi
00209 
00210   // construct rotation matrix
00211   TRotation rot;
00212   rot.RotateAxes( _u1, _u2, _u3 );
00213   rot.Invert();
00214 
00215   // rotate vector
00216   TVector3 vec3( vec.Vect() );
00217   vec3.Transform( rot );
00218 
00219   // retrieve parameters
00220   TMatrixD* tparams = new TMatrixD( _nPar, 1 );
00221   (*tparams)(0,0) = vec3(0)/_pini.P();
00222   (*tparams)(1,0) = vec3(1);
00223   (*tparams)(2,0) = vec3(2);
00224   (*tparams)(3,0) = vec.M()/_pini.M();
00225 
00226   return tparams;
00227 
00228 }