CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_6/src/PhysicsTools/KinFitter/src/TFitParticleMCMomDev.cc

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