checkmate is hosted by Hepforge, IPPP Durham
AnalysisBase Documentation  1.2.0
Check Object Conditions

Functions that remove particles checking different conditions. More...

Functions

template<class T >
std::vector< T * > AnalysisBase::filterPhaseSpace (std::vector< T * > unfiltered, double pTmin=0., double etamin=-100, double etamax=100, bool exclude_overlap=false)
 Require objects to have a certain ptmin and lie within a certain eta range. More...
 
template<class X , class Y >
std::vector< X * > AnalysisBase::overlapRemoval (std::vector< X * > candidates, std::vector< Y * > neighbours, double dR)
 Remove objects if they are to close. More...
 
template<class X >
std::vector< X * > AnalysisBase::overlapRemoval (std::vector< X * > candidates, double dR)
 Remove objects if they are to close to any other object in the same list. More...
 
std::vector< Electron * > AnalysisBase::filterIsolation (std::vector< Electron * > unfiltered, std::vector< int > relative_flags=std::vector< int >())
 Remove electrons that are not isolated. More...
 
std::vector< Electron * > AnalysisBase::filterIsolation (std::vector< Electron * > unfiltered, int relative_flag)
 Remove electrons that are not isolated (simplified function for exactly one condition given as one integer). More...
 
std::vector< Muon * > AnalysisBase::filterIsolation (std::vector< Muon * > unfiltered, std::vector< int > relative_flags=std::vector< int >())
 Remove muons that are not isolated. More...
 
std::vector< Muon * > AnalysisBase::filterIsolation (std::vector< Muon * > unfiltered, int relative_flag)
 Remove muons that are not isolated (simplified function for exactly one condition given as one integer) More...
 
std::vector< Photon * > AnalysisBase::filterIsolation (std::vector< Photon * > unfiltered, std::vector< int > relative_flags=std::vector< int >())
 Remove photons that are not isolated. More...
 
std::vector< Photon * > AnalysisBase::filterIsolation (std::vector< Photon * > unfiltered, int relative_flag)
 Remove photons that are not isolated (simplified function for exactly one condition given as one integer) More...
 
bool AnalysisBase::checkTauTag (Jet *candidate, std::string efficiency)
 Checks if candidate jet fulfills given tau identification cut. More...
 
bool AnalysisBase::checkBTag (Jet *candidate, int relative_flag=0, std::string option="")
 Checks if candidate jet fulfills given b-jet identification. More...
 

Detailed Description

Functions that remove particles checking different conditions.

These functions all take an object or a list of objects and check conditions as they appear in many analyses. Note that for functions that return vectors of objects, these vectors have the following properties:

  • The order of objects is preserved (a priori: high pt first)
  • The output does not contain copies but identically the same objects that went in.

Function Documentation

template<class T >
std::vector<T*> AnalysisBase::filterPhaseSpace ( std::vector< T * >  unfiltered,
double  pTmin = 0.,
double  etamin = -100,
double  etamax = 100,
bool  exclude_overlap = false 
)
inlineprotected

Require objects to have a certain ptmin and lie within a certain eta range.

A given set of objects (electrons, jets, ...) is filtered w.r.t minimum pt and a given eta range. the overlap region 1.37 <= |eta| <= 1.52 is common for many objects and hence it has a seperate parameter to be checked.

Parameters
unfilteredA vector containing objects of type Electron, Muon, Photon, Jet or FinalStateObject that should be filtered.
pTminThe required minimum transverse momentum of an object to pass the filter.
etaminThe required minimum pseudorapidity of an object to pass the filter.
etamaxThe required maximum pseudorapidity of an object to pass the filter.
exclude_overlapIf set to true, objects within 1.37 <= |pseudorapidity| <= 1.52 are always filtered.
Returns
A vector containing the objects that passed the filter.

Definition at line 137 of file AnalysisBase.h.

template<class X , class Y >
std::vector<X*> AnalysisBase::overlapRemoval ( std::vector< X * >  candidates,
std::vector< Y * >  neighbours,
double  dR 
)
inlineprotected

Remove objects if they are to close.

A given set X of objects (electrons, jets, ...) is compared to another set Y of objects (not necessarily of the same type). If any object in X is too close to any element in Y, it is not returned.

Parameters
candidatesA vector containing objects of type Electron, Muon, Photon, Jet or FinalStateObject that should be tested.
neighboursA vector containing objects of type Electron, Muon, Photon, Jet or FinalStateObject the candidates should be tested against. The two lists do not have to contain the same type of objects.
dRThe minimum separation a candidate must have to all neighbours in order to pass.
Returns
A vector containing candidates which all are at least dR away from all neighbours.

Definition at line 161 of file AnalysisBase.h.

template<class X >
std::vector<X*> AnalysisBase::overlapRemoval ( std::vector< X * >  candidates,
double  dR 
)
inlineprotected

Remove objects if they are to close to any other object in the same list.

