CMS 3D CMS Logo

PFRecHitCaloNavigatorWithTime.h
Go to the documentation of this file.
1 
2 #ifndef RecoParticleFlow_PFClusterProducer_PFRecHitCaloNavigatorWithTime_h
3 #define RecoParticleFlow_PFClusterProducer_PFRecHitCaloNavigatorWithTime_h
4 
8 
14 
19 
22 
24 
25 template <typename D, typename T, bool ownsTopo = true>
27 public:
29  sigmaCut2_ = pow(iConfig.getParameter<double>("sigmaCut"), 2);
30  const edm::ParameterSet& timeResConf = iConfig.getParameterSet("timeResolutionCalc");
31  _timeResolutionCalc.reset(new CaloRecHitResolutionProvider(timeResConf));
32  }
33 
35  if (!ownsTopo) {
36  topology_.release();
37  }
38  }
39 
41  std::unique_ptr<reco::PFRecHitCollection>& hits,
42  edm::RefProd<reco::PFRecHitCollection>& refProd) override {
43  DetId detid(hit.detId());
44 
45  CaloNavigator<D> navigator(detid, topology_.get());
46 
47  DetId N(0);
48  DetId E(0);
49  DetId S(0);
50  DetId W(0);
51  DetId NW(0);
52  DetId NE(0);
53  DetId SW(0);
54  DetId SE(0);
55 
56  N = navigator.north();
57  associateNeighbour(N, hit, hits, refProd, 0, 1);
58 
59  if (N != DetId(0)) {
60  NE = navigator.east();
61  } else {
62  navigator.home();
63  E = navigator.east();
64  NE = navigator.north();
65  }
66  associateNeighbour(NE, hit, hits, refProd, 1, 1);
67  navigator.home();
68 
69  S = navigator.south();
70  associateNeighbour(S, hit, hits, refProd, 0, -1);
71 
72  if (S != DetId(0)) {
73  SW = navigator.west();
74  } else {
75  navigator.home();
76  W = navigator.west();
77  SW = navigator.south();
78  }
79  associateNeighbour(SW, hit, hits, refProd, -1, -1);
80  navigator.home();
81 
82  E = navigator.east();
83  associateNeighbour(E, hit, hits, refProd, 1, 0);
84 
85  if (E != DetId(0)) {
86  SE = navigator.south();
87  } else {
88  navigator.home();
89  S = navigator.south();
90  SE = navigator.east();
91  }
92  associateNeighbour(SE, hit, hits, refProd, 1, -1);
93  navigator.home();
94 
95  W = navigator.west();
96  associateNeighbour(W, hit, hits, refProd, -1, 0);
97 
98  if (W != DetId(0)) {
99  NW = navigator.north();
100  } else {
101  navigator.home();
102  N = navigator.north();
103  NW = navigator.west();
104  }
105  associateNeighbour(NW, hit, hits, refProd, -1, 1);
106  }
107 
108 protected:
109  double sigmaCut2_;
110  std::unique_ptr<const T> topology_;
111  std::unique_ptr<CaloRecHitResolutionProvider> _timeResolutionCalc;
112 
113  void associateNeighbour(const DetId& id,
115  std::unique_ptr<reco::PFRecHitCollection>& hits,
117  short eta,
118  short phi) {
119  double sigma2 = 10000.0;
120 
121  auto found_hit = std::lower_bound(hits->begin(), hits->end(), id, [](const reco::PFRecHit& a, DetId b) {
122  if (DetId(a.detId()).det() == DetId::Hcal || b.det() == DetId::Hcal)
123  return (HcalDetId)(a.detId()) < (HcalDetId)(b);
124 
125  else
126  return a.detId() < b.rawId();
127  });
128 
129  if (found_hit != hits->end() &&
130  ((found_hit->detId() == id.rawId()) ||
131  ((id.det() == DetId::Hcal) && ((HcalDetId)(found_hit->detId()) == (HcalDetId)(id))))) {
132  sigma2 = _timeResolutionCalc->timeResolution2(hit.energy()) +
133  _timeResolutionCalc->timeResolution2(found_hit->energy());
134  const double deltaTime = hit.time() - found_hit->time();
135  if (deltaTime * deltaTime < sigmaCut2_ * sigma2) {
136  hit.addNeighbour(eta, phi, 0, found_hit - hits->begin());
137  }
138  }
139  }
140 };
141 
142 #endif
edm::RefProd
Definition: EDProductfwd.h:25
EcalPreshowerTopology.h
CaloRecHitResolutionProvider
Definition: CaloRecHitResolutionProvider.h:8
CaloNavigator.h
hfClusterShapes_cfi.hits
hits
Definition: hfClusterShapes_cfi.py:5
EBDetId.h
EEDetId.h
DetId::Hcal
Definition: DetId.h:28
EcalBarrelTopology.h
PFRecHitCaloNavigatorWithTime::~PFRecHitCaloNavigatorWithTime
~PFRecHitCaloNavigatorWithTime() override
Definition: PFRecHitCaloNavigatorWithTime.h:34
ESDetId.h
PFRecHitNavigatorBase.h
PFRecHitCaloNavigatorWithTime::associateNeighbour
void associateNeighbour(const DetId &id, reco::PFRecHit &hit, std::unique_ptr< reco::PFRecHitCollection > &hits, edm::RefProd< reco::PFRecHitCollection > &refProd, short eta, short phi)
Definition: PFRecHitCaloNavigatorWithTime.h:113
DetId
Definition: DetId.h:17
PFRecHitCaloNavigatorWithTime::_timeResolutionCalc
std::unique_ptr< CaloRecHitResolutionProvider > _timeResolutionCalc
Definition: PFRecHitCaloNavigatorWithTime.h:111
PVValHelper::eta
Definition: PVValidationHelpers.h:69
HLT_2018_cff.navigator
navigator
Definition: HLT_2018_cff.py:11734
PFRecHitCaloNavigatorWithTime::sigmaCut2_
double sigmaCut2_
Definition: PFRecHitCaloNavigatorWithTime.h:109
N
#define N
Definition: blowfish.cc:9
CaloTowerTopology.h
b
double b
Definition: hdecay.h:118
cuda_std::lower_bound
__host__ constexpr __device__ RandomIt lower_bound(RandomIt first, RandomIt last, const T &value, Compare comp={})
Definition: cudastdAlgorithm.h:27
S
double S(const TLorentzVector &, const TLorentzVector &)
Definition: Particle.cc:97
CaloSubdetectorGeometry.h
edm::ParameterSet
Definition: ParameterSet.h:36
a
double a
Definition: hdecay.h:119
EcalEndcapTopology.h
HcalDetId.h
PFRecHitNavigatorBase
Definition: PFRecHitNavigatorBase.h:26
HcalDetId
Definition: HcalDetId.h:12
CaloTowerDetId.h
DDAxes::phi
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
HcalTopology.h
CaloNavigator
Definition: CaloNavigator.h:7
triggerObjects_cff.id
id
Definition: triggerObjects_cff.py:31
CaloRecHitResolutionProvider.h
PFRecHitCaloNavigatorWithTime::PFRecHitCaloNavigatorWithTime
PFRecHitCaloNavigatorWithTime(const edm::ParameterSet &iConfig)
Definition: PFRecHitCaloNavigatorWithTime.h:28
CaloGeometry.h
PFRecHitCaloNavigatorWithTime::topology_
std::unique_ptr< const T > topology_
Definition: PFRecHitCaloNavigatorWithTime.h:110
S
Definition: CSCDBL1TPParametersExtended.h:16
reco::PFRecHit
Particle flow rechit (rechit + geometry and topology information). See clustering algorithm in PFClus...
Definition: PFRecHit.h:31
funct::pow
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:30
PFRecHitCaloNavigatorWithTime
Definition: PFRecHitCaloNavigatorWithTime.h:26
PFRecHitCaloNavigatorWithTime::associateNeighbours
void associateNeighbours(reco::PFRecHit &hit, std::unique_ptr< reco::PFRecHitCollection > &hits, edm::RefProd< reco::PFRecHitCollection > &refProd) override
Definition: PFRecHitCaloNavigatorWithTime.h:40
edm::ParameterSet::getParameterSet
ParameterSet const & getParameterSet(std::string const &) const
Definition: ParameterSet.cc:2121
hit
Definition: SiStripHitEffFromCalibTree.cc:88