CMS 3D CMS Logo

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

Generated on Tue Jun 9 17:41:16 2009 for CMSSW by  doxygen 1.5.4