Main Page
Namespaces
Classes
Package Documentation
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Properties
Friends
Macros
Groups
Pages
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
12
#include <list>
13
#include <algorithm>
14
15
#include "
DataFormats/Candidate/interface/Candidate.h
"
16
#include "
DataFormats/Math/interface/deltaR.h
"
17
#include <limits>
18
#include <vector>
19
#include <list>
20
21
#include "fastjet/PseudoJet.hh"
22
23
class
CMSInsideOutAlgorithm
{
24
public
:
25
typedef
reco::Particle::LorentzVector
LorentzVector
;
26
typedef
std::list<fastjet::PseudoJet>::iterator
inputListIter
;
27
// binary predicate to sort a std::list of std::list<InputItem> iterators by increasing deltaR
28
// from a eta-phi point specified in the ctor
29
class
ListIteratorLesserByDeltaR
{
30
public
:
31
ListIteratorLesserByDeltaR
(
const
double
&
eta
,
const
double
&
phi
) :
seedEta_
(eta),
seedPhi_
(phi) {}
32
bool
operator()
(
const
inputListIter
&
A
,
const
inputListIter
&
B
)
const
{
33
double
deltaR2A =
reco::deltaR2
((*A).eta(),
seedEta_
, (*A).phi(),
seedPhi_
);
34
double
deltaR2B =
reco::deltaR2
((*B).eta(),
seedEta_
, (*B).phi(),
seedPhi_
);
35
return
fabs(deltaR2A - deltaR2B) >
std::numeric_limits<double>::epsilon
()
36
? deltaR2A < deltaR2B
37
:
reco::deltaPhi
((*A).phi(),
seedPhi_
) <
reco::deltaPhi
((*B).phi(),
seedPhi_
);
38
};
39
40
private
:
41
double
seedEta_
,
seedPhi_
;
42
};
43
50
CMSInsideOutAlgorithm
(
double
seedObjectPt,
double
growthParameter,
double
maxSize
,
double
minSize)
51
:
seedThresholdPt_
(seedObjectPt),
52
growthParameterSquared_
(growthParameter * growthParameter),
53
maxSizeSquared_
(maxSize * maxSize),
54
minSizeSquared_
(minSize * minSize){};
55
57
void
run
(
const
std::vector<fastjet::PseudoJet>& fInput, std::vector<fastjet::PseudoJet>& fOutput);
58
59
private
:
60
double
seedThresholdPt_
;
61
double
growthParameterSquared_
;
62
double
maxSizeSquared_
;
63
double
minSizeSquared_
;
64
};
65
66
#endif
reco::deltaPhi
constexpr double deltaPhi(double phi1, double phi2)
Definition:
deltaPhi.h:26
CMSInsideOutAlgorithm::ListIteratorLesserByDeltaR::seedEta_
double seedEta_
Definition:
CMSInsideOutAlgorithm.h:38
deltaR.h
CMSInsideOutAlgorithm
Definition:
CMSInsideOutAlgorithm.h:23
CMSInsideOutAlgorithm::ListIteratorLesserByDeltaR::seedPhi_
double seedPhi_
Definition:
CMSInsideOutAlgorithm.h:38
CMSInsideOutAlgorithm::growthParameterSquared_
double growthParameterSquared_
Definition:
CMSInsideOutAlgorithm.h:61
PVValHelper::eta
Definition:
PVValidationHelpers.h:70
DDAxes::phi
CMSInsideOutAlgorithm::inputListIter
std::list< fastjet::PseudoJet >::iterator inputListIter
Definition:
CMSInsideOutAlgorithm.h:26
Candidate.h
CMSInsideOutAlgorithm::LorentzVector
reco::Particle::LorentzVector LorentzVector
Definition:
CMSInsideOutAlgorithm.h:25
CMSInsideOutAlgorithm::seedThresholdPt_
double seedThresholdPt_
Definition:
CMSInsideOutAlgorithm.h:60
reco_skim_cfg_mod.maxSize
tuple maxSize
'/store/data/Commissioning08/BeamHalo/RECO/StuffAlmostToP5_v1/000/061/642/10A0FE34-A67D-DD11-AD05-000...
Definition:
reco_skim_cfg_mod.py:154
CMSInsideOutAlgorithm::minSizeSquared_
double minSizeSquared_
Definition:
CMSInsideOutAlgorithm.h:63
CMSInsideOutAlgorithm::CMSInsideOutAlgorithm
CMSInsideOutAlgorithm(double seedObjectPt, double growthParameter, double maxSize, double minSize)
Definition:
CMSInsideOutAlgorithm.h:50
TtFullHadDaughter::B
static const std::string B
Definition:
TtFullHadronicEvent.h:9
CMSInsideOutAlgorithm::ListIteratorLesserByDeltaR
Definition:
CMSInsideOutAlgorithm.h:29
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:31
CMSInsideOutAlgorithm::ListIteratorLesserByDeltaR::operator()
bool operator()(const inputListIter &A, const inputListIter &B) const
Definition:
CMSInsideOutAlgorithm.h:32
CMSInsideOutAlgorithm::maxSizeSquared_
double maxSizeSquared_
Definition:
CMSInsideOutAlgorithm.h:62
CMSInsideOutAlgorithm::run
void run(const std::vector< fastjet::PseudoJet > &fInput, std::vector< fastjet::PseudoJet > &fOutput)
Build from input candidate collection.
Definition:
CMSInsideOutAlgorithm.cc:7
funct::A
A
Definition:
Factorize.h:45
geometryDiff.epsilon
int epsilon
Definition:
geometryDiff.py:26
reco::Particle::LorentzVector
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition:
Particle.h:21
Generated for CMSSW Reference Manual by
1.8.5