CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
EventShapeVariables Class Reference

Class for the calculation of several event shape variables. More...

#include "PhysicsTools/CandUtils/interface/EventShapeVariables.h"

Public Member Functions

double aplanarity ()
 
double C ()
 
double circularity (const unsigned int &numberOfSteps=1000) const
 
double D ()
 
 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...
 
const std::vector< double > & getEigenValues ()
 
const std::vector< double > & getEigenValuesNoNorm ()
 
const TMatrixD & getEigenVectors ()
 
double getFWmoment (unsigned l)
 
const std::vector< double > & getFWmoments ()
 
double isotropy (const unsigned int &numberOfSteps=1000) const
 
void set_r (double r)
 set exponent for computation of momentum tensor and related products More...
 
void setFWmax (unsigned m)
 set number of Fox-Wolfram moments to compute More...
 
double sphericity ()
 
 ~EventShapeVariables ()
 default destructor More...
 

Private Member Functions

void compTensorsAndVectors ()
 
void computeFWmoments ()
 

Private Attributes

std::vector< double > eigenValues_
 
std::vector< double > eigenValuesNoNorm_
 
TVectorD eigenValuesNoNormTmp_
 
TVectorD eigenValuesTmp_
 
TMatrixD eigenVectors_
 
std::vector< double > fwmom_
 
bool fwmom_computed_
 
unsigned fwmom_maxl_
 Owen ; save computed Fox-Wolfram moments. More...
 
std::vector< math::XYZVectorinputVectors_
 caching of input vectors More...
 
double r_
 caching of output More...
 
bool tensors_computed_
 

Detailed Description

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 https://arxiv.org/pdf/hep-ph/0603175v2.pdf#page=524 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.

Constructor & Destructor Documentation

EventShapeVariables::EventShapeVariables ( const edm::View< reco::Candidate > &  inputVectors)
explicit

constructor from reco::Candidates

Definition at line 6 of file EventShapeVariables.cc.

References inputVectors_, set_r(), setFWmax(), and edm::View< T >::size().

6  : eigenVectors_(3,3)
7 {
8  inputVectors_.reserve( inputVectors.size() );
9  for (const auto& vec : inputVectors){
10  inputVectors_.push_back(math::XYZVector(vec.px(), vec.py(), vec.pz()));
11  }
12  //default values
13  set_r(2.);
14  setFWmax(10);
15 }
size_type size() const
void setFWmax(unsigned m)
set number of Fox-Wolfram moments to compute
std::vector< math::XYZVector > inputVectors_
caching of input vectors
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
void set_r(double r)
set exponent for computation of momentum tensor and related products
EventShapeVariables::EventShapeVariables ( const std::vector< math::XYZVector > &  inputVectors)
explicit

constructor from XYZ coordinates

Definition at line 19 of file EventShapeVariables.cc.

References set_r(), and setFWmax().

19  : inputVectors_(inputVectors), eigenVectors_(3,3)
20 {
21  //default values
22  set_r(2.);
23  setFWmax(10);
24 }
void setFWmax(unsigned m)
set number of Fox-Wolfram moments to compute
std::vector< math::XYZVector > inputVectors_
caching of input vectors
void set_r(double r)
set exponent for computation of momentum tensor and related products
EventShapeVariables::EventShapeVariables ( const std::vector< math::RhoEtaPhiVector > &  inputVectors)
explicit

constructor from rho eta phi coordinates

Definition at line 27 of file EventShapeVariables.cc.

References inputVectors_, set_r(), and setFWmax().

27  : eigenVectors_(3,3)
28 {
29  inputVectors_.reserve( inputVectors.size() );
30  for (const auto& vec : inputVectors){
31  inputVectors_.push_back(math::XYZVector(vec.x(), vec.y(), vec.z()));
32  }
33  //default values
34  set_r(2.);
35  setFWmax(10);
36 }
void setFWmax(unsigned m)
set number of Fox-Wolfram moments to compute
std::vector< math::XYZVector > inputVectors_
caching of input vectors
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
void set_r(double r)
set exponent for computation of momentum tensor and related products
EventShapeVariables::EventShapeVariables ( const std::vector< math::RThetaPhiVector > &  inputVectors)
explicit

