Main Page
Namespaces
Classes
Package Documentation
RecoJets
JetAlgorithms
interface
CMSInsideOutAlgorithm.h
Go to the documentation of this file.
1
#ifndef RecoJets_JetAlgorithms_CMSInsideOutAlgorithm_h
2
#define RecoJets_JetAlgorithms_CMSInsideOutAlgorithm_h
3
13
#include <list>
14
#include <algorithm>
15
16
#include "
DataFormats/Candidate/interface/Candidate.h
"
17
#include "
DataFormats/Math/interface/deltaR.h
"
18
#include <limits>
19
#include <vector>
20
#include <list>
21
22
#include "fastjet/PseudoJet.hh"
23
24
25
26
27
class
CMSInsideOutAlgorithm
{
28
public
:
29
typedef
reco::Particle::LorentzVector
LorentzVector
;
30
typedef
std::list<fastjet::PseudoJet>::iterator
inputListIter
;
31
// binary predicate to sort a std::list of std::list<InputItem> iterators by increasing deltaR
32
// from a eta-phi point specified in the ctor
33
class
ListIteratorLesserByDeltaR
{
34
public
:
ListIteratorLesserByDeltaR
(
const
double
&
eta
,
const
double
&
phi
):
seedEta_
(eta),
seedPhi_
(phi){}
35
bool
operator()
(
const
inputListIter&
A
,
const
inputListIter&
B
)
const
{
36
double
deltaR2A =
reco::deltaR2
( (*A).eta(),
seedEta_
, (*A).phi(),
seedPhi_
);
37
double
deltaR2B =
reco::deltaR2
( (*B).eta(),
seedEta_
, (*B).phi(),
seedPhi_
);
38
return
39
fabs(deltaR2A - deltaR2B) >
std::numeric_limits<double>::epsilon
() ? deltaR2A < deltaR2B :
40
reco::deltaPhi
((*A).phi(),
seedPhi_
) <
reco::deltaPhi
((*B).phi(),
seedPhi_
);
41
};
42
private
:
43
double
seedEta_
,
seedPhi_
;
44
};
45
52
CMSInsideOutAlgorithm
(
double
seedObjectPt,
double
growthParameter,
double
maxSize
,
double
minSize):
53
seedThresholdPt_
(seedObjectPt),
54
growthParameterSquared_
(growthParameter*growthParameter),
55
maxSizeSquared_
(maxSize*maxSize),
56
minSizeSquared_
(minSize*minSize){};
57
58
59
61
void
run
(
const
std::vector<fastjet::PseudoJet>& fInput, std::vector<fastjet::PseudoJet> & fOutput);
62
63
private
:
64
double
seedThresholdPt_
;
65
double
growthParameterSquared_
;
66
double
maxSizeSquared_
;
67
double
minSizeSquared_
;
68
};
69
70
#endif
reco::deltaPhi
constexpr double deltaPhi(double phi1, double phi2)
Definition:
deltaPhi.h:22
CMSInsideOutAlgorithm::ListIteratorLesserByDeltaR::seedEta_
double seedEta_
Definition:
CMSInsideOutAlgorithm.h:41
CMSInsideOutAlgorithm
Definition:
CMSInsideOutAlgorithm.h:27
CMSInsideOutAlgorithm::ListIteratorLesserByDeltaR::seedPhi_
double seedPhi_
Definition:
CMSInsideOutAlgorithm.h:41
copyPickMerge_cfg.maxSize
maxSize
Definition:
copyPickMerge_cfg.py:39
CMSInsideOutAlgorithm::growthParameterSquared_
double growthParameterSquared_
Definition:
CMSInsideOutAlgorithm.h:65
PVValHelper::eta
Definition:
PVValidationHelpers.h:65
DDAxes::phi
geometryDiff.epsilon
int epsilon
Definition:
geometryDiff.py:25
CMSInsideOutAlgorithm::inputListIter
std::list< fastjet::PseudoJet >::iterator inputListIter
Definition:
CMSInsideOutAlgorithm.h:30
Candidate.h
CMSInsideOutAlgorithm::LorentzVector
reco::Particle::LorentzVector LorentzVector
Definition:
CMSInsideOutAlgorithm.h:29
CMSInsideOutAlgorithm::seedThresholdPt_
double seedThresholdPt_
Definition:
CMSInsideOutAlgorithm.h:64
deltaR.h
CMSInsideOutAlgorithm::minSizeSquared_
double minSizeSquared_
Definition:
CMSInsideOutAlgorithm.h:67
CMSInsideOutAlgorithm::CMSInsideOutAlgorithm
CMSInsideOutAlgorithm(double seedObjectPt, double growthParameter, double maxSize, double minSize)
Definition:
CMSInsideOutAlgorithm.h:52
TtFullHadDaughter::B
static const std::string B
Definition:
TtFullHadronicEvent.h:9
CMSInsideOutAlgorithm::ListIteratorLesserByDeltaR
Definition:
CMSInsideOutAlgorithm.h:33
reco::deltaR2
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition:
deltaR.h:16
CMSInsideOutAlgorithm::ListIteratorLesserByDeltaR::ListIteratorLesserByDeltaR
ListIteratorLesserByDeltaR(const double &eta, const double &phi)
Definition:
CMSInsideOutAlgorithm.h:34
CMSInsideOutAlgorithm::ListIteratorLesserByDeltaR::operator()
bool operator()(const inputListIter &A, const inputListIter &B) const
Definition:
CMSInsideOutAlgorithm.h:35
CMSInsideOutAlgorithm::maxSizeSquared_
double maxSizeSquared_
Definition:
CMSInsideOutAlgorithm.h:66
patCaloMETCorrections_cff.A
A
Definition:
patCaloMETCorrections_cff.py:43
CMSInsideOutAlgorithm::run
void run(const std::vector< fastjet::PseudoJet > &fInput, std::vector< fastjet::PseudoJet > &fOutput)
Build from input candidate collection.
Definition:
CMSInsideOutAlgorithm.cc:9
reco::Particle::LorentzVector
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition:
Particle.h:21
Generated for CMSSW Reference Manual by
1.8.11