CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_2_7_hltpatch2/src/RecoEgamma/EgammaPhotonAlgos/src/ConversionVertexFinder.cc

Go to the documentation of this file.
00001 #include <iostream>
00002 #include <memory>
00003 #include "RecoEgamma/EgammaPhotonAlgos/interface/ConversionVertexFinder.h"
00004 // Framework
00005 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00006 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00007 //
00008 //Kinematic constraint vertex fitter
00009 #include <RecoVertex/KinematicFitPrimitives/interface/KinematicParticleFactoryFromTransientTrack.h>
00010 
00011 #include "CommonTools/Statistics/interface/ChiSquaredProbability.h"
00012 
00013 // new templated one
00014 #include "RecoVertex/KinematicFit/interface/KinematicConstrainedVertexFitterT.h"
00015 #include "RecoVertex/KinematicFit/interface/ColinearityKinematicConstraintT.h"
00016 
00017 
00018 
00019 //
00020 
00021 
00022 
00023 
00024 
00025 ConversionVertexFinder::ConversionVertexFinder(const edm::ParameterSet& config ):
00026   conf_(config)
00027 { 
00028   LogDebug("ConversionVertexFinder") << "ConversionVertexFinder CTOR  " <<  "\n";  
00029   maxDelta_ = conf_.getParameter<double>("maxDelta");
00030   maxReducedChiSq_ = conf_.getParameter<double>("maxReducedChiSq");
00031   minChiSqImprovement_  = conf_.getParameter<double>("minChiSqImprovement");
00032   maxNbrOfIterations_ = conf_.getParameter<int>("maxNbrOfIterations");
00033   kcvFitter_ = new KinematicConstrainedVertexFitter();
00034   kcvFitter_->setParameters(conf_);
00035 
00036 }
00037 
00038 ConversionVertexFinder::~ConversionVertexFinder() {
00039 
00040   LogDebug("ConversionVertexFinder") << "ConversionVertexFinder DTOR " <<  "\n";  
00041   delete   kcvFitter_;
00042  
00043 }
00044 
00045 bool  ConversionVertexFinder::run(std::vector<reco::TransientTrack>  pair, reco::Vertex& the_vertex) {
00046   bool found= false;
00047 
00048   if ( pair.size() < 2) return found;
00049   
00050   float sigma = 0.00000001;
00051   float chi = 0.;
00052   float ndf = 0.;
00053   float mass = 0.000511;
00054   
00055 
00056   KinematicParticleFactoryFromTransientTrack pFactory;
00057   
00058   std::vector<RefCountedKinematicParticle> particles;
00059   
00060   particles.push_back(pFactory.particle (pair[0],mass,chi,ndf,sigma,*pair[0].innermostMeasurementState().freeState()));
00061   particles.push_back(pFactory.particle (pair[1],mass,chi,ndf,sigma,*pair[1].innermostMeasurementState().freeState()));
00062   
00063 
00064 #ifdef OldKineFit
00065   ColinearityKinematicConstraint constr(ColinearityKinematicConstraint::PhiTheta);
00066   
00067   RefCountedKinematicTree myTree = kcvFitter_->fit(particles, &constr);
00068 #else
00069 
00070   // bizzare way to the get the field...
00071   const MagneticField* mf = pair[0].field();
00072 
00073   ColinearityKinematicConstraintT<colinearityKinematic::PhiTheta> constr;
00074   KinematicConstrainedVertexFitterT<2,2> kcvFitter(mf);
00075   kcvFitter.setParameters(conf_);
00076   RefCountedKinematicTree myTree =  kcvFitter.fit(particles, &constr);
00077 
00078 #ifdef KineFitDebug
00079 
00080   ColinearityKinematicConstraint oldconstr(ColinearityKinematicConstraint::PhiTheta);
00081   
00082   RefCountedKinematicTree oldTree = kcvFitter_->fit(particles, &oldconstr);
00083 
00084 
00085   if( oldTree->isValid() ) {
00086     std::cout << "old " << kcvFitter_->getNit() << std::endl;
00087     RefCountedKinematicVertex gamma_dec_vertex = oldTree->currentDecayVertex();                                                          
00088     if( gamma_dec_vertex->vertexIsValid())
00089       std::cout << gamma_dec_vertex->chiSquared() <<  " " << gamma_dec_vertex->degreesOfFreedom() << std::endl;
00090     std::cout <<  oldTree->currentParticle()->currentState().globalMomentum() <<
00091       oldTree->currentParticle()->currentState().globalPosition()<< std::endl;
00092     std::vector<RefCountedKinematicParticle> fStates=oldTree->finalStateParticles();
00093     for (unsigned int kk=0; kk<fStates.size(); kk++) {
00094       std::cout <<  fStates[kk]->currentState().globalMomentum() << 
00095         fStates[kk]->currentState().globalPosition() << std::endl;
00096       std::cout <<  fStates[kk]->currentState().kinematicParametersError().matrix()<<std::endl;
00097     }
00098   } else       std::cout << "old invalid " << kcvFitter_->getNit() << std::endl;
00099   
00100   if( myTree->isValid() ) {
00101     std::cout << "new " << kcvFitter.getNit() << std::endl;
00102     RefCountedKinematicVertex gamma_dec_vertex = myTree->currentDecayVertex();                                                          
00103     if( gamma_dec_vertex->vertexIsValid())
00104       std::cout << gamma_dec_vertex->chiSquared() <<  " " << gamma_dec_vertex->degreesOfFreedom() << std::endl;
00105     std::cout <<  myTree->currentParticle()->currentState().globalMomentum() <<
00106       myTree->currentParticle()->currentState().globalPosition()<< std::endl;
00107     std::vector<RefCountedKinematicParticle> fStates=myTree->finalStateParticles();
00108     for (unsigned int kk=0; kk<fStates.size(); kk++) {
00109       std::cout <<  fStates[kk]->currentState().globalMomentum() << 
00110         fStates[kk]->currentState().globalPosition() << std::endl;
00111       std::cout <<  fStates[kk]->currentState().kinematicParametersError().matrix()<<std::endl;
00112     }
00113   } else       std::cout << "new invalid " << kcvFitter.getNit() << std::endl;
00114 
00115 #endif // TemplateKineFitDebug
00116 
00117 #endif
00118 
00119   if( myTree->isValid() ) {
00120     myTree->movePointerToTheTop();                                                                                
00121     RefCountedKinematicParticle the_photon = myTree->currentParticle();                                           
00122     if (the_photon->currentState().isValid()){                                                                    
00123       //const ParticleMass photon_mass = the_photon->currentState().mass();                                       
00124       RefCountedKinematicVertex gamma_dec_vertex;                                                               
00125       gamma_dec_vertex = myTree->currentDecayVertex();                                                          
00126       if( gamma_dec_vertex->vertexIsValid() ){                                                                  
00127         const float chi2Prob = ChiSquaredProbability(gamma_dec_vertex->chiSquared(), gamma_dec_vertex->degreesOfFreedom());
00128         if (chi2Prob>0.){// no longer cut here, only ask positive probability here 
00129           //const math::XYZPoint vtxPos(gamma_dec_vertex->position());                                           
00130           the_vertex = *gamma_dec_vertex;
00131           found = true;
00132         }
00133       }
00134     }
00135   }
00136  
00137   return found;
00138 }
00139 
00140 
00141 TransientVertex  ConversionVertexFinder::run(std::vector<reco::TransientTrack>  pair) {
00142   LogDebug("ConversionVertexFinder") << "ConversionVertexFinder run pair size " << pair.size() <<  "\n";  
00143   
00144   //for ( std::vector<reco::TransientTrack>::const_iterator iTk=pair.begin(); iTk!=pair.end(); ++iTk) {
00145   // LogDebug("ConversionVertexFinder") << "  ConversionVertexFinder  Tracks in the pair  charge " << iTk->charge() << " Num of RecHits " << iTk->recHitsSize() << " inner momentum " << iTk->track().innerMomentum() << "\n";  
00146   //}
00147 
00148 
00149   reco::Vertex theVertex;  
00150   KalmanVertexFitter fitter(true);
00151   TransientVertex transientVtx;
00152 
00153   const std::string metname =  "ConversionVertexFinder| ConversionVertexFinder";
00154   try{
00155 
00156     transientVtx = fitter.vertex(pair); 
00157 
00158   }  catch ( cms::Exception& e ) {
00159 
00160 
00161     edm::LogWarning(metname) << "cms::Exception caught in ConversionVertexFinder::run\n"
00162                              << e.explainSelf();
00163     
00164   }
00165   
00166 
00167   return transientVtx;
00168 
00169     
00170     
00171 }
00172