constructor from r theta phi coordinates

Definition at line 39 of file EventShapeVariables.cc.

References inputVectors_, set_r(), and setFWmax().

39  : eigenVectors_(3,3)
40 {
41  inputVectors_.reserve( inputVectors.size() );
42  for (const auto& vec : inputVectors){
43  inputVectors_.push_back(math::XYZVector(vec.x(), vec.y(), vec.z()));
44  }
45  //default values
46  set_r(2.);
47  setFWmax(10);
48 }
void setFWmax(unsigned m)
set number of Fox-Wolfram moments to compute
std::vector< math::XYZVector > inputVectors_
caching of input vectors
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
void set_r(double r)
set exponent for computation of momentum tensor and related products
EventShapeVariables::~EventShapeVariables ( )
inline

default destructor

Definition at line 44 of file EventShapeVariables.h.

References aplanarity(), C(), circularity(), D(), isotropy(), funct::m, alignCSCRings::r, set_r(), setFWmax(), and sphericity().

44 {};

Member Function Documentation

double EventShapeVariables::aplanarity ( )

1.5*q2 where q0>=q1>=q2>=0 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 177 of file EventShapeVariables.cc.

References compTensorsAndVectors(), eigenValues_, and tensors_computed_.

Referenced by TtFullHadSignalSel::TtFullHadSignalSel(), and ~EventShapeVariables().

178 {
180  return 1.5*eigenValues_[2];
181 }
std::vector< double > eigenValues_
double EventShapeVariables::C ( )

3.*(q0*q1+q0*q2+q1*q2) where q0>=q1>=q2>=0 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 187 of file EventShapeVariables.cc.

References compTensorsAndVectors(), eigenValues_, and tensors_computed_.

Referenced by TtFullHadSignalSel::TtFullHadSignalSel(), and ~EventShapeVariables().

188 {
191 }
std::vector< double > eigenValues_
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 76 of file EventShapeVariables.cc.

References Abs(), custom_jme_cff::area, hiPixelPairStep_cff::deltaPhi, mps_fire::i, inputVectors_, phi, Pi, and tmp.

Referenced by TtFullHadSignalSel::TtFullHadSignalSel(), and ~EventShapeVariables().

77 {
78  const double deltaPhi=2*TMath::Pi()/numberOfSteps;
79  double circularity=-1, phi=0, area = 0;
80  for(const auto& vec: inputVectors_) {
81  area+=TMath::Sqrt(vec.x()*vec.x()+vec.y()*vec.y());
82  }
83  for(unsigned int i=0; i<numberOfSteps; ++i){
84  phi+=deltaPhi;
85  double sum=0, tmp=0.;
86  double cosphi = TMath::Cos(phi);
87  double sinphi = TMath::Sin(phi);
88  for(const auto& vec: inputVectors_){
89  sum+=TMath::Abs(cosphi*vec.x()+sinphi*vec.y());
90  }
91  tmp=TMath::Pi()/2*sum/area;
92  if( circularity<0 || tmp<circularity ){
93  circularity=tmp;
94  }
95  }
96  return circularity;
97 }
const double Pi
double circularity(const unsigned int &numberOfSteps=1000) const
T Abs(T a)
Definition: MathUtil.h:49
std::vector< math::XYZVector > inputVectors_
caching of input vectors
std::vector< std::vector< double > > tmp
Definition: MVATrainer.cc:100
void EventShapeVariables::compTensorsAndVectors ( )
private

helper function to fill the 3 dimensional momentum tensor from the inputVectors where needed also fill the 3 dimensional vectors of eigen-values and eigen-vectors; the largest (smallest) eigen-value is stored at index position 0 (2)

Definition at line 114 of file EventShapeVariables.cc.

