CMS 3D CMS Logo

TFitConstraintM.cc

Go to the documentation of this file.
00001 // Classname: TFitConstraintM
00002 // Author: Jan E. Sundermann, Verena Klose (TU Dresden)      
00003 
00004 
00005 //________________________________________________________________
00006 // 
00007 // TFitConstraintM::
00008 // --------------------
00009 //
00010 // Fit constraint: mass conservation ( m_i - m_j - MassConstraint == 0 )
00011 //
00012 
00013 using namespace std;
00014 
00015 #include <iostream>
00016 #include <iomanip>
00017 #include "PhysicsTools/KinFitter/interface/TFitConstraintM.h"
00018 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00019 #include "TLorentzVector.h"
00020 #include "TClass.h"
00021 
00022 ClassImp(TFitConstraintM)
00023 
00024 //----------------
00025 // Constructor --
00026 //----------------
00027 TFitConstraintM::TFitConstraintM()
00028   : TAbsFitConstraint() 
00029   ,_ParList1(0)
00030   ,_ParList2(0)
00031   ,_TheMassConstraint(0)
00032 {
00033 
00034 }
00035 
00036 TFitConstraintM::TFitConstraintM(vector<TAbsFitParticle*>* ParList1,
00037                                  vector<TAbsFitParticle*>* ParList2, Double_t Mass)
00038   : TAbsFitConstraint() 
00039   ,_ParList1(0)
00040   ,_ParList2(0)
00041 {
00042   // ParList1: Vector containing first list of constrained particles 
00043   //           ( sum[ m_i ] - sum[ m_j ] - MassConstraint == 0 )
00044   // ParList2: Vector containing second list of constrained particles 
00045   //           ( sum[ m_i ] - sum[ m_j ]  - MassConstraint == 0 )
00046 
00047   if (ParList1) {
00048     _ParList1 = (*ParList1);
00049   }
00050   if (ParList2) {
00051     _ParList2 = (*ParList2);
00052   }
00053   if (Mass >= 0) {
00054     _TheMassConstraint = Mass;
00055   }
00056   else if(Mass < 0) {
00057     edm::LogWarning ("NegativeMassConstr")
00058       << "Mass constraint in TFitConstraintM cannot be set to a negative value, will be set to 0.";
00059     _TheMassConstraint = 0.;
00060   }
00061 }
00062 
00063 TFitConstraintM::TFitConstraintM(const TString &name, const TString &title,
00064                                  vector<TAbsFitParticle*>* ParList1,
00065                                  vector<TAbsFitParticle*>* ParList2, Double_t Mass)
00066   : TAbsFitConstraint(name, title) 
00067   ,_ParList1(0)
00068   ,_ParList2(0) 
00069 {
00070   // ParList1: Vector containing first list of constrained particles 
00071   //           ( sum[ m_i ] - sum[ m_j ] - MassConstraint == 0 )
00072   // ParList2: Vector containing second list of constrained particles 
00073   //           ( sum[ m_i ] - sum[ m_j ] - MassConstraint == 0 )
00074 
00075   if (ParList1) {
00076     _ParList1 = (*ParList1);
00077   }
00078   if (ParList2) {
00079     _ParList2 = (*ParList2);
00080   }  
00081   if (Mass >= 0) {
00082     _TheMassConstraint = Mass;
00083   }
00084   else if(Mass < 0) {
00085     edm::LogWarning ("NegativeMassConstr")
00086       << "Mass constraint in TFitConstraintM cannot be set to a negative value, will be set to 0.";
00087     _TheMassConstraint = 0.;
00088   }
00089 }
00090 void TFitConstraintM::addParticle1( TAbsFitParticle* particle ) {
00091   // Add one constrained particle to first list of particles
00092   // ( sum[ m_i ] - sum[ m_j ] - MassConstraint == 0 )
00093 
00094   _ParList1.push_back( particle );
00095 
00096 }
00097 
00098 void TFitConstraintM::addParticle2( TAbsFitParticle* particle ) {
00099   // Add one constrained particle to second list of particles
00100   // ( sum[ m_i ] - sum[ m_j ] - MassConstraint == 0 )
00101 
00102   _ParList2.push_back( particle );
00103 
00104 }
00105 
00106 void TFitConstraintM::addParticles1( TAbsFitParticle* p1, TAbsFitParticle* p2, TAbsFitParticle* p3, TAbsFitParticle* p4,
00107                                      TAbsFitParticle* p5, TAbsFitParticle* p6, TAbsFitParticle* p7, TAbsFitParticle* p8,
00108                                      TAbsFitParticle* p9, TAbsFitParticle* p10) {
00109   // Add many constrained particle to first list of particles
00110   // ( sum[ m_i ] - sum[ m_j ] - MassConstraint == 0 )
00111 
00112   if (p1) addParticle1( p1 );
00113   if (p2) addParticle1( p2 );
00114   if (p3) addParticle1( p3 );
00115   if (p4) addParticle1( p4 );
00116   if (p5) addParticle1( p5 );
00117   if (p6) addParticle1( p6 );
00118   if (p7) addParticle1( p7 );
00119   if (p8) addParticle1( p8 );
00120   if (p9) addParticle1( p9 );
00121   if (p10) addParticle1( p10 );
00122 
00123 }
00124 
00125 void TFitConstraintM::addParticles2( TAbsFitParticle* p1, TAbsFitParticle* p2, TAbsFitParticle* p3, TAbsFitParticle* p4,
00126                                      TAbsFitParticle* p5, TAbsFitParticle* p6, TAbsFitParticle* p7, TAbsFitParticle* p8,
00127                                      TAbsFitParticle* p9, TAbsFitParticle* p10) {
00128   // Add many constrained particle to second list of particles
00129   // ( sum[ m_i ] - sum[ m_j ] - MassConstraint == 0 )
00130 
00131   if (p1) addParticle2( p1 );
00132   if (p2) addParticle2( p2 );
00133   if (p3) addParticle2( p3 );
00134   if (p4) addParticle2( p4 );
00135   if (p5) addParticle2( p5 );
00136   if (p6) addParticle2( p6 );
00137   if (p7) addParticle2( p7 );
00138   if (p8) addParticle2( p8 );
00139   if (p9) addParticle2( p9 );
00140   if (p10) addParticle2( p10 );
00141 
00142 }
00143 
00144 
00145 //--------------
00146 // Destructor --
00147 //--------------
00148 TFitConstraintM::~TFitConstraintM() {
00149 
00150 }
00151 
00152 //--------------
00153 // Operations --
00154 //--------------
00155 TMatrixD* TFitConstraintM::getDerivative( TAbsFitParticle* particle ) {
00156   // returns derivative df/dP with P=(p,E) and f the constraint (f=0).
00157   // The matrix contains one row (df/dp, df/dE).
00158 
00159   TMatrixD* DerivativeMatrix = new TMatrixD(1,4);
00160   (*DerivativeMatrix) *= 0.;
00161 
00162   // Pf[4] is the 4-Mom (p,E) of the sum of particles on 
00163   // the list particle is part of 
00164   
00165   Double_t Factor = 0.;
00166   TLorentzVector Pf(0., 0., 0., 0.);
00167 
00168   if( OnList( &_ParList1, particle) ) {
00169     UInt_t Npart = _ParList1.size();
00170     for (unsigned int i=0; i<Npart; i++) {
00171       const TLorentzVector* FourVec = (_ParList1[i])->getCurr4Vec();
00172       Pf += (*FourVec);
00173     }
00174     if( Pf.M() == 0. ) {
00175       edm::LogInfo ("KinFitter")
00176         << "Division by zero in "
00177         << IsA()->GetName() << " (named " << GetName() << ", titled " << GetTitle()
00178         << ") will lead to Inf in derivative matrix for particle "
00179         << particle->GetName() << ".";
00180     }
00181     Factor = 1./ Pf.M();
00182   } else if (OnList( &_ParList2, particle) ) {
00183     UInt_t Npart = _ParList2.size();
00184     for (unsigned int i=0; i<Npart; i++) {
00185       const TLorentzVector* FourVec = (_ParList2[i])->getCurr4Vec();
00186       Pf += (*FourVec);
00187     }
00188     if( Pf.M() == 0. ) {
00189       edm::LogInfo ("KinFitter")
00190         << "Division by zero in "
00191         << IsA()->GetName() << " (named " << GetName() << ", titled " << GetTitle()
00192         << ") will lead to Inf in derivative matrix for particle "
00193         << particle->GetName() << ".";
00194     }
00195     Factor = -1./Pf.M();
00196   } else {
00197     Factor = 0.; 
00198   }
00199   
00200   (*DerivativeMatrix)(0,0) = -Pf[0] ;
00201   (*DerivativeMatrix)(0,1) = -Pf[1];
00202   (*DerivativeMatrix)(0,2) = -Pf[2];
00203   (*DerivativeMatrix)(0,3) = +Pf[3];
00204   (*DerivativeMatrix) *= Factor;
00205 
00206   return DerivativeMatrix;
00207 
00208 }
00209 
00210 
00211 Double_t TFitConstraintM::getInitValue() {
00212   // Get initial value of constraint (before the fit)
00213 
00214   Double_t InitValue(0) ;   
00215   InitValue = CalcMass(&_ParList1,true) - CalcMass(&_ParList2,true) - _TheMassConstraint ;
00216   return InitValue;
00217 }
00218 
00219 Double_t TFitConstraintM::getCurrentValue() {
00220   // Get value of constraint after the fit
00221 
00222   Double_t CurrentValue(0);
00223   CurrentValue= CalcMass(&_ParList1,false) - CalcMass(&_ParList2,false) - _TheMassConstraint;
00224   return CurrentValue;
00225 }
00226  
00227 
00228 Bool_t TFitConstraintM::OnList(vector<TAbsFitParticle*>* List,
00229                                TAbsFitParticle* particle) {
00230   // Checks whether list contains given particle
00231 
00232   Bool_t ok(false);
00233   UInt_t Npart = List->size();
00234   for (unsigned int i=0;i<Npart;i++) {
00235     ok = (particle == (*List)[i]);
00236     if (ok) break;
00237   }
00238   return ok;
00239 }
00240 
00241 Double_t TFitConstraintM::CalcMass(vector<TAbsFitParticle*>* List, Bool_t IniVal) {
00242   // Calculates initial/current invariant mass of provided list of particles
00243 
00244   TLorentzVector P(0., 0., 0., 0.);
00245   UInt_t Npart = List->size();
00246   for (unsigned int i=0;i<Npart;i++) {
00247     const TLorentzVector* FourVec = 0;
00248     if (IniVal)
00249       FourVec = ((*List)[i])->getIni4Vec();
00250     else      
00251       FourVec = ((*List)[i])->getCurr4Vec();
00252     P += (*FourVec);
00253   }
00254   return P.M();
00255 }
00256 
00257 TString TFitConstraintM::getInfoString() {
00258   // Collect information to be used for printout
00259 
00260   stringstream info;
00261   info << scientific << setprecision(6);
00262 
00263   info << "__________________________" << endl
00264        << endl;
00265   info <<"OBJ: " << IsA()->GetName() << "\t" << GetName() << "\t" << GetTitle() << endl;
00266 
00267   info << "initial value: " << getInitValue() << endl;
00268   info << "current value: " << getCurrentValue() << endl;
00269   info << "mass: " << _TheMassConstraint << endl;
00270 
00271   return info.str();
00272 
00273 }
00274 
00275 void TFitConstraintM::print() {
00276   // Print constraint contents
00277 
00278   edm::LogVerbatim("KinFitter") << this->getInfoString();
00279 
00280 }

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