CMS 3D CMS Logo

PFTauDiscriminants::DiscriminantBase< T > Class Template Reference

#include <RecoTauTag/TauTagTools/interface/PFTauDiscriminantBase.h>

Inheritance diagram for PFTauDiscriminants::DiscriminantBase< T >:

PFTauDiscriminants::Discriminant

List of all members.

Public Types

typedef std::vector< T >
::const_iterator 
myVectorIterator

Public Member Functions

void branchTree (TTree *theTree)
 add a branch to a ttree corresponding to this variable
void compute (PFTauDiscriminantManager *input)
 computes the associated quanity for the tau object that is loaded in the PFTauDiscriminantManager implemented in derived implementation class
 DiscriminantBase (string name, string rootTypeName, bool branchAsSimpleDataType, bool isMultiple, T defaultValue)
void fillMVA (std::vector< PhysicsTools::Variable::Value > &mvaHolder) const
virtual void setNullResult (PFTauDiscriminantManager *input)
virtual ~DiscriminantBase ()

Protected Member Functions

virtual void doComputation (PFTauDiscriminantManager *input, vector< T > &result)=0

Private Attributes

defaultValue_
bool isMultiple_
std::vector< T > result_
std::vector< T > * resultPtr_
singleResult_


Detailed Description

template<class T>
class PFTauDiscriminants::DiscriminantBase< T >

Definition at line 70 of file PFTauDiscriminantBase.h.


Member Typedef Documentation

template<class T>
typedef std::vector<T>::const_iterator PFTauDiscriminants::DiscriminantBase< T >::myVectorIterator

Definition at line 77 of file PFTauDiscriminantBase.h.


Constructor & Destructor Documentation

template<class T>
PFTauDiscriminants::DiscriminantBase< T >::DiscriminantBase ( string  name,
string  rootTypeName,
bool  branchAsSimpleDataType,
bool  isMultiple,
defaultValue 
) [inline, explicit]

Definition at line 72 of file PFTauDiscriminantBase.h.

00073                                                                                      :Discriminant(name, rootTypeName, branchAsSimpleDataType),isMultiple_(isMultiple),defaultValue_(defaultValue){
00074          resultPtr_ = &result_;
00075       };

template<class T>
virtual PFTauDiscriminants::DiscriminantBase< T >::~DiscriminantBase (  )  [inline, virtual]

Definition at line 77 of file PFTauDiscriminantBase.h.

00077 {};


Member Function Documentation

template<class T>
void PFTauDiscriminants::DiscriminantBase< T >::branchTree ( TTree *  theTree  )  [inline, virtual]

add a branch to a ttree corresponding to this variable

Implements PFTauDiscriminants::Discriminant.

Definition at line 112 of file PFTauDiscriminantBase.h.

00112                                       {
00113          if (!this->branchSimply())
00114          {
00115             edm::LogInfo("PFTauDiscriminantBase") << "Branching TTree: " << theTree->GetName() << " with full class name (bronch)";
00116             theTree->Branch(name().c_str(), rootTypeName().c_str(), &resultPtr_); 
00117          }
00118          else
00119          {
00120             edm::LogInfo("PFTauDiscriminantBase") << "Branching TTree: " << theTree->GetName() << " with struct style branch (leaflist)";
00121             stringstream branchType;
00122             branchType << name() << "/" << rootTypeName(); //eg D, F, I, etc
00123             theTree->Branch(this->name().c_str(), &singleResult_, branchType.str().c_str());
00124          }
00125       }

template<class T>
void PFTauDiscriminants::DiscriminantBase< T >::compute ( PFTauDiscriminantManager input  )  [inline, virtual]

computes the associated quanity for the tau object that is loaded in the PFTauDiscriminantManager implemented in derived implementation class

Implements PFTauDiscriminants::Discriminant.

Definition at line 88 of file PFTauDiscriminantBase.h.

00089       {
00090          result_.clear();
00091 
00092          if (input)
00093             doComputation(input, result_); 
00094          else
00095             edm::LogError("DiscriminantBase") << "Error in DiscriminantBase - trying to compute discriminants on null PFTauDecayMode pointer!";
00096 
00097          size_t numberOfResultsReturned = result_.size();
00098          if(!numberOfResultsReturned) //if there are no results, ROOT branches of simple variables must be filled w/ the default value
00099          {
00100             singleResult_ = defaultValue_;
00101          } else
00102          {
00103             if(!isMultiple_ && numberOfResultsReturned > 1)
00104             {
00105                edm::LogWarning("PFTauDiscriminants::DiscriminantBase") << "Warning, multiple discriminant values recieved for a non-multiple branch, taking only the first!"; 
00106             }
00107             singleResult_ = result_[0];
00108          }
00109       }

