CMS 3D CMS Logo

LegacyPFClusterProducer.cc
Go to the documentation of this file.
1 #include <Eigen/Core>
2 #include <Eigen/Dense>
3 #include <algorithm>
4 #include <array>
5 #include <cmath>
6 #include <iostream>
7 #include <memory>
8 #include <string>
9 #include <type_traits>
10 #include <unordered_map>
11 #include <utility>
12 #include <vector>
13 
24 
27 
35 
37 public:
39  : pfClusterSoAToken_(consumes(config.getParameter<edm::InputTag>("src"))),
40  pfRecHitFractionSoAToken_(consumes(config.getParameter<edm::InputTag>("src"))),
41  InputPFRecHitSoA_Token_{consumes(config.getParameter<edm::InputTag>("PFRecHitsLabelIn"))},
42  legacyPfClustersToken_(produces()),
43  recHitsLabel_(consumes(config.getParameter<edm::InputTag>("recHitsSource"))),
44  hcalCutsToken_(esConsumes<HcalPFCuts, HcalPFCutsRcd>(edm::ESInputTag("", "withTopo"))),
45  cutsFromDB_(config.getParameter<bool>("usePFThresholdsFromDB")) {
46  edm::ConsumesCollector cc = consumesCollector();
47 
48  //setup pf cluster builder if requested
49  const edm::ParameterSet& pfcConf = config.getParameterSet("pfClusterBuilder");
50  if (!pfcConf.empty()) {
51  const auto& acConf = pfcConf.getParameterSet("positionCalc");
52  if (!acConf.empty()) {
53  const std::string& algoac = acConf.getParameter<std::string>("algoName");
54  if (!algoac.empty())
55  positionCalc_ = PFCPositionCalculatorFactory::get()->create(algoac, acConf, cc);
56  }
57 
58  const auto& acConf2 = pfcConf.getParameterSet("allCellsPositionCalc");
59  if (!acConf2.empty()) {
60  const std::string& algoac = acConf2.getParameter<std::string>("algoName");
61  if (!algoac.empty())
62  allCellsPositionCalc_ = PFCPositionCalculatorFactory::get()->create(algoac, acConf2, cc);
63  }
64  }
65  }
66 
67  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
69  desc.add<edm::InputTag>("src", edm::InputTag("pfClusterSoAProducer"));
70  desc.add<edm::InputTag>("PFRecHitsLabelIn", edm::InputTag("pfRecHitSoAProducerHCAL"));
71  desc.add<edm::InputTag>("recHitsSource", edm::InputTag("legacyPFRecHitProducer"));
72  desc.add<bool>("usePFThresholdsFromDB", true);
73  {
75  pfClusterBuilder.add<unsigned int>("maxIterations", 5);
76  pfClusterBuilder.add<double>("minFracTot", 1e-20);
77  pfClusterBuilder.add<double>("minFractionToKeep", 1e-7);
78  pfClusterBuilder.add<bool>("excludeOtherSeeds", true);
79  pfClusterBuilder.add<double>("showerSigma", 10.);
80  pfClusterBuilder.add<double>("stoppingTolerance", 1e-8);
81  pfClusterBuilder.add<double>("timeSigmaEB", 10.);
82  pfClusterBuilder.add<double>("timeSigmaEE", 10.);
83  pfClusterBuilder.add<double>("maxNSigmaTime", 10.);
84  pfClusterBuilder.add<double>("minChi2Prob", 0.);
85  pfClusterBuilder.add<bool>("clusterTimeResFromSeed", false);
86  pfClusterBuilder.add<std::string>("algoName", "");
87  pfClusterBuilder.add<edm::ParameterSetDescription>("positionCalcForConvergence", {});
88  {
90  validator.add<std::string>("detector", "");
91  validator.add<std::vector<int>>("depths", {});
92  validator.add<std::vector<double>>("recHitEnergyNorm", {});
93  std::vector<edm::ParameterSet> vDefaults(2);
94  vDefaults[0].addParameter<std::string>("detector", "HCAL_BARREL1");
95  vDefaults[0].addParameter<std::vector<int>>("depths", {1, 2, 3, 4});
96  vDefaults[0].addParameter<std::vector<double>>("recHitEnergyNorm", {0.1, 0.2, 0.3, 0.3});
97  vDefaults[1].addParameter<std::string>("detector", "HCAL_ENDCAP");
98  vDefaults[1].addParameter<std::vector<int>>("depths", {1, 2, 3, 4, 5, 6, 7});
99  vDefaults[1].addParameter<std::vector<double>>("recHitEnergyNorm", {0.1, 0.2, 0.2, 0.2, 0.2, 0.2, 0.2});
100  pfClusterBuilder.addVPSet("recHitEnergyNorms", validator, vDefaults);
101  }
102  {
104  bar.add<std::string>("algoName", "Basic2DGenericPFlowPositionCalc");
105  bar.add<double>("minFractionInCalc", 1e-9);
106  bar.add<int>("posCalcNCrystals", 5);
107  {
109  validator.add<std::string>("detector", "");
110  validator.add<std::vector<int>>("depths", {});
111  validator.add<std::vector<double>>("logWeightDenominator", {});
112  std::vector<edm::ParameterSet> vDefaults(2);
113  vDefaults[0].addParameter<std::string>("detector", "HCAL_BARREL1");
114  vDefaults[0].addParameter<std::vector<int>>("depths", {1, 2, 3, 4});
115  vDefaults[0].addParameter<std::vector<double>>("logWeightDenominator", {0.1, 0.2, 0.3, 0.3});
116  vDefaults[1].addParameter<std::string>("detector", "HCAL_ENDCAP");
117  vDefaults[1].addParameter<std::vector<int>>("depths", {1, 2, 3, 4, 5, 6, 7});
118  vDefaults[1].addParameter<std::vector<double>>("logWeightDenominator", {0.1, 0.2, 0.2, 0.2, 0.2, 0.2, 0.2});
119  bar.addVPSet("logWeightDenominatorByDetector", validator, vDefaults);
120  }
121  bar.add<double>("minAllowedNormalization", 1e-9);
122  bar.add<edm::ParameterSetDescription>("timeResolutionCalcBarrel", {});
123  bar.add<edm::ParameterSetDescription>("timeResolutionCalcEndcap", {});
124  pfClusterBuilder.add("positionCalc", bar);
125  }
126  {
128  bar.add<std::string>("algoName", "Basic2DGenericPFlowPositionCalc");
129  bar.add<double>("minFractionInCalc", 1e-9);
130  bar.add<int>("posCalcNCrystals", -1);
131  {
133  validator.add<std::string>("detector", "");
134  validator.add<std::vector<int>>("depths", {});
135  validator.add<std::vector<double>>("logWeightDenominator", {});
136  std::vector<edm::ParameterSet> vDefaults(2);
137  vDefaults[0].addParameter<std::string>("detector", "HCAL_BARREL1");
138  vDefaults[0].addParameter<std::vector<int>>("depths", {1, 2, 3, 4});
139  vDefaults[0].addParameter<std::vector<double>>("logWeightDenominator", {0.1, 0.2, 0.3, 0.3});
140  vDefaults[1].addParameter<std::string>("detector", "HCAL_ENDCAP");
141  vDefaults[1].addParameter<std::vector<int>>("depths", {1, 2, 3, 4, 5, 6, 7});
142  vDefaults[1].addParameter<std::vector<double>>("logWeightDenominator", {0.1, 0.2, 0.2, 0.2, 0.2, 0.2, 0.2});
143  bar.addVPSet("logWeightDenominatorByDetector", validator, vDefaults);
144  }
145  bar.add<double>("minAllowedNormalization", 1e-9);
146  bar.add<edm::ParameterSetDescription>("timeResolutionCalcBarrel", {});
147  bar.add<edm::ParameterSetDescription>("timeResolutionCalcEndcap", {});
148  pfClusterBuilder.add("allCellsPositionCalc", bar);
149  }
150  {
152  bar.add<double>("corrTermLowE", 0.);
153  bar.add<double>("threshLowE", 6.);
154  bar.add<double>("noiseTerm", 21.86);
155  bar.add<double>("constantTermLowE", 4.24);
156  bar.add<double>("noiseTermLowE", 8.);
157  bar.add<double>("threshHighE", 15.);
158  bar.add<double>("constantTerm", 2.82);
159  pfClusterBuilder.add("timeResolutionCalcBarrel", bar);
160  }
161  {
163  bar.add<double>("corrTermLowE", 0.);
164  bar.add<double>("threshLowE", 6.);
165  bar.add<double>("noiseTerm", 21.86);
166  bar.add<double>("constantTermLowE", 4.24);
167  bar.add<double>("noiseTermLowE", 8.);
168  bar.add<double>("threshHighE", 15.);
169  bar.add<double>("constantTerm", 2.82);
170  pfClusterBuilder.add("timeResolutionCalcEndcap", bar);
171  }
172  {
174  pfClusterBuilder.add("positionReCalc", bar);
175  }
176  {
178  pfClusterBuilder.add("energyCorrector", bar);
179  }
180  desc.add("pfClusterBuilder", pfClusterBuilder);
181  }
182  descriptions.addWithDefaultLabel(desc);
183  }
184 
185 private:
186  void produce(edm::Event&, const edm::EventSetup&) override;
193  const bool cutsFromDB_;
194  // the actual algorithm
195  std::unique_ptr<PFCPositionCalculatorBase> positionCalc_;
196  std::unique_ptr<PFCPositionCalculatorBase> allCellsPositionCalc_;
197 };
198 
201 
202  HcalPFCuts const* paramPF = cutsFromDB_ ? &setup.getData(hcalCutsToken_) : nullptr;
203 
204  auto const& pfClusterSoA = event.get(pfClusterSoAToken_).const_view();
205  auto const& pfRecHitFractionSoA = event.get(pfRecHitFractionSoAToken_).const_view();
206 
207  int nRH = 0;
208  if (pfRecHits->metadata().size() != 0)
209  nRH = pfRecHits.view().size();
211  out.reserve(nRH);
212 
213  auto const rechitsHandle = event.getHandle(recHitsLabel_);
214 
215  if (nRH != 0) {
216  // Build PFClusters in legacy format
217  std::vector<int> nTopoSeeds(nRH, 0);
218 
219  for (int i = 0; i < pfClusterSoA.nSeeds(); i++) {
220  nTopoSeeds[pfClusterSoA[i].topoId()]++;
221  }
222 
223  // Looping over SoA PFClusters to produce legacy PFCluster collection
224  for (int i = 0; i < pfClusterSoA.nSeeds(); i++) {
225  unsigned int n = pfClusterSoA[i].seedRHIdx();
227  temp.setSeed((*rechitsHandle)[n].detId()); // Pulling the detId of this PFRecHit from the legacy format input
228  int offset = pfClusterSoA[i].rhfracOffset();
229  for (int k = offset; k < (offset + pfClusterSoA[i].rhfracSize()) && k >= 0;
230  k++) { // Looping over PFRecHits in the same topo cluster
231  if (pfRecHitFractionSoA[k].pfrhIdx() < nRH && pfRecHitFractionSoA[k].pfrhIdx() > -1 &&
232  pfRecHitFractionSoA[k].frac() > 0.0) {
233  const reco::PFRecHitRef& refhit = reco::PFRecHitRef(rechitsHandle, pfRecHitFractionSoA[k].pfrhIdx());
234  temp.addRecHitFraction(reco::PFRecHitFraction(refhit, pfRecHitFractionSoA[k].frac()));
235  }
236  }
237 
238  // Now PFRecHitFraction of this PFCluster is set. Now compute calculateAndSetPosition (energy, position etc)
239  if (nTopoSeeds[pfClusterSoA[i].topoId()] == 1 && allCellsPositionCalc_) {
240  allCellsPositionCalc_->calculateAndSetPosition(temp, paramPF);
241  } else {
242  positionCalc_->calculateAndSetPosition(temp, paramPF);
243  }
244  out.emplace_back(std::move(temp));
245  }
246  }
247 
248  event.emplace(legacyPfClustersToken_, std::move(out));
249 }
250 
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
std::unique_ptr< PFCPositionCalculatorBase > positionCalc_
const edm::EDGetTokenT< reco::PFClusterHostCollection > pfClusterSoAToken_
Particle flow cluster, see clustering algorithm in PFClusterAlgo.
Definition: PFCluster.h:42
const edm::EDGetTokenT< reco::PFRecHitFractionHostCollection > pfRecHitFractionSoAToken_
const edm::EDGetTokenT< reco::PFRecHitHostCollection > InputPFRecHitSoA_Token_
uint32_t cc[maxCellsPerHit]
Definition: gpuFishbone.h:49
ParameterSet const & getParameterSet(std::string const &) const
Fraction of a PFRecHit (rechits can be shared between several PFCluster&#39;s)
Definition: config.py:1
const edm::EDPutTokenT< reco::PFClusterCollection > legacyPfClustersToken_
std::unique_ptr< PFCPositionCalculatorBase > allCellsPositionCalc_
bool empty() const
Definition: ParameterSet.h:202
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
edm::Ref< PFRecHitCollection > PFRecHitRef
persistent reference to PFRecHit objects
Definition: PFRecHitFwd.h:15
LegacyPFClusterProducer(edm::ParameterSet const &config)
const edm::ESGetToken< HcalPFCuts, HcalPFCutsRcd > hcalCutsToken_
HLT enums.
const edm::EDGetTokenT< reco::PFRecHitCollection > recHitsLabel_
std::vector< PFCluster > PFClusterCollection
collection of PFCluster objects
Definition: PFClusterFwd.h:9
#define get
void produce(edm::Event &, const edm::EventSetup &) override
def move(src, dest)
Definition: eostools.py:511
Definition: event.py:1