References eigenValues_, eigenValuesNoNorm_, eigenValuesNoNormTmp_, eigenValuesTmp_, eigenVectors_, inputVectors_, p2, r_, and tensors_computed_.

Referenced by aplanarity(), C(), D(), getEigenValues(), getEigenValuesNoNorm(), getEigenVectors(), and sphericity().

115 {
116  if(tensors_computed_) return;
117 
118  if ( inputVectors_.size() < 2 ){
119  tensors_computed_ = true;
120  return;
121  }
122 
123  TMatrixDSym momentumTensor(3);
124  momentumTensor.Zero();
125 
126  // fill momentumTensor from inputVectors
127  double norm = 0.;
128  for (const auto& vec: inputVectors_){
129  double p2 = vec.Dot(vec);
130  double pR = ( r_ == 2. ) ? p2 : TMath::Power(p2, 0.5*r_);
131  norm += pR;
132  double pRminus2 = ( r_ == 2. ) ? 1. : TMath::Power(p2, 0.5*r_ - 1.);
133  momentumTensor(0,0) += pRminus2*vec.x()*vec.x();
134  momentumTensor(0,1) += pRminus2*vec.x()*vec.y();
135  momentumTensor(0,2) += pRminus2*vec.x()*vec.z();
136  momentumTensor(1,0) += pRminus2*vec.y()*vec.x();
137  momentumTensor(1,1) += pRminus2*vec.y()*vec.y();
138  momentumTensor(1,2) += pRminus2*vec.y()*vec.z();
139  momentumTensor(2,0) += pRminus2*vec.z()*vec.x();
140  momentumTensor(2,1) += pRminus2*vec.z()*vec.y();
141  momentumTensor(2,2) += pRminus2*vec.z()*vec.z();
142  }
143 
144  if( momentumTensor.IsSymmetric() && ( momentumTensor.NonZeros() != 0 ) ){
145  momentumTensor.EigenVectors(eigenValuesNoNormTmp_);
146  }
150 
151  // momentumTensor normalized to determinant 1
152  momentumTensor *= (1./norm);
153 
154  // now get eigens
155  if( momentumTensor.IsSymmetric() && ( momentumTensor.NonZeros() != 0 ) ){
156  eigenVectors_ = momentumTensor.EigenVectors(eigenValuesTmp_);
157  }
161 
162  tensors_computed_ = true;
163 }
std::vector< double > eigenValuesNoNorm_
std::vector< double > eigenValues_
double r_
caching of output
double p2[4]
Definition: TauolaWrapper.h:90
std::vector< math::XYZVector > inputVectors_
caching of input vectors
void EventShapeVariables::computeFWmoments ( )
private

reduce computation by exploiting symmetry: all off-diagonal elements appear twice in the sum

compute higher legendre polynomials recursively need to keep track of two previous values

initial cases

store new value

Definition at line 230 of file EventShapeVariables.cc.

References fwmom_, fwmom_computed_, fwmom_maxl_, mps_fire::i, inputVectors_, and gen::n.

Referenced by getEigenVectors(), getFWmoment(), and getFWmoments().

