CMS 3D CMS Logo

RealisticSimClusterMapper.cc
Go to the documentation of this file.
1 // Author: Felice Pantaleo
3 // Date: 30/06/2017
4 // Email: felice@cern.ch
6 
7 #include "RealisticCluster.h"
8 
19 
20 #include <unordered_map>
21 
23 public:
25  : InitialClusteringStepBase(conf, cc),
26  invisibleFraction_(conf.getParameter<double>("invisibleFraction")),
27  exclusiveFraction_(conf.getParameter<double>("exclusiveFraction")),
28  maxDistanceFilter_(conf.getParameter<bool>("maxDistanceFilter")),
29  maxDistance_(conf.getParameter<double>("maxDistance")),
30  maxDforTimingSquared_(conf.getParameter<double>("maxDforTimingSquared")),
31  timeOffset_(conf.getParameter<double>("timeOffset")),
32  minNHitsforTiming_(conf.getParameter<unsigned int>("minNHitsforTiming")),
33  useMCFractionsForExclEnergy_(conf.getParameter<bool>("useMCFractionsForExclEnergy")),
34  calibMinEta_(conf.getParameter<double>("calibMinEta")),
35  calibMaxEta_(conf.getParameter<double>("calibMaxEta")),
36  hadronCalib_(conf.getParameter<std::vector<double> >("hadronCalib")),
37  egammaCalib_(conf.getParameter<std::vector<double> >("egammaCalib")),
38  simClusterToken_(cc.consumes<SimClusterCollection>(conf.getParameter<edm::InputTag>("simClusterSrc"))),
39  geomToken_(cc.esConsumes<edm::Transition::BeginLuminosityBlock>()) {}
40 
44 
45  void updateEvent(const edm::Event&) final;
46  void update(const edm::EventSetup&) final;
47 
49  const std::vector<bool>&,
50  const std::vector<bool>&,
51  reco::PFClusterCollection&) override;
52 
53 private:
55  const float invisibleFraction_ = 0.3f;
56  const float exclusiveFraction_ = 0.7f;
57  const bool maxDistanceFilter_ = false;
58  const float maxDistance_ = 10.f;
59  const float maxDforTimingSquared_ = 4.0f;
60  const float timeOffset_;
61  const unsigned int minNHitsforTiming_ = 3;
62  const bool useMCFractionsForExclEnergy_ = false;
63  const float calibMinEta_ = 1.4;
64  const float calibMaxEta_ = 3.0;
65  std::vector<double> hadronCalib_;
66  std::vector<double> egammaCalib_;
67 
70 
72 };
73 
75 
76 #ifdef PFLOW_DEBUG
77 #define LOGVERB(x) edm::LogVerbatim(x)
78 #define LOGWARN(x) edm::LogWarning(x)
79 #define LOGERR(x) edm::LogError(x)
80 #define LOGDRESSED(x) edm::LogInfo(x)
81 #else
82 #define LOGVERB(x) LogTrace(x)
83 #define LOGWARN(x) edm::LogWarning(x)
84 #define LOGERR(x) edm::LogError(x)
85 #define LOGDRESSED(x) LogDebug(x)
86 #endif
87 
88 namespace {
89 
90  inline bool isPi0(int pdgId) { return pdgId == 111; }
91 
92  inline bool isEGamma(int pdgId) {
93  pdgId = std::abs(pdgId);
94  return (pdgId == 11) or (pdgId == 22);
95  }
96 
97  inline bool isHadron(int pdgId) {
98  pdgId = std::abs(pdgId) % 10000;
99  return ((pdgId > 100 and pdgId < 900) or (pdgId > 1000 and pdgId < 9000));
100  }
101 } // namespace
102 
104 
106 
108  const std::vector<bool>& rechitMask,
109  const std::vector<bool>& seedable,
111  const SimClusterCollection& simClusters = *simClusterH_;
112  auto const& hits = *input;
113  RealisticHitToClusterAssociator realisticAssociator;
116  realisticAssociator.init(hits.size(), simClusters.size(), numberOfLayers + 1);
117  // for quick indexing back to hit energy
118  std::unordered_map<uint32_t, size_t> detIdToIndex(hits.size());
119  for (uint32_t i = 0; i < hits.size(); ++i) {
120  detIdToIndex[hits[i].detId()] = i;
121  auto ref = makeRefhit(input, i);
122  const auto& hitPos = rhtools_.getPosition(ref->detId());
123 
124  realisticAssociator.insertHitPosition(hitPos.x(), hitPos.y(), hitPos.z(), i);
125  realisticAssociator.insertHitEnergy(ref->energy(), i);
126  realisticAssociator.insertLayerId(rhtools_.getLayerWithOffset(ref->detId()), i);
127  }
128 
129  for (unsigned int ic = 0; ic < simClusters.size(); ++ic) {
130  const auto& sc = simClusters[ic];
131  const auto& hitsAndFractions = sc.hits_and_fractions();
132  for (const auto& hAndF : hitsAndFractions) {
133  auto itr = detIdToIndex.find(hAndF.first);
134  if (itr == detIdToIndex.end()) {
135  continue; // hit wasn't saved in reco or did not pass the SNR threshold
136  }
137  auto hitId = itr->second;
138  auto ref = makeRefhit(input, hitId);
139  float fraction = hAndF.second;
140  float associatedEnergy = fraction * ref->energy();
141  realisticAssociator.insertSimClusterIdAndFraction(ic, fraction, hitId, associatedEnergy);
142  }
143  }
144  realisticAssociator.computeAssociation(
147  realisticAssociator.findCentersOfGravity();
148  if (maxDistanceFilter_)
149  realisticAssociator.filterHitsByDistance(maxDistance_);
150 
151  const auto& realisticClusters = realisticAssociator.realisticClusters();
152  unsigned int nClusters = realisticClusters.size();
153  float bin_norm = 1. / (calibMaxEta_ - calibMinEta_);
154  float egamma_bin_norm = egammaCalib_.size() * bin_norm;
155  float hadron_bin_norm = hadronCalib_.size() * bin_norm;
156  for (unsigned ic = 0; ic < nClusters; ++ic) {
157  float highest_energy = 0.0f;
158  output.emplace_back();
159  reco::PFCluster& back = output.back();
161  float energyCorrection = 1.f;
162  float timeRealisticSC = -99.;
163  if (realisticClusters[ic].isVisible()) {
164  int pdgId = simClusters[ic].pdgId();
165  auto abseta = std::abs(simClusters[ic].eta());
166  if ((abseta >= calibMinEta_) and (abseta <= calibMaxEta_)) //protecting range
167  {
168  if ((isEGamma(pdgId) or isPi0(pdgId)) and !egammaCalib_.empty()) {
169  unsigned int etabin = std::floor(((abseta - calibMinEta_) * egamma_bin_norm));
170  energyCorrection = egammaCalib_[etabin];
171  } else if (isHadron(pdgId) and !(isPi0(pdgId)) and
172  !hadronCalib_
173  .empty()) // this function is expensive.. should we treat as hadron everything which is not egamma?
174  {
175  unsigned int etabin = std::floor(((abseta - calibMinEta_) * hadron_bin_norm));
176  energyCorrection = hadronCalib_[etabin];
177  }
178  }
179  std::vector<float> timeHits;
180  const auto& hitsIdsAndFractions = realisticClusters[ic].hitsIdsAndFractions();
181  for (const auto& idAndF : hitsIdsAndFractions) {
182  auto fraction = idAndF.second;
183  if (fraction > 0.f) {
184  auto ref = makeRefhit(input, idAndF.first);
186  const float hit_energy = fraction * ref->energy();
187  if (hit_energy > highest_energy || highest_energy == 0.0) {
188  highest_energy = hit_energy;
189  seed = ref;
190  }
191  //select hits good for timing
192  if (ref->time() > -1.) {
193  int rhLayer = rhtools_.getLayerWithOffset(ref->detId());
194  std::array<float, 3> scPosition = realisticClusters[ic].getCenterOfGravity(rhLayer);
195  float distanceSquared =
196  std::pow((ref->position().x() - scPosition[0]), 2) + std::pow((ref->position().y() - scPosition[1]), 2);
197  if (distanceSquared < maxDforTimingSquared_) {
198  timeHits.push_back(ref->time());
199  }
200  }
201  }
202  }
203  //assign time if minimum number of hits
205  timeRealisticSC = (timeEstimator.fixSizeHighestDensity(timeHits)).first;
206  }
207  if (!back.hitsAndFractions().empty()) {
208  back.setSeed(seed->detId());
209  back.setEnergy(realisticClusters[ic].getEnergy());
210  back.setCorrectedEnergy(energyCorrection * realisticClusters[ic].getEnergy()); //applying energy correction
211  back.setTime(timeRealisticSC);
212  } else {
213  back.setSeed(-1);
214  back.setEnergy(0.f);
215  }
216  }
217 }
void insertLayerId(unsigned int layerId, unsigned int hitIndex)
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
const std::vector< std::pair< DetId, float > > & hitsAndFractions() const
Definition: CaloCluster.h:210
reco::PFRecHitRef makeRefhit(const edm::Handle< reco::PFRecHitCollection > &h, const unsigned i) const
unsigned int getLayer(DetId::Detector type, bool nose=false) const
Definition: RecHitTools.cc:307
Particle flow cluster, see clustering algorithm in PFClusterAlgo.
Definition: PFCluster.h:42
void buildClusters(const edm::Handle< reco::PFRecHitCollection > &, const std::vector< bool > &, const std::vector< bool > &, reco::PFClusterCollection &) override
constexpr uint32_t numberOfLayers
double getEnergy(HBHERecHitCollection::const_iterator hit, int useRaw=0, bool debug=false)
void setTime(float time, float timeError=0)
Definition: PFCluster.h:84
Fraction of a PFRecHit (rechits can be shared between several PFCluster&#39;s)
void setEnergy(double energy)
Definition: CaloCluster.h:136
int getGeometryType() const
Definition: RecHitTools.h:85
static std::string const input
Definition: EdmProvDump.cc:50
void setSeed(const DetId &id)
Definition: CaloCluster.h:146
RealisticSimClusterMapper(const edm::ParameterSet &conf, edm::ConsumesCollector &cc)
std::pair< float, float > fixSizeHighestDensity(std::vector< float > &time, std::vector< float > weight=std::vector< float >(), unsigned int minNhits=3, float deltaT=0.210, float timeWidthBy=0.5)
void init(std::size_t numberOfHits, std::size_t numberOfSimClusters, std::size_t numberOfLayers)
unsigned int lastLayerFH() const
Definition: RecHitTools.h:76
GlobalPoint getPosition(const DetId &id) const
Definition: RecHitTools.cc:129
const std::vector< RealisticCluster > & realisticClusters() const
void setCorrectedEnergy(double cenergy)
Definition: CaloCluster.h:137
void insertSimClusterIdAndFraction(unsigned int scIdx, float fraction, unsigned int hitIndex, float associatedEnergy)
void computeAssociation(float exclusiveFraction, bool useMCFractionsForExclEnergy, unsigned int fhOffset, unsigned int bhOffset)
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e< void, edm::EventID const &, edm::Timestamp const & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
Transition
Definition: Transition.h:12
void findAndMergeInvisibleClusters(float invisibleFraction, float exclusiveFraction)
void insertHitPosition(float x, float y, float z, unsigned int hitIndex)
double f[11][100]
bool getData(T &iHolder) const
Definition: EventSetup.h:122
RealisticSimClusterMapper & operator=(const RealisticSimClusterMapper &)=delete
void update(const edm::EventSetup &) final
edm::EDGetTokenT< SimClusterCollection > simClusterToken_
void insertHitEnergy(float energy, unsigned int hitIndex)
void updateEvent(const edm::Event &) final
void addRecHitFraction(const reco::PFRecHitFraction &frac)
add a given fraction of the rechit
Definition: PFCluster.cc:33
void setGeometry(CaloGeometry const &)
Definition: RecHitTools.cc:68
HLT enums.
std::vector< PFCluster > PFClusterCollection
collection of PFCluster objects
Definition: PFClusterFwd.h:9
#define DEFINE_EDM_PLUGIN(factory, type, name)
edm::Handle< SimClusterCollection > simClusterH_
std::vector< SimCluster > SimClusterCollection
Definition: SimClusterFwd.h:8
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
unsigned int lastLayerEE(bool nose=false) const
Definition: RecHitTools.h:75
edm::ESGetToken< CaloGeometry, CaloGeometryRecord > geomToken_
unsigned int getLayerWithOffset(const DetId &) const
Definition: RecHitTools.cc:365