#include <AlphaT.h>
Public Member Functions | |
template<class T > | |
AlphaT (std::vector< T const * > const &p4, bool use_et=true) | |
template<class T > | |
AlphaT (std::vector< T > const &p4, bool use_et=true) | |
double | value (void) const |
double | value (std::vector< bool > &jet_sign) const |
Private Member Functions | |
double | value_ (std::vector< bool > *jet_sign) const |
Private Attributes | |
std::vector< double > | et_ |
std::vector< double > | px_ |
std::vector< double > | py_ |
AlphaT::AlphaT | ( | std::vector< T const * > const & | p4, |
bool | use_et = true |
||
) |
Definition at line 33 of file AlphaT.h.
References et_, reco::tau::disc::Pt(), px_, py_, and create_public_pileup_plots::transform.
{ std::transform( p4.begin(), p4.end(), back_inserter(et_), ( use_et ? std::mem_fun(&T::Et) : std::mem_fun(&T::Pt) ) ); std::transform( p4.begin(), p4.end(), back_inserter(px_), std::mem_fun(&T::Px) ); std::transform( p4.begin(), p4.end(), back_inserter(py_), std::mem_fun(&T::Py) ); }
AlphaT::AlphaT | ( | std::vector< T > const & | p4, |
bool | use_et = true |
||
) |
Definition at line 41 of file AlphaT.h.
References et_, reco::tau::disc::Pt(), px_, py_, and create_public_pileup_plots::transform.
{ std::transform( p4.begin(), p4.end(), back_inserter(et_), std::mem_fun_ref( use_et ? &T::Et : &T::Pt ) ); std::transform( p4.begin(), p4.end(), back_inserter(px_), std::mem_fun_ref(&T::Px) ); std::transform( p4.begin(), p4.end(), back_inserter(py_), std::mem_fun_ref(&T::Py) ); }
double AlphaT::value | ( | void | ) | const [inline] |
Definition at line 50 of file AlphaT.h.
References value_().
Referenced by HLTAlphaTFilter< T >::hltFilter().
{ return value_(0); }
double AlphaT::value | ( | std::vector< bool > & | jet_sign | ) | const [inline] |
double AlphaT::value_ | ( | std::vector< bool > * | jet_sign | ) | const [private] |
Definition at line 3 of file AlphaT.cc.
References abs, et_, i, j, max(), px_, py_, and mathSSE::sqrt().
Referenced by value().
{ // Clear pseudo-jet container if (jet_sign) { jet_sign->clear(); jet_sign->resize(et_.size()); } // check the size of the input collection if (et_.size() == 0) // empty jet collection, return AlphaT = 0 return 0.; if (et_.size() > (unsigned int) std::numeric_limits<unsigned int>::digits) // too many jets, return AlphaT = a very large number return std::numeric_limits<double>::max(); // Momentum sums in transverse plane const double sum_et = std::accumulate( et_.begin(), et_.end(), 0. ); const double sum_px = std::accumulate( px_.begin(), px_.end(), 0. ); const double sum_py = std::accumulate( py_.begin(), py_.end(), 0. ); // Minimum Delta Et for two pseudo-jets double min_delta_sum_et = sum_et; for (unsigned int i = 0; i < (1U << (et_.size() - 1)); i++) { //@@ iterate through different combinations double delta_sum_et = 0.; for (unsigned int j = 0; j < et_.size(); ++j) { //@@ iterate through jets if (i & (1U << j)) delta_sum_et -= et_[j]; else delta_sum_et += et_[j]; } delta_sum_et = std::abs(delta_sum_et); if (delta_sum_et < min_delta_sum_et) { min_delta_sum_et = delta_sum_et; if (jet_sign) { for (unsigned int j = 0; j < et_.size(); ++j) (*jet_sign)[j] = ((i & (1U << j)) == 0); } } } // Alpha_T return (0.5 * (sum_et - min_delta_sum_et) / sqrt( sum_et*sum_et - (sum_px*sum_px+sum_py*sum_py) )); }
std::vector<double> AlphaT::et_ [private] |
std::vector<double> AlphaT::px_ [private] |
std::vector<double> AlphaT::py_ [private] |