CMS 3D CMS Logo

Classes | Public Member Functions | Static Public Attributes | Private Attributes | Friends

PhysicsTools::TreeReader Class Reference

#include <TreeReader.h>

List of all members.

Classes

struct  Bool
class  Value

Public Member Functions

void addBranch (const std::string &expression, AtomicId name=AtomicId(), bool opt=true)
void addBranch (TBranch *branch, AtomicId name=AtomicId(), bool opt=true)
template<typename T >
void addMulti (AtomicId name, const std::vector< T > *value)
template<typename T >
void addSingle (AtomicId name, const T *value, bool opt=false)
void addTypeMulti (AtomicId name, const void *value, char type)
void addTypeSingle (AtomicId name, const void *value, char type, bool opt)
void automaticAdd (bool skipTarget=false, bool skipWeight=false)
Variable::ValueList fill ()
double fill (const MVAComputer *mva)
uint64_t loop (const MVAComputer *mva)
TreeReaderoperator= (const TreeReader &orig)
void reset ()
void setOptional (AtomicId name, bool opt, double optVal=kOptVal)
void setTree (TTree *tree)
 TreeReader ()
 TreeReader (const TreeReader &orig)
 TreeReader (TTree *tree, bool skipTarget=false, bool skipWeight=false)
void update ()
std::vector< AtomicIdvariables () const
virtual ~TreeReader ()

Static Public Attributes

static const double kOptVal = -999.0

Private Attributes

std::vector< std::pair< void
*, std::vector< Bool_t > > > 
multiBool
std::vector< std::pair< void
*, std::vector< Double_t > > > 
multiDouble
std::vector< std::pair< void
*, std::vector< Float_t > > > 
multiFloat
std::vector< std::pair< void
*, std::vector< Int_t > > > 
multiInt
std::vector< BoolsingleBool
std::vector< Double_t > singleDouble
std::vector< Float_t > singleFloat
std::vector< Int_t > singleInt
TTree * tree
bool upToDate
std::map< AtomicId, ValuevalueMap
Variable::ValueList values

Friends

class Value

Detailed Description

Definition at line 19 of file TreeReader.h.


Constructor & Destructor Documentation

PhysicsTools::TreeReader::TreeReader ( )

Definition at line 25 of file TreeReader.cc.

                       :
        tree(0), upToDate(false)
{
}
PhysicsTools::TreeReader::TreeReader ( const TreeReader orig)

Definition at line 30 of file TreeReader.cc.

References operator=().

{
        this->operator = (orig);
}
PhysicsTools::TreeReader::TreeReader ( TTree *  tree,
bool  skipTarget = false,
bool  skipWeight = false 
)

Definition at line 35 of file TreeReader.cc.

References automaticAdd().

                                                                    :
        tree(tree), upToDate(false)
{
        automaticAdd(skipTarget, skipWeight);
}
PhysicsTools::TreeReader::~TreeReader ( ) [virtual]

Definition at line 41 of file TreeReader.cc.

{
}

Member Function Documentation

void PhysicsTools::TreeReader::addBranch ( const std::string &  expression,
AtomicId  name = AtomicId(),
bool  opt = true 
)

Definition at line 72 of file TreeReader.cc.

References Exception, and tree.

Referenced by automaticAdd().

{
        if (!tree)
                throw cms::Exception("NoTreeAvailable")
                        << "No TTree set in TreeReader::addBranch."
                        << std::endl;

        TBranch *branch = tree->GetBranch(expression.c_str());
        if (!branch)
                throw cms::Exception("BranchMissing")
                        << "Tree branch \"" << expression << "\" missing."
                        << std::endl;

        addBranch(branch, name, opt);
}
void PhysicsTools::TreeReader::addBranch ( TBranch *  branch,
AtomicId  name = AtomicId(),
bool  opt = true 
)

Definition at line 89 of file TreeReader.cc.

References addTypeMulti(), addTypeSingle(), Exception, python::Node::leaf, AlCaRecoCosmics_cfg::name, and valueMap.