230  {
231  if(fwmom_computed_) return;
232 
233  double esum_total(0.) ;
234  for ( unsigned int i = 0 ; i < inputVectors_.size() ; i ++ ) {
235  esum_total += inputVectors_[i].R() ;
236  } // i
237  double esum_total_sq = esum_total * esum_total ;
238 
239  for ( unsigned int i = 0 ; i < inputVectors_.size() ; i ++ ) {
240 
241  double p_i = inputVectors_[i].R() ;
242  if ( p_i <= 0 ) continue ;
243 
244  for ( unsigned int j = 0 ; j <= i ; j ++ ) {
245 
246  double p_j = inputVectors_[j].R() ;
247  if ( p_j <= 0 ) continue ;
248 
251  int symmetry_factor = 2;
252  if( j == i ) symmetry_factor = 1;
253  double p_ij = p_i * p_j;
254  double cosTheta = inputVectors_[i].Dot( inputVectors_[j] ) / (p_ij) ;
255  double pi_pj_over_etot2 = p_ij / esum_total_sq ;
256 
259  double Pn1 = 0;
260  double Pn2 = 0;
261  for ( unsigned n = 0; n < fwmom_maxl_; n ++ ) {
263  if ( n == 0 ) {
264  Pn2 = pi_pj_over_etot2;
265  fwmom_[0] += Pn2*symmetry_factor;
266  }
267  else if ( n == 1 ) {
268  Pn1 = pi_pj_over_etot2 * cosTheta;
269  fwmom_[1] += Pn1*symmetry_factor;
270  }
271  else {
272  double Pn = ((2*n-1)*cosTheta*Pn1 - (n-1)*Pn2)/n;
273  fwmom_[n] += Pn*symmetry_factor;
275  Pn2 = Pn1;
276  Pn1 = Pn;
277  }
278  }
279 
280  } // j
281  } // i
282 
283  fwmom_computed_ = true;
284 
285 } // computeFWmoments
std::vector< math::XYZVector > inputVectors_
caching of input vectors
unsigned fwmom_maxl_
Owen ; save computed Fox-Wolfram moments.
std::vector< double > fwmom_
double EventShapeVariables::D ( )

27.*(q0*q1*q2) where q0>=q1>=q2>=0 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 4-jet structure of the event (D vanishes for a planar event)

27.*(q0*q1*q2) where q0>=q1>=q2>=0 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 197 of file EventShapeVariables.cc.

References compTensorsAndVectors(), eigenValues_, and tensors_computed_.

Referenced by TtFullHadSignalSel::TtFullHadSignalSel(), and ~EventShapeVariables().

198 {
200  return 27.*eigenValues_[0]*eigenValues_[1]*eigenValues_[2];
201 }
std::vector< double > eigenValues_
const std::vector<double>& EventShapeVariables::getEigenValues ( )
inline
const std::vector<double>& EventShapeVariables::getEigenValuesNoNorm ( )
inline
const TMatrixD& EventShapeVariables::getEigenVectors ( )
inline
double EventShapeVariables::getFWmoment ( unsigned  l)

Definition at line 214 of file EventShapeVariables.cc.

References computeFWmoments(), fwmom_, fwmom_computed_, fwmom_maxl_, and checklumidiff::l.

Referenced by getEigenVectors().

214  {
215 
216  if ( l > fwmom_maxl_ ) return 0. ;
217 
219 
220  return fwmom_[l] ;
221 
222 } // getFWmoment
unsigned fwmom_maxl_
Owen ; save computed Fox-Wolfram moments.
std::vector< double > fwmom_
const std::vector< double > & EventShapeVariables::getFWmoments ( )

Definition at line 224 of file EventShapeVariables.cc.

References computeFWmoments(), fwmom_, and fwmom_computed_.

Referenced by getEigenVectors().

224  {
226 
227  return fwmom_;
228 }
std::vector< double > fwmom_
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 54 of file EventShapeVariables.cc.

References Abs(), hiPixelPairStep_cff::deltaPhi, mps_fire::i, inputVectors_, phi, and Pi.

Referenced by TtFullHadSignalSel::TtFullHadSignalSel(), and ~EventShapeVariables().

55 {
56  const double deltaPhi=2*TMath::Pi()/numberOfSteps;
57  double phi = 0, eIn =-1., eOut=-1.;
58  for(unsigned int i=0; i<numberOfSteps; ++i){
59  phi+=deltaPhi;
60  double sum=0;
61  double cosphi = TMath::Cos(phi);
62  double sinphi = TMath::Sin(phi);
63  for(const auto& vec : inputVectors_){
64  // sum over inner product of unit vectors and momenta
65  sum+=TMath::Abs(cosphi*vec.x()+sinphi*vec.y());
66  }
67  if( eOut<0. || sum<eOut ) eOut=sum;
68  if( eIn <0. || sum>eIn ) eIn =sum;
69  }
70  return (eIn-eOut)/eIn;
71 }
const double Pi
T Abs(T a)
Definition: MathUtil.h:49
std::vector< math::XYZVector > inputVectors_
caching of input vectors
void EventShapeVariables::set_r ( double  r)