Similar to overlapRemoval(std::vector<X*> candidates, std::vector<Y*> neighbours, double dR), but instead comparing the candidates to themselves. The differences are:

  • A candidate is never compared to itself
  • If a candidate overlaps with another candidate, only the one with lower pt is removed.
See Also
overlapRemoval(std::vector<X*> candidates, std::vector<Y*> neighbours, double dR)

Definition at line 191 of file AnalysisBase.h.

std::vector<Electron*> AnalysisBase::filterIsolation ( std::vector< Electron * >  unfiltered,
std::vector< int >  relative_flags = std::vector<int>() 
)
inlineprotected

Remove electrons that are not isolated.

The input list of Electrons is checked against the electron isolation criteria defined by the user in the analysis manager and only those that pass all criteria are returned.

Important If the user enters N conditions in the AnalysisManager, there will internally be N+1 conditions (i.e. the user defined plus one CheckMATE internal). This function only has access to the N conditions defined by the user, labelled with indices 0 to N-1.

Parameters
unfilteredVector of electrons that should be tested
relative_flagsVector of integers of the flags to be tested. The first condition entered into the AnalysisManager corresponds to 0, the second (if it exists) to 1, etc. If no second parameter is provided, all objects are compared to all existing conditions.

Definition at line 227 of file AnalysisBase.h.

std::vector<Electron*> AnalysisBase::filterIsolation ( std::vector< Electron * >  unfiltered,
int  relative_flag 
)
inlineprotected

Remove electrons that are not isolated (simplified function for exactly one condition given as one integer).

See Also
filterIsolation(std::vector<Electron*> unfiltered, std::vector<int> relative_flags = std::vector<int>())

Definition at line 246 of file AnalysisBase.h.

std::vector<Muon*> AnalysisBase::filterIsolation ( std::vector< Muon * >  unfiltered,
std::vector< int >  relative_flags = std::vector<int>() 
)
inlineprotected

Remove muons that are not isolated.

See Also
filterIsolation(std::vector<Electron*> unfiltered, std::vector<int> relative_flags = std::vector<int>())

Definition at line 258 of file AnalysisBase.h.

std::vector<Muon*> AnalysisBase::filterIsolation ( std::vector< Muon * >  unfiltered,
int  relative_flag 
)
inlineprotected

Remove muons that are not isolated (simplified function for exactly one condition given as one integer)

See Also
filterIsolation(std::vector<Electron*> unfiltered, int relative_flag)

Definition at line 277 of file AnalysisBase.h.

std::vector<Photon*> AnalysisBase::filterIsolation ( std::vector< Photon * >  unfiltered,
std::vector< int >  relative_flags = std::vector<int>() 
)
inlineprotected

Remove photons that are not isolated.

See Also
filterIsolation(std::vector<Electron*> unfiltered, std::vector<int> relative_flags = std::vector<int>())

Definition at line 291 of file AnalysisBase.h.

std::vector<Photon*> AnalysisBase::filterIsolation ( std::vector< Photon * >  unfiltered,
int  relative_flag 
)
inlineprotected

Remove photons that are not isolated (simplified function for exactly one condition given as one integer)

See Also
filterIsolation(std::vector<Electron*> unfiltered, int relative_flag)

Definition at line 310 of file AnalysisBase.h.

bool AnalysisBase::checkTauTag ( Jet candidate,
std::string  efficiency 
)
protected

Checks if candidate jet fulfills given tau identification cut.

If tau tagging was activated in the AnalysisManager, a given jet candidate can be tested for tau tags according to the three main working points 'loose', 'medium' and 'tight'. Note that every 'tight' is always also a 'medium' tau jet and similarly every 'medium' is always also 'loose'.

Parameters
candidateThe jet candidate to be tested.
efficiencyThe to be tested efficiency 'loose', 'medium' or 'tight'.
Returns
'true' if the candidate was tagged as a tau jet, otherwise 'false.

Definition at line 258 of file AnalysisBase.cc.

bool AnalysisBase::checkBTag ( Jet candidate,
int  relative_flag = 0,
std::string  option = "" 
)
protected

Checks if candidate jet fulfills given b-jet identification.

If b tagging was activated in the AnalysisManager, a given jet candidate can be tested for b-jets giving the respective index of the defined working point, where the first b-efficiency defined in CheckMATE is tested with index 0, the second one with 1 etc. If the candidate belongs to the 'second jet' type, the user has to give the option 'secondJet' and consequently can only use the btags that have been defined for the 2nd jet in the AnalysisManager. The tags are set such that each jet that passed the working point with signal efficiency X will always also pass those with working points with effiencies larger than X. (E.g. every 40% jet will always be a 60% jet too).

Parameters
candidateThe jet candidate to be tested.
relative_flagThe condition to be tested with 0 being the first defined in the AnalysisManager.
optionIf option=="secondJet", the btags defined for the second jet type are checked.
Returns
'true' if the candidate was tagged as a b jet, otherwise 'false.

Definition at line 269 of file AnalysisBase.cc.