{
        TString branchName = branch->GetName();
        if (!name)
                name = (const char*)branchName;

        TLeaf *leaf = dynamic_cast<TLeaf*>(branch->GetLeaf(branchName));
        if (!leaf)
                throw cms::Exception("InvalidBranch")
                        << "Tree branch \"" << branchName << "\" has no leaf."
                        << std::endl;

        TString typeName = leaf->GetTypeName();
        char typeId = 0;
        bool multi = false;
        if (typeName == "Double_t" || typeName == "double")
                typeId = 'D';
        else if (typeName == "Float_t" || typeName == "float")
                typeId = 'F';
        else if (typeName == "Int_t" || typeName == "int")
                typeId = 'I';
        else if (typeName == "Bool_t" || typeName == "bool")
                typeId = 'B';
        else {
                multi = true;
                if (typeName == "vector<double>" ||
                    typeName == "Vector<Double_t>")
                        typeId = 'D';
                else if (typeName == "vector<float>" ||
                         typeName == "Vector<Float_t>")
                        typeId = 'F';
                else if (typeName == "vector<int>" ||
                         typeName == "Vector<Int_t>")
                        typeId = 'I';
                else if (typeName == "vector<bool>" ||
                         typeName == "Vector<Bool_t>")
                        typeId = 'B';
        }

        if (!typeId)
                throw cms::Exception("InvalidBranch")
                        << "Tree branch \"" << branchName << "\" is of "
                           "unsupported type \"" << typeName << "\"."
                        << std::endl;

        if (multi)
                addTypeMulti(name, 0, typeId);
        else
                addTypeSingle(name, 0, typeId, opt);

        valueMap[name].setBranchName(branch->GetName());
}
template<typename T >
void PhysicsTools::TreeReader::addMulti ( AtomicId  name,
const std::vector< T > *  value 
)
template<typename T >
void PhysicsTools::TreeReader::addSingle ( AtomicId  name,
const T *  value,
bool  opt = false 
)
void PhysicsTools::TreeReader::addTypeMulti ( AtomicId  name,
const void *  value,
char  type 
)

Definition at line 199 of file TreeReader.cc.

References Exception, getHLTprescales::index, multiBool, multiDouble, multiFloat, multiInt, AlCaRecoCosmics_cfg::name, pos, upToDate, Value, and valueMap.

Referenced by addBranch().

{
        std::map<AtomicId, Value>::const_iterator pos = valueMap.find(name);
        if (pos != valueMap.end())
                throw cms::Exception("DuplicateVariable")
                        << "Duplicate Variable \"" << name << "\"."
                        << std::endl;

        if (type != 'D' && type != 'F' && type != 'I' && type != 'B')
                throw cms::Exception("InvalidType")
                        << "Unsupported type '" << type << "' in call to"
                           "TreeReader::addTypeMulti." << std::endl;

        int index = -1;
        if (!value) {
                switch(type) {
                    case 'D':
                        index = (int)multiDouble.size();
                        multiDouble.push_back(makeMulti<Double_t>());
                        break;
                    case 'F':
                        index = (int)multiFloat.size();
                        multiFloat.push_back(makeMulti<Float_t>());
                        break;
                    case 'I':
                        index = (int)multiInt.size();
                        multiInt.push_back(makeMulti<Int_t>());
                        break;
                    case 'B':
                        index = (int)multiBool.size();
                        multiBool.push_back(makeMulti<Bool_t>());
                        break;
                }
        }

        valueMap[name] = Value(index, true, false, type);
        if (value)
                valueMap[name].setPtr(value);

        upToDate = false;
}
void PhysicsTools::TreeReader::addTypeSingle ( AtomicId  name,
const void *  value,
char  type,
bool  opt 
)

Definition at line 153 of file TreeReader.cc.

References Exception, getHLTprescales::index, AlCaRecoCosmics_cfg::name, pos, singleBool, singleDouble, singleFloat, singleInt, upToDate, Value, and valueMap.

Referenced by addBranch().

