CMS 3D CMS Logo

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>

List of all members.

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 std::vector< math::RThetaPhiVector > &inputVectors)
 constructor from r theta phi coordinates
 EventShapeVariables (const std::vector< math::RhoEtaPhiVector > &inputVectors)
 constructor from rho eta phi coordinates
 EventShapeVariables (const std::vector< math::XYZVector > &inputVectors)
 constructor from XYZ coordinates
 EventShapeVariables (const edm::View< reco::Candidate > &inputVectors)
 constructor from reco::Candidates
double isotropy (const unsigned int &numberOfSteps=1000) const
double sphericity (double=2.) const
 ~EventShapeVariables ()
 default destructor

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

Private Attributes

std::vector< math::XYZVectorinputVectors_
 cashing of input vectors

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 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.


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 edm::View< T >::begin(), gather_cfg::cout, edm::View< T >::end(), inputVectors_, and edm::View< T >::size().

{
  std::cout << "inputVectors.size = " << inputVectors.size() << std::endl;
  inputVectors_.reserve( inputVectors.size() );
  for ( edm::View<reco::Candidate>::const_iterator vec = inputVectors.begin(); vec != inputVectors.end(); ++vec){
    inputVectors_.push_back(math::XYZVector(vec->px(), vec->py(), vec->pz()));
  }
}
EventShapeVariables::EventShapeVariables ( const std::vector< math::XYZVector > &  inputVectors) [explicit]

constructor from XYZ coordinates

Definition at line 16 of file EventShapeVariables.cc.

  : inputVectors_(inputVectors)
{}
EventShapeVariables::EventShapeVariables ( const std::vector< math::RhoEtaPhiVector > &  inputVectors) [explicit]

constructor from rho eta phi coordinates

Definition at line 21 of file EventShapeVariables.cc.

References inputVectors_.

{
  inputVectors_.reserve( inputVectors.size() );
  for ( std::vector<math::RhoEtaPhiVector>::const_iterator vec = inputVectors.begin(); vec != inputVectors.end(); ++vec ){
    inputVectors_.push_back(math::XYZVector(vec->x(), vec->y(), vec->z()));
  }
}
EventShapeVariables::EventShapeVariables ( const std::vector< math::RThetaPhiVector > &  inputVectors) [explicit]

constructor from r theta phi coordinates

Definition at line 30 of file EventShapeVariables.cc.

References inputVectors_.

{
  inputVectors_.reserve( inputVectors.size() );
  for(std::vector<math::RThetaPhiVector>::const_iterator vec = inputVectors.begin(); vec != inputVectors.end(); ++vec ){
    inputVectors_.push_back(math::XYZVector(vec->x(), vec->y(), vec->z()));
  }
}
EventShapeVariables::~EventShapeVariables ( ) [inline]

default destructor

Definition at line 44 of file EventShapeVariables.h.

{};

Member Function Documentation

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().

{
  TVectorD eigenValues = compEigenValues(r);
  return 1.5*eigenValues(2);
}
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().

{
  TVectorD eigenValues = compEigenValues(r);
  return 3.*(eigenValues(0)*eigenValues(1) + eigenValues(0)*eigenValues(2) + eigenValues(1)*eigenValues(2));
}
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 SiPixelRawToDigiRegional_cfi::deltaPhi, i, inputVectors_, j, phi, Pi, tmp, x, and detailsBasic3DVector::y.

Referenced by TtFullHadSignalSel::TtFullHadSignalSel().

