-
Notifications
You must be signed in to change notification settings - Fork 4.3k
/
ECALAndBREMLinker.cc
41 lines (36 loc) · 1.83 KB
/
ECALAndBREMLinker.cc
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
#include "RecoParticleFlow/PFProducer/interface/BlockElementLinkerBase.h"
#include "DataFormats/ParticleFlowReco/interface/PFCluster.h"
#include "DataFormats/ParticleFlowReco/interface/PFBlockElementCluster.h"
#include "DataFormats/ParticleFlowReco/interface/PFBlockElementBrem.h"
#include "RecoParticleFlow/PFClusterTools/interface/LinkByRecHit.h"
class ECALAndBREMLinker : public BlockElementLinkerBase {
public:
ECALAndBREMLinker(const edm::ParameterSet& conf)
: BlockElementLinkerBase(conf),
useKDTree_(conf.getParameter<bool>("useKDTree")),
debug_(conf.getUntrackedParameter<bool>("debug", false)) {}
double testLink(const reco::PFBlockElement*, const reco::PFBlockElement*) const override;
private:
bool useKDTree_, debug_;
};
DEFINE_EDM_PLUGIN(BlockElementLinkerFactory, ECALAndBREMLinker, "ECALAndBREMLinker");
double ECALAndBREMLinker::testLink(const reco::PFBlockElement* elem1, const reco::PFBlockElement* elem2) const {
constexpr reco::PFTrajectoryPoint::LayerType ECALShowerMax = reco::PFTrajectoryPoint::ECALShowerMax;
const reco::PFBlockElementCluster* ecalelem(nullptr);
const reco::PFBlockElementBrem* bremelem(nullptr);
double dist(-1.0);
if (elem1->type() < elem2->type()) {
ecalelem = static_cast<const reco::PFBlockElementCluster*>(elem1);
bremelem = static_cast<const reco::PFBlockElementBrem*>(elem2);
} else {
ecalelem = static_cast<const reco::PFBlockElementCluster*>(elem2);
bremelem = static_cast<const reco::PFBlockElementBrem*>(elem1);
}
const reco::PFClusterRef& clusterref = ecalelem->clusterRef();
const reco::PFRecTrack& track = bremelem->trackPF();
const reco::PFTrajectoryPoint& tkAtECAL = track.extrapolatedPoint(ECALShowerMax);
if (tkAtECAL.isValid()) {
dist = LinkByRecHit::testTrackAndClusterByRecHit(track, *clusterref, true, debug_);
}
return dist;
}