ATLAS Offline Software
BoostToCenterOfMass.cxx
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2020 CERN for the benefit of the ATLAS collaboration
3 */
4 
6 #include "TLorentzVector.h"
7 using namespace std;
8 
9 vector<fastjet::PseudoJet> JetSubStructureUtils::boostToCenterOfMass(const fastjet::PseudoJet& jet,
10  vector<fastjet::PseudoJet> constit_pseudojets) {
11  vector<fastjet::PseudoJet> clusters;
12  if(jet.e() < 1e-20) { // FPE
13  return clusters;
14  }
15 
16  double bx = jet.px()/jet.e();
17  double by = jet.py()/jet.e();
18  double bz = jet.pz()/jet.e();
19 
20  if(bx*bx + by*by + bz*bz >= 1) { // Faster than light
21  return clusters;
22  }
23 
24  for(unsigned int i1=0; i1 < constit_pseudojets.size(); i1++) {
25  TLorentzVector v;
26  v.SetPxPyPzE(constit_pseudojets.at(i1).px(), constit_pseudojets.at(i1).py(),constit_pseudojets.at(i1).pz(),constit_pseudojets.at(i1).e());
27  v.Boost(-bx,-by,-bz);
28  fastjet::PseudoJet v2(v.Px(), v.Py(), v.Pz(), v.E());
29  clusters.push_back(v2);
30  }
31 
32  return clusters;
33 }
JetSubStructureUtils::boostToCenterOfMass
std::vector< fastjet::PseudoJet > boostToCenterOfMass(const fastjet::PseudoJet &jet, std::vector< fastjet::PseudoJet > constituents)
Definition: BoostToCenterOfMass.cxx:9
jet
Definition: JetCalibTools_PlotJESFactors.cxx:23
fitman.bz
bz
Definition: fitman.py:412
fitman.bx
bx
Definition: fitman.py:410
fitman.by
by
Definition: fitman.py:411
BoostToCenterOfMass.h
ReadCellNoiseFromCoolCompare.v2
v2
Definition: ReadCellNoiseFromCoolCompare.py:364
python.PyAthena.v
v
Definition: PyAthena.py:157
DiTauMassTools::MaxHistStrategyV2::e
e
Definition: PhysicsAnalysis/TauID/DiTauMassTools/DiTauMassTools/HelperFunctions.h:26
RunTileMonitoring.clusters
clusters
Definition: RunTileMonitoring.py:133