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 
5 
9 
15 
20 
23 
25 
26 template <typename D,typename T,bool ownsTopo=true>
28  public:
30  sigmaCut2_ = pow(iConfig.getParameter<double>("sigmaCut"), 2);
31  const edm::ParameterSet& timeResConf = iConfig.getParameterSet("timeResolutionCalc");
32  _timeResolutionCalc.reset(new CaloRecHitResolutionProvider(timeResConf));
33  }
34 
35 
36  virtual ~PFRecHitCaloNavigatorWithTime() { if(!ownsTopo) { topology_.release(); } }
37 
38 
39  void associateNeighbours(reco::PFRecHit& hit,std::unique_ptr<reco::PFRecHitCollection>& hits,edm::RefProd<reco::PFRecHitCollection>& refProd) {
40  DetId detid( hit.detId() );
41 
42  CaloNavigator<D> navigator(detid, topology_.get());
43 
44  DetId N(0);
45  DetId E(0);
46  DetId S(0);
47  DetId W(0);
48  DetId NW(0);
49  DetId NE(0);
50  DetId SW(0);
51  DetId SE(0);
52 
53 
54  N=navigator.north();
55  associateNeighbour(N,hit,hits,refProd,0,1);
56 
57 
58  if (N !=DetId(0)) {
59  NE=navigator.east();
60  }
61  else
62  {
63  navigator.home();
64  E=navigator.east();
65  NE=navigator.north();
66  }
67  associateNeighbour(NE,hit,hits,refProd,1,1);
68  navigator.home();
69 
70  S = navigator.south();
71  associateNeighbour(S,hit,hits,refProd,0,-1);
72 
73  if (S !=DetId(0)) {
74  SW = navigator.west();
75  } else {
76  navigator.home();
77  W=navigator.west();
78  SW=navigator.south();
79  }
80  associateNeighbour(SW,hit,hits,refProd,-1,-1);
81  navigator.home();
82 
83  E = navigator.east();
84  associateNeighbour(E,hit,hits,refProd,1,0);
85 
86  if (E !=DetId(0)) {
87  SE = navigator.south();
88  } else {
89  navigator.home();
90  S=navigator.south();
91  SE=navigator.east();
92  }
93  associateNeighbour(SE,hit,hits,refProd,1,-1);
94  navigator.home();
95 
96 
97  W = navigator.west();
98  associateNeighbour(W,hit,hits,refProd,-1,0);
99 
100  if (W !=DetId(0)) {
101  NW = navigator.north();
102  } else {
103  navigator.home();
104  N=navigator.north();
105  NW=navigator.west();
106  }
107  associateNeighbour(NW,hit,hits,refProd,-1,1);
108  }
109 
110 
111 
112  protected:
113 
114  double sigmaCut2_;
115  std::unique_ptr<const T> topology_;
116  std::unique_ptr<CaloRecHitResolutionProvider> _timeResolutionCalc;
117 
118 
119  void associateNeighbour(const DetId& id, reco::PFRecHit& hit,std::unique_ptr<reco::PFRecHitCollection>& hits,edm::RefProd<reco::PFRecHitCollection>& refProd,short eta, short phi) {
120  double sigma2=10000.0;
121 
122 
123  auto found_hit = std::lower_bound(hits->begin(),hits->end(),
124  id,
125  [](const reco::PFRecHit& a,
126  DetId b){
127  if (DetId(a.detId()).det() == DetId::Hcal || b.det() == DetId::Hcal) return (HcalDetId)(a.detId()) < (HcalDetId)(b);
128 
129  else return a.detId() < b.rawId();
130  });
131 
132 
133  if (found_hit != hits->end() && ((found_hit->detId() == id.rawId()) ||
134  ((id.det()==DetId::Hcal) &&
135  ((HcalDetId)(found_hit->detId()) == (HcalDetId)(id))))) {
136  sigma2 = _timeResolutionCalc->timeResolution2(hit.energy()) + _timeResolutionCalc->timeResolution2(found_hit->energy());
137  const double deltaTime = hit.time()-found_hit->time();
138  if(deltaTime*deltaTime<sigmaCut2_*sigma2) {
139  hit.addNeighbour(eta,phi,0,found_hit-hits->begin());
140  }
141  }
142 
143  }
144 
145 
146 
147 };
148 
149 #endif
150 
151 
T getParameter(std::string const &) const
float time() const
timing for cleaned hits
Definition: PFRecHit.h:118
unsigned detId() const
rechit detId
Definition: PFRecHit.h:108
void addNeighbour(short x, short y, short z, unsigned int)
Definition: PFRecHit.cc:5
void associateNeighbour(const DetId &id, reco::PFRecHit &hit, std::unique_ptr< reco::PFRecHitCollection > &hits, edm::RefProd< reco::PFRecHitCollection > &refProd, short eta, short phi)
Particle flow rechit (rechit + geometry and topology information). See clustering algorithm in PFClus...
Definition: PFRecHit.h:31
std::unique_ptr< CaloRecHitResolutionProvider > _timeResolutionCalc
float energy() const
rechit energy
Definition: PFRecHit.h:114
Definition: DetId.h:18
#define N
Definition: blowfish.cc:9
ParameterSet const & getParameterSet(std::string const &) const
void associateNeighbours(reco::PFRecHit &hit, std::unique_ptr< reco::PFRecHitCollection > &hits, edm::RefProd< reco::PFRecHitCollection > &refProd)
double b
Definition: hdecay.h:120
double S(const TLorentzVector &, const TLorentzVector &)
Definition: Particle.cc:99
double a
Definition: hdecay.h:121
PFRecHitCaloNavigatorWithTime(const edm::ParameterSet &iConfig)
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40