template<class T>
virtual void PFTauDiscriminants::DiscriminantBase< T >::doComputation ( PFTauDiscriminantManager input,
vector< T > &  result 
) [protected, pure virtual]

Implemented in PFTauDiscriminants::DecayMode, PFTauDiscriminants::OutlierNCharged, PFTauDiscriminants::Pt, PFTauDiscriminants::Eta, PFTauDiscriminants::MainTrackPt, PFTauDiscriminants::MainTrackAngle, PFTauDiscriminants::TrackPt, PFTauDiscriminants::PiZeroPt, PFTauDiscriminants::TrackAngle, PFTauDiscriminants::PiZeroAngle, PFTauDiscriminants::Dalitz, PFTauDiscriminants::InvariantMassOfSignal, PFTauDiscriminants::InvariantMass, PFTauDiscriminants::OutlierPt, PFTauDiscriminants::OutlierAngle, PFTauDiscriminants::ChargedOutlierPt, PFTauDiscriminants::ChargedOutlierAngle, PFTauDiscriminants::NeutralOutlierPt, and PFTauDiscriminants::NeutralOutlierAngle.

Referenced by PFTauDiscriminants::DiscriminantBase< int >::compute().

template<class T>
void PFTauDiscriminants::DiscriminantBase< T >::fillMVA ( std::vector< PhysicsTools::Variable::Value > &  mvaHolder  )  const [inline, virtual]

Implements PFTauDiscriminants::Discriminant.

Definition at line 127 of file PFTauDiscriminantBase.h.

00127                                                                           {
00128          if (isMultiple_)
00129          {
00130             for(myVectorIterator aResult = result_.begin(); aResult != result_.end(); ++aResult)
00131             {
00132                mvaHolder.push_back(PhysicsTools::Variable::Value(theAtomicId(), static_cast<double>(*aResult)));
00133             }
00134          }
00135          else
00136          {
00137             mvaHolder.push_back(PhysicsTools::Variable::Value(theAtomicId(), static_cast<double>(singleResult_)));
00138          }
00139       }

template<class T>
virtual void PFTauDiscriminants::DiscriminantBase< T >::setNullResult ( PFTauDiscriminantManager input  )  [inline, virtual]

Implements PFTauDiscriminants::Discriminant.

Definition at line 80 of file PFTauDiscriminantBase.h.

00081       {
00082          result_.clear();
00083          singleResult_ = defaultValue_;
00084       }


Member Data Documentation

template<class T>
T PFTauDiscriminants::DiscriminantBase< T >::defaultValue_ [private]

Definition at line 145 of file PFTauDiscriminantBase.h.

Referenced by PFTauDiscriminants::DiscriminantBase< int >::compute(), and PFTauDiscriminants::DiscriminantBase< int >::setNullResult().

template<class T>
bool PFTauDiscriminants::DiscriminantBase< T >::isMultiple_ [private]

Definition at line 144 of file PFTauDiscriminantBase.h.

Referenced by PFTauDiscriminants::DiscriminantBase< int >::compute(), and PFTauDiscriminants::DiscriminantBase< int >::fillMVA().

template<class T>
std::vector<T> PFTauDiscriminants::DiscriminantBase< T >::result_ [private]

Definition at line 147 of file PFTauDiscriminantBase.h.

Referenced by PFTauDiscriminants::DiscriminantBase< int >::compute(), PFTauDiscriminants::DiscriminantBase< int >::DiscriminantBase(), PFTauDiscriminants::DiscriminantBase< int >::fillMVA(), and PFTauDiscriminants::DiscriminantBase< int >::setNullResult().

template<class T>
std::vector<T>* PFTauDiscriminants::DiscriminantBase< T >::resultPtr_ [private]

Definition at line 148 of file PFTauDiscriminantBase.h.

Referenced by PFTauDiscriminants::DiscriminantBase< int >::branchTree(), and PFTauDiscriminants::DiscriminantBase< int >::DiscriminantBase().

template<class T>
T PFTauDiscriminants::DiscriminantBase< T >::singleResult_ [private]

Definition at line 146 of file PFTauDiscriminantBase.h.

Referenced by PFTauDiscriminants::DiscriminantBase< int >::branchTree(), PFTauDiscriminants::DiscriminantBase< int >::compute(), PFTauDiscriminants::DiscriminantBase< int >::fillMVA(), and PFTauDiscriminants::DiscriminantBase< int >::setNullResult().


The documentation for this class was generated from the following file:
Generated on Tue Jun 9 18:50:00 2009 for CMSSW by  doxygen 1.5.4