CMS 3D CMS Logo

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

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