ATLAS Offline Software
Public Member Functions | List of all members
JetSubStructureUtils::SphericityTensor Class Reference

#include <SphericityTensor.h>

Inheritance diagram for JetSubStructureUtils::SphericityTensor:
Collaboration diagram for JetSubStructureUtils::SphericityTensor:

Public Member Functions

virtual std::map< std::string, double > result (const fastjet::PseudoJet &jet) const
 
virtual TOut result (const xAOD::Jet &jet) const
 
virtual std::map< std::string, double > result (const xAOD::Jet &jet) const
 

Detailed Description

Definition at line 11 of file SphericityTensor.h.

Member Function Documentation

◆ result() [1/3]

map< string, double > SphericityTensor::result ( const fastjet::PseudoJet &  jet) const
virtual

Definition at line 15 of file SphericityTensor.cxx.

16 {
17  map<string, double> Variables;
18  Variables["Sphericity"] = -999.*1000.;
19  Variables["Aplanarity"] = -999.*1000.;
20 
21  vector<fastjet::PseudoJet> clusters = boostToCenterOfMass(jet, jet.constituents());
22  if(clusters.size() < 2) return Variables;
23 
24  TMatrixD MomentumTensor(3,3);
25  double P2Sum = 0;
26 
27  for(std::vector<fastjet::PseudoJet>::const_iterator Itr=clusters.begin(); Itr!=clusters.end(); ++Itr) {
28  MomentumTensor(0,0) += (*Itr).px()*(*Itr).px();
29  MomentumTensor(0,1) += (*Itr).px()*(*Itr).py();
30  MomentumTensor(0,2) += (*Itr).px()*(*Itr).pz();
31  MomentumTensor(1,0) += (*Itr).py()*(*Itr).px();
32  MomentumTensor(1,1) += (*Itr).py()*(*Itr).py();
33  MomentumTensor(1,2) += (*Itr).py()*(*Itr).pz();
34  MomentumTensor(2,0) += (*Itr).pz()*(*Itr).px();
35  MomentumTensor(2,1) += (*Itr).pz()*(*Itr).py();
36  MomentumTensor(2,2) += (*Itr).pz()*(*Itr).pz();
37 
38  P2Sum += (*Itr).px()*(*Itr).px()+(*Itr).py()*(*Itr).py()+(*Itr).pz()*(*Itr).pz();
39  }
40 
41  double Aplanarity = -1;
42  double Sphericity = -1;
43 
44  if(P2Sum > 0) {
45  const double inv_P2Sum = 1. / P2Sum;
46  for(int i=0; i<3; i++) {
47  for(int j=0; j<3; j++) {
48  MomentumTensor(i,j) *= inv_P2Sum;
49  }
50  }
51 
52  TDecompSVD * aSVD = new TDecompSVD(MomentumTensor);
53  TVectorD Lambda = aSVD->GetSig();
54 
55  Aplanarity = 1.5*Lambda[2];
56  Sphericity = 1.5*(Lambda[1]+Lambda[2]);
57 
58  delete aSVD;
59  }
60  else {
61  return Variables;
62  }
63 
64  Variables["Aplanarity"] = Aplanarity;
65  Variables["Sphericity"] = Sphericity;
66  return Variables;
67 }

◆ result() [2/3]

virtual TOut JetSubStructureUtils::SubstructureCalculator< TOut >::result
inline

Definition at line 25 of file SubstructureCalculator.h.

25  {
26  // PS 4/18 master developent
27  // std::vector<fastjet::PseudoJet> constit_pseudojets =
28  // jet::JetConstituentFiller::constituentPseudoJets(jet);
29 
30  std::vector<fastjet::PseudoJet> constit_pseudojets;
31  std::transform(jet.getConstituents().begin(),
32  jet.getConstituents().end(),
33  std::back_inserter(constit_pseudojets),
34  [](const auto& c){
35  const xAOD::IParticle* ip = c->rawConstituent();
36  return
37  // fastjet::PseudoJet((c->rawConstituent())->p4());
38  fastjet::PseudoJet(ip->p4());
39  });
40 
41  fastjet::PseudoJet pjet = fastjet::join(constit_pseudojets);
42 
43  return result(pjet);
44  }

◆ result() [3/3]

virtual std::map< std::string, double > JetSubStructureUtils::SubstructureCalculator< std::map< std::string, double > >::result ( const xAOD::Jet jet) const
inlinevirtualinherited

Definition at line 25 of file SubstructureCalculator.h.

25  {
26  // PS 4/18 master developent
27  // std::vector<fastjet::PseudoJet> constit_pseudojets =
28  // jet::JetConstituentFiller::constituentPseudoJets(jet);
29 
30  std::vector<fastjet::PseudoJet> constit_pseudojets;
31  std::transform(jet.getConstituents().begin(),
32  jet.getConstituents().end(),
33  std::back_inserter(constit_pseudojets),
34  [](const auto& c){
35  const xAOD::IParticle* ip = c->rawConstituent();
36  return
37  // fastjet::PseudoJet((c->rawConstituent())->p4());
38  fastjet::PseudoJet(ip->p4());
39  });
40 
41  fastjet::PseudoJet pjet = fastjet::join(constit_pseudojets);
42 
43  return result(pjet);
44  }

The documentation for this class was generated from the following files:
JetSubStructureUtils::boostToCenterOfMass
std::vector< fastjet::PseudoJet > boostToCenterOfMass(const fastjet::PseudoJet &jet, std::vector< fastjet::PseudoJet > constituents)
Definition: BoostToCenterOfMass.cxx:9
JetSubStructureUtils::SubstructureCalculator< std::map< std::string, double > >::result
virtual std::map< std::string, double > result(const xAOD::Jet &jet) const
Definition: SubstructureCalculator.h:25
JetSubStructureUtils::SphericityTensor::result
virtual std::map< std::string, double > result(const fastjet::PseudoJet &jet) const
Definition: SphericityTensor.cxx:15
jet
Definition: JetCalibTools_PlotJESFactors.cxx:23
lumiFormat.i
int i
Definition: lumiFormat.py:92
Amg::transform
Amg::Vector3D transform(Amg::Vector3D &v, Amg::Transform3D &tr)
Transform a point from a Trasformation3D.
Definition: GeoPrimitivesHelpers.h:156
TCS::join
std::string join(const std::vector< std::string > &v, const char c=',')
Definition: Trigger/TrigT1/L1Topo/L1TopoCommon/Root/StringUtils.cxx:10
GlobalVariables.Variables
Variables
Definition: GlobalVariables.py:276
RunTileMonitoring.clusters
clusters
Definition: RunTileMonitoring.py:133
FlavorTagDiscriminants::TruthDecoratorHelpers::TruthType::Lambda
@ Lambda
Definition: TruthDecoratorHelpers.h:21
python.compressB64.c
def c
Definition: compressB64.py:93