{
  const double deltaPhi=2*TMath::Pi()/numberOfSteps;
  double circularity=-1, phi=0, area = 0;
  for(unsigned int i=0;i<inputVectors_.size();i++) {
    area+=TMath::Sqrt(inputVectors_[i].x()*inputVectors_[i].x()+inputVectors_[i].y()*inputVectors_[i].y());
  }
  for(unsigned int i=0; i<numberOfSteps; ++i){
    phi+=deltaPhi;
    double sum=0, tmp=0.;
    for(unsigned int j=0; j<inputVectors_.size(); ++j){
      sum+=TMath::Abs(TMath::Cos(phi)*inputVectors_[j].x()+TMath::Sin(phi)*inputVectors_[j].y());
    }
    tmp=TMath::Pi()/2*sum/area;
    if( circularity<0 || tmp<circularity ){
      circularity=tmp;
    }
  }
  return circularity;
}
TVectorD EventShapeVariables::compEigenValues ( double  r = 2.) const [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(), and sphericity().

{
  TVectorD eigenValues(3);
  TMatrixDSym myTensor = compMomentumTensor(r);
  if( myTensor.IsSymmetric() ){
    if( myTensor.NonZeros() != 0 ) myTensor.EigenVectors(eigenValues);
  }

  // CV: TMatrixDSym::EigenVectors returns eigen-values and eigen-vectors
  //     ordered by descending eigen-values, so no need to do any sorting here...
  //std::cout << "eigenValues(0) = " << eigenValues(0) << ","
  //          << " eigenValues(1) = " << eigenValues(1) << ","
  //          << " eigenValues(2) = " << eigenValues(2) << std::endl;

  return eigenValues;
}
TMatrixDSym EventShapeVariables::compMomentumTensor ( double  r = 2.) const [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 i, inputVectors_, and p2.

Referenced by compEigenValues().

{
  TMatrixDSym momentumTensor(3);
  momentumTensor.Zero();

  if ( inputVectors_.size() < 2 ){
    return momentumTensor;
  }

  // fill momentumTensor from inputVectors
  double norm = 1.;
  for ( int i = 0; i < (int)inputVectors_.size(); ++i ){
    double p2 = inputVectors_[i].Dot(inputVectors_[i]);
    double pR = ( r == 2. ) ? p2 : TMath::Power(p2, 0.5*r);
    norm += pR;
    double pRminus2 = ( r == 2. ) ? 1. : TMath::Power(p2, 0.5*r - 1.);
    momentumTensor(0,0) += pRminus2*inputVectors_[i].x()*inputVectors_[i].x();
    momentumTensor(0,1) += pRminus2*inputVectors_[i].x()*inputVectors_[i].y();
    momentumTensor(0,2) += pRminus2*inputVectors_[i].x()*inputVectors_[i].z();
    momentumTensor(1,0) += pRminus2*inputVectors_[i].y()*inputVectors_[i].x();
    momentumTensor(1,1) += pRminus2*inputVectors_[i].y()*inputVectors_[i].y();
    momentumTensor(1,2) += pRminus2*inputVectors_[i].y()*inputVectors_[i].z();
    momentumTensor(2,0) += pRminus2*inputVectors_[i].z()*inputVectors_[i].x();
    momentumTensor(2,1) += pRminus2*inputVectors_[i].z()*inputVectors_[i].y();
    momentumTensor(2,2) += pRminus2*inputVectors_[i].z()*inputVectors_[i].z();
  }

  //std::cout << "momentumTensor:" << std::endl;
  //std::cout << momentumTensor(0,0) << " " << momentumTensor(0,1) << " " << momentumTensor(0,2) 
  //          << momentumTensor(1,0) << " " << momentumTensor(1,1) << " " << momentumTensor(1,2) 
  //          << momentumTensor(2,0) << " " << momentumTensor(2,1) << " " << momentumTensor(2,2) << std::endl;

  // return momentumTensor normalized to determinant 1
  return (1./norm)*momentumTensor;
}
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().

{
  TVectorD eigenValues = compEigenValues(r);
  return 27.*eigenValues(0)*eigenValues(1)*eigenValues(2);
}
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 SiPixelRawToDigiRegional_cfi::deltaPhi, i, inputVectors_, j, phi, Pi, x, and detailsBasic3DVector::y.

Referenced by TtFullHadSignalSel::TtFullHadSignalSel().

{
  const double deltaPhi=2*TMath::Pi()/numberOfSteps;
  double phi = 0, eIn =-1., eOut=-1.;
  for(unsigned int i=0; i<numberOfSteps; ++i){
    phi+=deltaPhi;
    double sum=0;
    for(unsigned int j=0; j<inputVectors_.size(); ++j){
      // sum over inner product of unit vectors and momenta
      sum+=TMath::Abs(TMath::Cos(phi)*inputVectors_[j].x()+TMath::Sin(phi)*inputVectors_[j].y());
    }
    if( eOut<0. || sum<eOut ) eOut=sum;
    if( eIn <0. || sum>eIn  ) eIn =sum;
  }
  return (eIn-eOut)/eIn;
}
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().

{
  TVectorD eigenValues = compEigenValues(r);
  return 1.5*(eigenValues(1) + eigenValues(2));
}

Member Data Documentation

cashing of input vectors

Definition at line 79 of file EventShapeVariables.h.

Referenced by circularity(), compMomentumTensor(), EventShapeVariables(), and isotropy().