Class for the calculation of several event shape variables. More...
#include "PhysicsTools/CandUtils/interface/EventShapeVariables.h"
Public Member Functions | |
double | aplanarity (double=2.) const |
double | C (double=2.) const |
double | circularity (const unsigned int &numberOfSteps=1000) const |
double | D (double=2.) const |
EventShapeVariables (const edm::View< reco::Candidate > &inputVectors) | |
constructor from reco::Candidates More... | |
EventShapeVariables (const std::vector< math::XYZVector > &inputVectors) | |
constructor from XYZ coordinates More... | |
EventShapeVariables (const std::vector< math::RhoEtaPhiVector > &inputVectors) | |
constructor from rho eta phi coordinates More... | |
EventShapeVariables (const std::vector< math::RThetaPhiVector > &inputVectors) | |
constructor from r theta phi coordinates More... | |
double | isotropy (const unsigned int &numberOfSteps=1000) const |
double | sphericity (double=2.) const |
~EventShapeVariables () | |
default destructor More... | |
Private Member Functions | |
TVectorD | compEigenValues (double=2.) const |
TMatrixDSym | compMomentumTensor (double=2.) const |
helper function to fill the 3 dimensional momentum tensor from the inputVecotrs where needed More... | |
Private Attributes | |
std::vector< math::XYZVector > | inputVectors_ |
cashing of input vectors More... | |
Class for the calculation of several event shape variables.
Class for the calculation of several event shape variables. Isotropy, sphericity, aplanarity and circularity are supported. The class supports vectors of 3d vectors and edm::Views of reco::Candidates as input. The 3d vectors can be given in cartesian, cylindrical or polar coordinates. It exploits the ROOT::TMatrixDSym for the calculation of the sphericity and aplanarity.
See http://cepa.fnal.gov/psm/simulation/mcgen/lund/pythia_manual/pythia6.3/pythia6301/node213.html for an explanation of sphericity, aplanarity and the quantities C and D.
Author: Sebastian Naumann-Emme, University of Hamburg Roger Wolf, University of Hamburg Christian Veelken, UC Davis
Definition at line 32 of file EventShapeVariables.h.
|
explicit |
constructor from reco::Candidates
Definition at line 6 of file EventShapeVariables.cc.
References edm::View< T >::begin(), edm::View< T >::end(), inputVectors_, and edm::View< T >::size().
|
explicit |
constructor from XYZ coordinates
Definition at line 16 of file EventShapeVariables.cc.
|
explicit |
constructor from rho eta phi coordinates
Definition at line 21 of file EventShapeVariables.cc.
References inputVectors_.
|
explicit |
constructor from r theta phi coordinates
Definition at line 30 of file EventShapeVariables.cc.
References inputVectors_.
|
inline |
default destructor
Definition at line 44 of file EventShapeVariables.h.
References aplanarity(), C(), circularity(), compEigenValues(), compMomentumTensor(), D(), isotropy(), and sphericity().
double EventShapeVariables::aplanarity | ( | double | r = 2. | ) | const |
1.5*q1 where 0<=q1<=q2<=q3 are the eigenvalues of the momemtum tensor sum{p_j[a]*p_j[b]}/sum{p_j**2} normalized to 1. Return values are 0.5 for spherical and 0 for plane and linear events
1.5*q1 where 0<=q1<=q2<=q3 are the eigenvalues of the momentum tensor sum{p_j[a]*p_j[b]}/sum{p_j**2} normalized to 1. Return values are 0.5 for spherical and 0 for plane and linear events
Definition at line 153 of file EventShapeVariables.cc.
References compEigenValues().
Referenced by TtFullHadSignalSel::TtFullHadSignalSel(), and ~EventShapeVariables().
double EventShapeVariables::C | ( | double | r = 2. | ) | const |
3.*(q1*q2+q1*q3+q2*q3) where 0<=q1<=q2<=q3 are the eigenvalues of the momemtum tensor sum{p_j[a]*p_j[b]}/sum{p_j**2} normalized to 1. Return value is between 0 and 1 and measures the 3-jet structure of the event (C vanishes for a "perfect" 2-jet event)
3.*(q1*q2+q1*q3+q2*q3) where 0<=q1<=q2<=q3 are the eigenvalues of the momentum tensor sum{p_j[a]*p_j[b]}/sum{p_j**2} normalized to 1. Return value is between 0 and 1 and measures the 3-jet structure of the event (C vanishes for a "perfect" 2-jet event)
Definition at line 163 of file EventShapeVariables.cc.
References compEigenValues().
Referenced by TtFullHadSignalSel::TtFullHadSignalSel(), and ~EventShapeVariables().
double EventShapeVariables::circularity | ( | const unsigned int & | numberOfSteps = 1000 | ) | const |
the return value is 1 for spherical and 0 linear events in r-phi. This function needs the number of steps to determine how fine the granularity of the algorithm in phi should be
Definition at line 62 of file EventShapeVariables.cc.
References Abs(), MuonCkfTrajectoryBuilder_cfi::deltaPhi, mps_fire::i, inputVectors_, phi, Pi, tmp, x, and y.
Referenced by TtFullHadSignalSel::TtFullHadSignalSel(), and ~EventShapeVariables().
|
private |
helper function to fill the 3 dimensional vector of eigen-values; the largest (smallest) eigen-value is stored at index position 0 (2)
Definition at line 124 of file EventShapeVariables.cc.
References compMomentumTensor().
Referenced by aplanarity(), C(), D(), sphericity(), and ~EventShapeVariables().
|
private |
helper function to fill the 3 dimensional momentum tensor from the inputVecotrs where needed
helper function to fill the 3 dimensional momentum tensor from the inputVectors where needed
Definition at line 85 of file EventShapeVariables.cc.
References mps_fire::i, inputVectors_, createfilelist::int, and p2.
Referenced by compEigenValues(), and ~EventShapeVariables().
double EventShapeVariables::D | ( | double | r = 2. | ) | const |
27.*(q1*q2*q3) where 0<=q1<=q2<=q3 are the eigenvalues of the momemtum tensor sum{p_j[a]*p_j[b]}/sum{p_j**2} normalized to 1. Return value is between 0 and 1 and measures the 4-jet structure of the event (D vanishes for a planar event)
Definition at line 173 of file EventShapeVariables.cc.
References compEigenValues().
Referenced by TtFullHadSignalSel::TtFullHadSignalSel(), and ~EventShapeVariables().
double EventShapeVariables::isotropy | ( | const unsigned int & | numberOfSteps = 1000 | ) | const |
the return value is 1 for spherical events and 0 for events linear in r-phi. This function needs the number of steps to determine how fine the granularity of the algorithm in phi should be
Definition at line 42 of file EventShapeVariables.cc.
References Abs(), MuonCkfTrajectoryBuilder_cfi::deltaPhi, mps_fire::i, inputVectors_, phi, Pi, x, and y.
Referenced by TtFullHadSignalSel::TtFullHadSignalSel(), and ~EventShapeVariables().
double EventShapeVariables::sphericity | ( | double | r = 2. | ) | const |
1.5*(q1+q2) where 0<=q1<=q2<=q3 are the eigenvalues of the momemtum tensor sum{p_j[a]*p_j[b]}/sum{p_j**2} normalized to 1. Return values are 1 for spherical, 3/4 for plane and 0 for linear events
1.5*(q1+q2) where 0<=q1<=q2<=q3 are the eigenvalues of the momentum tensor sum{p_j[a]*p_j[b]}/sum{p_j**2} normalized to 1. Return values are 1 for spherical, 3/4 for plane and 0 for linear events
Definition at line 144 of file EventShapeVariables.cc.
References compEigenValues().
Referenced by TtFullHadSignalSel::TtFullHadSignalSel(), and ~EventShapeVariables().
|
private |
cashing of input vectors
Definition at line 79 of file EventShapeVariables.h.
Referenced by circularity(), compMomentumTensor(), EventShapeVariables(), and isotropy().