set exponent for computation of momentum tensor and related products

invalidate previous cached computations

Definition at line 101 of file EventShapeVariables.cc.

References eigenValues_, eigenValuesNoNorm_, alignCSCRings::r, r_, and tensors_computed_.

Referenced by EventShapeVariables(), and ~EventShapeVariables().

102 {
103  r_ = r;
105  tensors_computed_ = false;
106  eigenValues_ = std::vector<double>(3,0);
107  eigenValuesNoNorm_ = std::vector<double>(3,0);
108 }
std::vector< double > eigenValuesNoNorm_
std::vector< double > eigenValues_
double r_
caching of output
void EventShapeVariables::setFWmax ( unsigned  m)

set number of Fox-Wolfram moments to compute

Definition at line 207 of file EventShapeVariables.cc.

References fwmom_, fwmom_computed_, fwmom_maxl_, and funct::m.

Referenced by EventShapeVariables(), and ~EventShapeVariables().

208 {
209  fwmom_maxl_ = m;
210  fwmom_computed_ = false;
211  fwmom_ = std::vector<double>(fwmom_maxl_,0.);
212 }
unsigned fwmom_maxl_
Owen ; save computed Fox-Wolfram moments.
std::vector< double > fwmom_
double EventShapeVariables::sphericity ( )

1.5*(q1+q2) where q0>=q1>=q2>=0 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 168 of file EventShapeVariables.cc.

References compTensorsAndVectors(), eigenValues_, and tensors_computed_.

Referenced by TtFullHadSignalSel::TtFullHadSignalSel(), and ~EventShapeVariables().

169 {
171  return 1.5*(eigenValues_[1] + eigenValues_[2]);
172 }
std::vector< double > eigenValues_

Member Data Documentation

std::vector<double> EventShapeVariables::eigenValues_
private
std::vector<double> EventShapeVariables::eigenValuesNoNorm_
private

Definition at line 101 of file EventShapeVariables.h.

Referenced by compTensorsAndVectors(), getEigenValuesNoNorm(), and set_r().

TVectorD EventShapeVariables::eigenValuesNoNormTmp_
private

Definition at line 100 of file EventShapeVariables.h.

Referenced by compTensorsAndVectors().

TVectorD EventShapeVariables::eigenValuesTmp_
private

Definition at line 100 of file EventShapeVariables.h.

Referenced by compTensorsAndVectors().

TMatrixD EventShapeVariables::eigenVectors_
private

Definition at line 99 of file EventShapeVariables.h.

Referenced by compTensorsAndVectors(), and getEigenVectors().

std::vector<double> EventShapeVariables::fwmom_
private

Definition at line 105 of file EventShapeVariables.h.

Referenced by computeFWmoments(), getFWmoment(), getFWmoments(), and setFWmax().

bool EventShapeVariables::fwmom_computed_
private

Definition at line 106 of file EventShapeVariables.h.

Referenced by computeFWmoments(), getFWmoment(), getFWmoments(), and setFWmax().

unsigned EventShapeVariables::fwmom_maxl_
private

Owen ; save computed Fox-Wolfram moments.

Definition at line 104 of file EventShapeVariables.h.

Referenced by computeFWmoments(), getFWmoment(), and setFWmax().

std::vector<math::XYZVector> EventShapeVariables::inputVectors_
private

caching of input vectors

Definition at line 94 of file EventShapeVariables.h.

Referenced by circularity(), compTensorsAndVectors(), computeFWmoments(), EventShapeVariables(), and isotropy().

double EventShapeVariables::r_
private

caching of output

Definition at line 97 of file EventShapeVariables.h.

Referenced by compTensorsAndVectors(), and set_r().

bool EventShapeVariables::tensors_computed_
private