RecoParticleFlow
PFProducer
plugins
linkers
SCAndHGCalLinker.cc
Go to the documentation of this file.
1
#include "
RecoParticleFlow/PFProducer/interface/BlockElementLinkerBase.h
"
2
#include "
DataFormats/ParticleFlowReco/interface/PFCluster.h
"
3
#include "
DataFormats/ParticleFlowReco/interface/PFBlockElementCluster.h
"
4
#include "
DataFormats/ParticleFlowReco/interface/PFBlockElementSuperCluster.h
"
5
#include "
RecoParticleFlow/PFClusterTools/interface/LinkByRecHit.h
"
6
#include "
RecoParticleFlow/PFClusterTools/interface/ClusterClusterMapping.h
"
7
8
class
SCAndHGCalLinker
:
public
BlockElementLinkerBase
{
9
public
:
10
SCAndHGCalLinker
(
const
edm::ParameterSet
& conf)
11
:
BlockElementLinkerBase
(conf),
12
useKDTree_
(conf.getParameter<
bool
>(
"useKDTree"
)),
13
debug_
(conf.getUntrackedParameter<
bool
>(
"debug"
,
false
)),
14
superClusterMatchByRef_
(conf.getParameter<
bool
>(
"SuperClusterMatchByRef"
)) {}
15
16
double
testLink
(
const
reco::PFBlockElement
*,
const
reco::PFBlockElement
*)
const override
;
17
18
private
:
19
bool
useKDTree_
,
debug_
,
superClusterMatchByRef_
;
20
};
21
22
DEFINE_EDM_PLUGIN
(
BlockElementLinkerFactory
,
SCAndHGCalLinker
,
"SCAndHGCalLinker"
);
23
24
double
SCAndHGCalLinker::testLink
(
const
reco::PFBlockElement
* elem1,
const
reco::PFBlockElement
* elem2)
const
{
25
double
dist = -1.0;
26
const
reco::PFBlockElementCluster
* ecalelem(
nullptr
);
27
const
reco::PFBlockElementSuperCluster
* scelem(
nullptr
);
28
if
(elem1->
type
() > elem2->
type
()) {
29
ecalelem = static_cast<const reco::PFBlockElementCluster*>(elem1);
30
scelem = static_cast<const reco::PFBlockElementSuperCluster*>(elem2);
31
}
else
{
32
ecalelem = static_cast<const reco::PFBlockElementCluster*>(elem2);
33
scelem = static_cast<const reco::PFBlockElementSuperCluster*>(elem1);
34
}
35
const
reco::PFClusterRef
& clus = ecalelem->
clusterRef
();
36
const
reco::SuperClusterRef
& sclus = scelem->
superClusterRef
();
37
if
(sclus.
isNull
()) {
38
throw
cms::Exception
(
"BadRef"
) <<
"SuperClusterRef is invalid!"
;
39
}
40
41
if
(
superClusterMatchByRef_
) {
42
if
(sclus == ecalelem->
superClusterRef
())
43
dist = 0.001;
44
}
else
{
45
if
(
ClusterClusterMapping::overlap
(*sclus, *clus)) {
46
dist =
LinkByRecHit::computeDist
(
47
sclus->position().eta(), sclus->position().phi(), clus->positionREP().Eta(), clus->positionREP().Phi());
48
}
49
}
50
51
return
dist;
52
}
SCAndHGCalLinker::useKDTree_
bool useKDTree_
Definition:
SCAndHGCalLinker.cc:19
electrons_cff.bool
bool
Definition:
electrons_cff.py:393
funct::false
false
Definition:
Factorize.h:29
PFBlockElementCluster.h
edm::Ref::isNull
bool isNull() const
Checks for null.
Definition:
Ref.h:235
reco::PFBlockElementSuperCluster
Cluster Element.
Definition:
PFBlockElementSuperCluster.h:15
BlockElementLinkerBase
Definition:
BlockElementLinkerBase.h:10
edm::Ref< PFClusterCollection >
SCAndHGCalLinker
Definition:
SCAndHGCalLinker.cc:8
PFCluster.h
reco::PFBlockElementCluster::superClusterRef
const SuperClusterRef & superClusterRef() const
Definition:
PFBlockElementCluster.h:30
DEFINE_EDM_PLUGIN
#define DEFINE_EDM_PLUGIN(factory, type, name)
Definition:
PluginFactory.h:124
ClusterClusterMapping::overlap
static bool overlap(const reco::CaloCluster &sc1, const reco::CaloCluster &sc, float minfrac=0.01, bool debug=false)
Definition:
ClusterClusterMapping.cc:4
edm::ParameterSet
Definition:
ParameterSet.h:47
edmplugin::PluginFactory
Definition:
PluginFactory.h:34
SCAndHGCalLinker::SCAndHGCalLinker
SCAndHGCalLinker(const edm::ParameterSet &conf)
Definition:
SCAndHGCalLinker.cc:10
SCAndHGCalLinker::testLink
double testLink(const reco::PFBlockElement *, const reco::PFBlockElement *) const override
Definition:
SCAndHGCalLinker.cc:24
LinkByRecHit.h
PFBlockElementSuperCluster.h
reco::PFBlockElement
Abstract base class for a PFBlock element (track, cluster...)
Definition:
PFBlockElement.h:26
reco::PFBlockElementSuperCluster::superClusterRef
const SuperClusterRef & superClusterRef() const
Definition:
PFBlockElementSuperCluster.h:36
ClusterClusterMapping.h
SCAndHGCalLinker::debug_
bool debug_
Definition:
SCAndHGCalLinker.cc:19
reco::PFBlockElementCluster
Cluster Element.
Definition:
PFBlockElementCluster.h:16
Exception
Definition:
hltDiff.cc:246
BlockElementLinkerBase.h
reco::PFBlockElementCluster::clusterRef
const PFClusterRef & clusterRef() const override
Definition:
PFBlockElementCluster.h:29
reco::PFBlockElement::type
Type type() const
Definition:
PFBlockElement.h:69
SCAndHGCalLinker::superClusterMatchByRef_
bool superClusterMatchByRef_
Definition:
SCAndHGCalLinker.cc:19
LinkByRecHit::computeDist
static double computeDist(double eta1, double phi1, double eta2, double phi2, bool etaPhi=true)
computes a chisquare
Definition:
LinkByRecHit.cc:519
Generated for CMSSW Reference Manual by
1.8.16