{
        std::map<AtomicId, Value>::const_iterator pos = valueMap.find(name);
        if (pos != valueMap.end())
                throw cms::Exception("DuplicateVariable")
                        << "Duplicate Variable \"" << name << "\"."
                        << std::endl;

        if (type != 'D' && type != 'F' && type != 'I' && type != 'B')
                throw cms::Exception("InvalidType")
                        << "Unsupported type '" << type << "' in call to"
                           "TreeReader::addTypeSingle." << std::endl;

        int index = -1;
        if (!value) {
                switch(type) {
                    case 'D':
                        index = (int)singleDouble.size();
                        singleDouble.push_back(Double_t());
                        break;
                    case 'F':
                        index = (int)singleFloat.size();
                        singleFloat.push_back(Float_t());
                        break;
                    case 'I':
                        index = (int)singleInt.size();
                        singleInt.push_back(Int_t());
                        break;
                    case 'B':
                        index = (int)singleBool.size();
                        singleBool.push_back(Bool());
                        break;
                }
        }

        valueMap[name] = Value(index, false, opt, type);
        if (value)
                valueMap[name].setPtr(value);

        upToDate = false;
}
void PhysicsTools::TreeReader::automaticAdd ( bool  skipTarget = false,
bool  skipWeight = false 
)

Definition at line 241 of file TreeReader.cc.

References addBranch(), Exception, VarParsing::obj, and tree.

Referenced by TreeReader().

{
        if (!tree)
                throw cms::Exception("NoTreeAvailable")
                        << "No TTree set in TreeReader::automaticAdd."
                        << std::endl;

        TIter iter(tree->GetListOfBranches());
        TObject *obj;
        while((obj = iter())) {
                TBranch *branch = dynamic_cast<TBranch*>(obj);
                if (!branch)
                        continue;

                if (skipTarget &&
                    !std::strcmp(branch->GetName(), "__TARGET__"))
                        continue;

                if (skipWeight &&
                    !std::strcmp(branch->GetName(), "__WEIGHT__"))
                        continue;

                addBranch(branch);
        }
}
Variable::ValueList PhysicsTools::TreeReader::fill ( void  )

Definition at line 330 of file TreeReader.cc.

References PhysicsTools::Variable::ValueList::clear(), query::result, valueMap, and values.

Referenced by loop().

{
        for(std::map<AtomicId, Value>::const_iterator iter = valueMap.begin();
            iter != valueMap.end(); iter++)
                iter->second.fill(iter->first, this);

        Variable::ValueList result = values;
        values.clear();

        return result;
}
double PhysicsTools::TreeReader::fill ( const MVAComputer mva)

Definition at line 318 of file TreeReader.cc.

References PhysicsTools::Variable::ValueList::clear(), PhysicsTools::MVAComputer::eval(), query::result, valueMap, and values.

{
        for(std::map<AtomicId, Value>::const_iterator iter = valueMap.begin();
            iter != valueMap.end(); iter++)
                iter->second.fill(iter->first, this);

        double result = mva->eval(values);
        values.clear();

        return result;
}
uint64_t PhysicsTools::TreeReader::loop ( const MVAComputer mva)

Definition at line 298 of file TreeReader.cc.

References python::tagInventory::entries, Exception, fill(), tree, update(), and upToDate.

Referenced by PhysicsTools::TreeTrainer::iteration().

{
        if (!tree)
                throw cms::Exception("NoTreeAvailable")
                        << "No TTree set in TreeReader::automaticAdd."
                        << std::endl;

        if (!upToDate)
                update();

        Long64_t entries = tree->GetEntries();
        for(Long64_t entry = 0; entry < entries; entry++)
        {
                tree->GetEntry(entry);
                fill(mva);
        }

        return entries;
}
TreeReader & PhysicsTools::TreeReader::operator= ( const TreeReader orig)

Definition at line 45 of file TreeReader.cc.

References multiBool, multiDouble, multiFloat, multiInt, reset(), singleBool, singleDouble, singleFloat, singleInt, tree, and valueMap.

Referenced by TreeReader().

