CMS 3D CMS Logo

PFHFRecHitCreator.h
Go to the documentation of this file.
1 #ifndef RecoParticleFlow_PFClusterProducer_PFHFRecHitCreator_h
2 #define RecoParticleFlow_PFClusterProducer_PFHFRecHitCreator_h
3 
5 
11 
18 
20 
21 class PFHFRecHitCreator final : public PFRecHitCreatorBase {
22 public:
25  EM_Depth_ = iConfig.getParameter<double>("EMDepthCorrection");
26  HAD_Depth_ = iConfig.getParameter<double>("HADDepthCorrection");
27  shortFibre_Cut = iConfig.getParameter<double>("ShortFibre_Cut");
28  longFibre_Fraction = iConfig.getParameter<double>("LongFibre_Fraction");
29  longFibre_Cut = iConfig.getParameter<double>("LongFibre_Cut");
30  shortFibre_Fraction = iConfig.getParameter<double>("ShortFibre_Fraction");
31  thresh_HF_ = iConfig.getParameter<double>("thresh_HF");
32  HFCalib_ = iConfig.getParameter<double>("HFCalib29");
33  }
34 
35  void importRecHits(std::unique_ptr<reco::PFRecHitCollection>& out,
36  std::unique_ptr<reco::PFRecHitCollection>& cleaned,
37  const edm::Event& iEvent,
38  const edm::EventSetup& iSetup) override {
40 
41  beginEvent(iEvent, iSetup);
42 
44 
46  iSetup.get<CaloGeometryRecord>().get(geoHandle);
47 
48  // get the ecal geometry
50 
51  iEvent.getByToken(recHitToken_, recHitHandle);
52  for (const auto& erh : *recHitHandle) {
53  const HcalDetId& detid = (HcalDetId)erh.detid();
54  auto depth = detid.depth();
55 
56  // ATTN: skip dual anode in HF for now (should be fixed in upstream changes)
57  if (depth > 2)
58  continue;
59 
60  auto energy = erh.energy();
61  auto time = erh.time();
62 
63  std::shared_ptr<const CaloCellGeometry> thisCell = hcalGeo->getGeometry(detid);
64  auto zp = dynamic_cast<IdealZPrism const*>(thisCell.get());
65  assert(zp);
66  thisCell = zp->forPF();
67 
68  // find rechit geometry
69  if (!thisCell) {
70  edm::LogError("PFHFRecHitCreator")
71  << "warning detid " << detid.rawId() << " not found in geometry" << std::endl;
72  continue;
73  }
74 
76 
77  reco::PFRecHit rh(thisCell, detid.rawId(), layer, energy);
78  rh.setTime(time);
79  rh.setDepth(depth);
80 
81  bool rcleaned = false;
82  bool keep = true;
83 
84  //Apply Q tests
85  for (const auto& qtest : qualityTests_) {
86  if (!qtest->test(rh, erh, rcleaned)) {
87  keep = false;
88  }
89  }
90 
91  if (keep) {
92  tmpOut.push_back(std::move(rh));
93  } else if (rcleaned)
94  cleaned->push_back(std::move(rh));
95  }
96  //Sort by DetID the collection
98  if (!tmpOut.empty())
99  std::sort(tmpOut.begin(), tmpOut.end(), sorter);
100 
102 
103  for (auto& hit : tmpOut) {
104  reco::PFRecHit newHit = hit;
105  const HcalDetId& detid = (HcalDetId)hit.detId();
106  if (detid.depth() == 1) {
107  double lONG = hit.energy();
108  //find the short hit
109  HcalDetId shortID(HcalForward, detid.ieta(), detid.iphi(), 2);
110  auto found_hit =
111  std::lower_bound(tmpOut.begin(), tmpOut.end(), shortID, [](const reco::PFRecHit& a, HcalDetId b) {
112  return a.detId() < b.rawId();
113  });
114  if (found_hit != tmpOut.end() && found_hit->detId() == shortID.rawId()) {
115  double sHORT = found_hit->energy();
116  //Ask for fraction
117  double energy = lONG - sHORT;
118 
119  if (abs(detid.ieta()) <= 32)
120  energy *= HFCalib_;
121  newHit.setEnergy(energy);
122  if (!(lONG > longFibre_Cut && (sHORT / lONG < shortFibre_Fraction)))
123  if (energy > thresh_HF_)
124  out->push_back(newHit);
125  } else {
126  //make only long hit
127  double energy = lONG;
128  if (abs(detid.ieta()) <= 32)
129  energy *= HFCalib_;
130  newHit.setEnergy(energy);
131 
132  if (energy > thresh_HF_)
133  out->push_back(newHit);
134  }
135 
136  } else {
137  double sHORT = hit.energy();
138  HcalDetId longID(HcalForward, detid.ieta(), detid.iphi(), 1);
139  auto found_hit =
140  std::lower_bound(tmpOut.begin(), tmpOut.end(), longID, [](const reco::PFRecHit& a, HcalDetId b) {
141  return a.detId() < b.rawId();
142  });
143  double energy = 2 * sHORT;
144  if (found_hit != tmpOut.end() && found_hit->detId() == longID.rawId()) {
145  double lONG = found_hit->energy();
146  //Ask for fraction
147 
148  //If in this case lONG-sHORT<0 add the energy to the sHORT
149  if ((lONG - sHORT) < thresh_HF_)
150  energy = lONG + sHORT;
151 
152  if (abs(detid.ieta()) <= 32)
153  energy *= HFCalib_;
154 
155  newHit.setEnergy(energy);
156  if (!(sHORT > shortFibre_Cut && (lONG / sHORT < longFibre_Fraction)))
157  if (energy > thresh_HF_)
158  out->push_back(newHit);
159 
160  } else {
161  //only short hit!
162  if (abs(detid.ieta()) <= 32)
163  energy *= HFCalib_;
164  newHit.setEnergy(energy);
165  if (energy > thresh_HF_)
166  out->push_back(newHit);
167  }
168  }
169  }
170  }
171 
172 protected:
174  double EM_Depth_;
175  double HAD_Depth_;
176  // Don't allow large energy in short fibres if there is no energy in long fibres
179 
180  // Don't allow large energy in long fibres if there is no energy in short fibres
183  double thresh_HF_;
184  double HFCalib_;
185 
186  class DetIDSorter {
187  public:
188  DetIDSorter() = default;
189  ~DetIDSorter() = default;
190 
191  bool operator()(const reco::PFRecHit& a, const reco::PFRecHit& b) {
192  if (DetId(a.detId()).det() == DetId::Hcal || DetId(b.detId()).det() == DetId::Hcal)
193  return (HcalDetId)(a.detId()) < (HcalDetId)(b.detId());
194  else
195  return a.detId() < b.detId();
196  }
197  };
198 };
199 #endif
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
T getParameter(std::string const &) const
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
access the subdetector geometry for the given subdetector directly
Definition: CaloGeometry.cc:34
void importRecHits(std::unique_ptr< reco::PFRecHitCollection > &out, std::unique_ptr< reco::PFRecHitCollection > &cleaned, const edm::Event &iEvent, const edm::EventSetup &iSetup) override
std::vector< std::unique_ptr< PFRecHitQTestBase > > qualityTests_
PFHFRecHitCreator(const edm::ParameterSet &iConfig, edm::ConsumesCollector &iC)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
unsigned detId() const
rechit detId
Definition: PFRecHit.h:89
edm::EDGetTokenT< edm::SortedCollection< HFRecHit > > recHitToken_
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
std::vector< PFRecHit > PFRecHitCollection
collection of PFRecHit objects
Definition: PFRecHitFwd.h:9
bool operator()(const reco::PFRecHit &a, const reco::PFRecHit &b)
const int keep
int depth() const
get the tower depth
Definition: HcalDetId.h:164
int iEvent
Definition: GenABIO.cc:224
Particle flow rechit (rechit + geometry and topology information). See clustering algorithm in PFClus...
Definition: PFRecHit.h:31
int ieta() const
get the cell ieta
Definition: HcalDetId.h:155
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
void beginEvent(const edm::Event &event, const edm::EventSetup &setup)
void setTime(double time)
Definition: PFRecHit.h:69
Layer
layer definition
Definition: PFLayer.h:29
int iphi() const
get the cell iphi
Definition: HcalDetId.h:157
Definition: DetId.h:17
virtual std::shared_ptr< const CaloCellGeometry > getGeometry(const DetId &id) const
Get the cell geometry of a given detector id. Should return false if not found.
double b
Definition: hdecay.h:118
double a
Definition: hdecay.h:119
void setEnergy(float energy)
Definition: PFRecHit.h:65
T get() const
Definition: EventSetup.h:73
def move(src, dest)
Definition: eostools.py:511