CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/DataFormats/BTauReco/interface/TrackProbabilityTagInfo.h

Go to the documentation of this file.
00001 #ifndef BTauReco_BJetTagTrackProbability_h
00002 #define BTauReco_BJetTagTrackProbability_h
00003 
00004 #include "DataFormats/BTauReco/interface/RefMacros.h"
00005 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00006 #include "DataFormats/BTauReco/interface/JTATagInfo.h"
00007 #include "DataFormats/JetReco/interface/JetTracksAssociation.h"
00008 
00009 namespace reco {
00010  
00011 class TrackProbabilityTagInfo : public JTATagInfo
00012  {
00013   public:
00014 
00015  TrackProbabilityTagInfo(
00016    std::vector<double> probability2d,
00017    std::vector<double> probability3d,
00018    std::vector<int> trackOrder2d,
00019    std::vector<int> trackOrder3d,const JetTracksAssociationRef & jtaRef) : JTATagInfo(jtaRef),
00020      m_probability2d(probability2d),
00021      m_probability3d(probability3d),
00022      m_trackOrder2d(trackOrder2d),
00023      m_trackOrder3d(trackOrder3d)     {}
00024 
00025   TrackProbabilityTagInfo() {}
00026   
00027   virtual ~TrackProbabilityTagInfo() {}
00028 
00029 int factorial(int n) const
00030 {
00031  if(n<2) return 1;
00032  else return n*factorial(n-1);
00033 }
00034   
00035   virtual float probability(size_t n,int ip) const 
00036    {
00037     if(ip == 0)
00038     {
00039      if(n <m_probability3d.size())
00040       return m_probability3d[n];  
00041     }
00042     else
00043     {
00044      if(n <m_probability2d.size())
00045       return m_probability2d[n];  
00046     }
00047     return -10.; 
00048    }
00049 
00050   virtual float jetProbability(int ip, float minTrackProb) const
00051 {
00052 
00053  const std::vector<double> * vp;
00054    if(ip==0) vp= &m_probability3d;
00055    else vp= &m_probability2d;
00056    const std::vector<double> & v =*vp;
00057 
00058    int ngoodtracks=v.size();
00059    double SumJet=0.;
00060 
00061   for(std::vector<double>::const_iterator q = v.begin(); q != v.end(); q++){
00062     SumJet+=(*q>minTrackProb)?log(*q):log(minTrackProb);
00063   }
00064 
00065   double ProbJet;
00066   double Loginvlog=0;
00067 
00068   if(SumJet<0.){
00069     if(ngoodtracks>=2){
00070       Loginvlog=log(-SumJet);
00071     }
00072     double Prob=1.;
00073     for(int l=1; l!=ngoodtracks; l++){
00074 
00075       Prob+=exp(l*Loginvlog-log(1.*factorial(l)));
00076     }
00077     double LogProb=log(Prob);
00078     ProbJet=
00079       std::min(exp(std::max(LogProb+SumJet,-30.)),1.);
00080   }else{
00081     ProbJet=1.;
00082   }
00083   if(ProbJet>1)
00084    std::cout << "ProbJet too high: "  << ProbJet << std::endl;
00085 
00086   //double LogProbJet=-log(ProbJet);
00087   //return 1.-ProbJet;
00088   return -log10(ProbJet)/4.;
00089 }
00090 
00091    
00100   virtual float discriminator(int ipType, float minProb) const { return jetProbability(ipType,minProb); }
00101 
00102   virtual int selectedTracks(int ipType) const
00103   {
00104    if(ipType == 0) return m_probability3d.size();
00105    else return m_probability2d.size();
00106   }
00107      virtual int trackIndex(size_t n,int ip) const
00108    {
00109     if(ip == 0)
00110     {
00111      if(n <m_probability3d.size())
00112       return m_trackOrder3d[n];
00113     }
00114     else
00115     {
00116      if(n <m_probability2d.size())
00117       return m_trackOrder2d[n];
00118     }
00119     return 0;
00120    }
00121  
00122   virtual const Track & track(size_t n,int ipType) const
00123   {
00124     return *tracks()[trackIndex(n,ipType)];
00125   }
00126  
00127   virtual TrackProbabilityTagInfo* clone() const { return new TrackProbabilityTagInfo( * this ); }
00128   
00129   private:
00130    std::vector<double> m_probability2d;     //
00131    std::vector<double> m_probability3d;     // create a smarter container instead of 
00132    std::vector<int> m_trackOrder2d;         // this pair of vectors. 
00133    std::vector<int> m_trackOrder3d;         //
00134  };
00135 
00136 //typedef edm::ExtCollection< TrackProbabilityTagInfo,JetTagCollection> TrackProbabilityExtCollection;
00137 //typedef edm::OneToOneAssociation<JetTagCollection, TrackProbabilityTagInfo> TrackProbabilityExtCollection;
00138  
00139 DECLARE_EDM_REFS( TrackProbabilityTagInfo )
00140 
00141 }
00142 #endif