ATLAS Offline Software
Loading...
Searching...
No Matches
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

Detailed Description

Definition at line 11 of file SphericityTensor.h.

Member Function Documentation

◆ result() [1/2]

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}
std::vector< fastjet::PseudoJet > boostToCenterOfMass(const fastjet::PseudoJet &jet, std::vector< fastjet::PseudoJet > constituents)

◆ result() [2/2]

virtual TOut JetSubStructureUtils::SubstructureCalculator< TOut >::result ( const xAOD::Jet & jet) const
inlinevirtual

Reimplemented from JetSubStructureUtils::SubstructureCalculator< std::map< std::string, double > >.

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 }
virtual std::map< std::string, double > result(const fastjet::PseudoJet &jet) const
iterator begin() const
iterator on the first constituent
iterator end() const
iterator after the last constituent
JetConstituentVector getConstituents() const
Return a vector of consituents. The object behaves like vector<const IParticle*>. See JetConstituentV...
Definition Jet_v1.cxx:147

The documentation for this class was generated from the following files: