CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_5/src/PhysicsTools/KinFitter/src/TFitParticleEScaledMomDev.cc

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