{
        reset();

        tree = orig.tree;

        multiDouble.resize(orig.multiDouble.size());
        multiFloat.resize(orig.multiFloat.size());
        multiInt.resize(orig.multiInt.size());
        multiBool.resize(orig.multiBool.size());

        singleDouble.resize(orig.singleDouble.size());
        singleFloat.resize(orig.singleFloat.size());
        singleInt.resize(orig.singleInt.size());
        singleBool.resize(orig.singleBool.size());

        valueMap = orig.valueMap;

        return *this;
}
void PhysicsTools::TreeReader::reset ( void  )

Definition at line 267 of file TreeReader.cc.

References multiBool, multiDouble, multiFloat, multiInt, singleBool, singleDouble, singleFloat, singleInt, upToDate, and valueMap.

Referenced by operator=().

{
        multiDouble.clear();
        multiFloat.clear();
        multiInt.clear();
        multiBool.clear();

        singleDouble.clear();
        singleFloat.clear();
        singleInt.clear();
        singleBool.clear();

        valueMap.clear();

        upToDate = false;
}
void PhysicsTools::TreeReader::setOptional ( AtomicId  name,
bool  opt,
double  optVal = kOptVal 
)

Definition at line 142 of file TreeReader.cc.

References pos, and valueMap.

{
        std::map<AtomicId, Value>::iterator pos = valueMap.find(name);
        if (pos == valueMap.end())
                throw cms::Exception("UnknownVariable")
                        << "Variable \"" <<name << "\" is not known to the "
                           "TreeReader." << std::endl;

        pos->second.setOpt(opt, optVal);
}
void PhysicsTools::TreeReader::setTree ( TTree *  tree)

Definition at line 66 of file TreeReader.cc.

References tree, and upToDate.

{
        this->tree = tree;
        upToDate = false;
}
void PhysicsTools::TreeReader::update ( void  )

Definition at line 284 of file TreeReader.cc.

References Exception, tree, upToDate, and valueMap.

Referenced by loop().

{
        if (!tree)
                throw cms::Exception("NoTreeAvailable")
                        << "No TTree set in TreeReader::automaticAdd."
                        << std::endl;

        for(std::map<AtomicId, Value>::iterator iter = valueMap.begin();
            iter != valueMap.end(); iter++)
                iter->second.update(this);

        upToDate = true;
}
std::vector< AtomicId > PhysicsTools::TreeReader::variables ( ) const

Definition at line 342 of file TreeReader.cc.

References query::result, and valueMap.

{
        std::vector<AtomicId> result;
        for(std::map<AtomicId, Value>::const_iterator iter = valueMap.begin();
            iter != valueMap.end(); iter++)
                result.push_back(iter->first);

        return result;
}

Friends And Related Function Documentation

friend class Value [friend]

Definition at line 106 of file TreeReader.h.

Referenced by addTypeMulti(), and addTypeSingle().


Member Data Documentation

const double PhysicsTools::TreeReader::kOptVal = -999.0 [static]

Definition at line 57 of file TreeReader.h.

std::vector<std::pair<void*, std::vector<Bool_t> > > PhysicsTools::TreeReader::multiBool [private]
std::vector<std::pair<void*, std::vector<Double_t> > > PhysicsTools::TreeReader::multiDouble [private]
std::vector<std::pair<void*, std::vector<Float_t> > > PhysicsTools::TreeReader::multiFloat [private]
std::vector<std::pair<void*, std::vector<Int_t> > > PhysicsTools::TreeReader::multiInt [private]
std::vector<Double_t> PhysicsTools::TreeReader::singleDouble [private]
std::vector<Float_t> PhysicsTools::TreeReader::singleFloat [private]
std::vector<Int_t> PhysicsTools::TreeReader::singleInt [private]

Definition at line 110 of file TreeReader.h.

Referenced by addTypeMulti(), addTypeSingle(), loop(), reset(), setTree(), and update().

Definition at line 109 of file TreeReader.h.

Referenced by fill(), and PhysicsTools::TreeReader::Value::fill().