CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_0/src/PhysicsTools/KinFitter/src/TFitParticleEMomDev.cc

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