CMS 3D CMS Logo

SiPixelDigisClustersFromSoA.cc
Go to the documentation of this file.
19 
20 // local include(s)
21 #include "PixelClusterizerBase.h"
23 
24 template <typename TrackerTraits>
26 public:
27  explicit SiPixelDigisClustersFromSoAT(const edm::ParameterSet& iConfig);
28  ~SiPixelDigisClustersFromSoAT() override = default;
29 
30  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
31 
32 private:
33  void produce(edm::StreamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const override;
34 
36 
38 
41 
42  const SiPixelClusterThresholds clusterThresholds_; // Cluster threshold in electrons
43 
44  const bool produceDigis_;
45  const bool storeDigis_;
46 };
47 
48 template <typename TrackerTraits>
50  : topoToken_(esConsumes()),
51  digiGetToken_(consumes<SiPixelDigisSoA>(iConfig.getParameter<edm::InputTag>("src"))),
52  clusterPutToken_(produces<SiPixelClusterCollectionNew>()),
53  clusterThresholds_(iConfig.getParameter<int>("clusterThreshold_layer1"),
54  iConfig.getParameter<int>("clusterThreshold_otherLayers")),
55  produceDigis_(iConfig.getParameter<bool>("produceDigis")),
56  storeDigis_(iConfig.getParameter<bool>("produceDigis") && iConfig.getParameter<bool>("storeDigis")) {
57  if (produceDigis_)
58  digiPutToken_ = produces<edm::DetSetVector<PixelDigi>>();
59 }
60 
61 template <typename TrackerTraits>
64  desc.add<edm::InputTag>("src", edm::InputTag("siPixelDigisSoA"));
65  desc.add<int>("clusterThreshold_layer1", gpuClustering::clusterThresholdLayerOne);
66  desc.add<int>("clusterThreshold_otherLayers", gpuClustering::clusterThresholdOtherLayers);
67  desc.add<bool>("produceDigis", true);
68  desc.add<bool>("storeDigis", true);
69 
70  descriptions.addWithDefaultLabel(desc);
71 }
72 
73 template <typename TrackerTraits>
76  const edm::EventSetup& iSetup) const {
77  const auto& digis = iEvent.get(digiGetToken_);
78  const uint32_t nDigis = digis.size();
79  const auto& ttopo = iSetup.getData(topoToken_);
80  constexpr auto maxModules = TrackerTraits::numberOfModules;
81 
82  std::unique_ptr<edm::DetSetVector<PixelDigi>> collection;
83  if (produceDigis_)
84  collection = std::make_unique<edm::DetSetVector<PixelDigi>>();
85  if (storeDigis_)
86  collection->reserve(maxModules);
87  auto outputClusters = std::make_unique<SiPixelClusterCollectionNew>();
88  outputClusters->reserve(maxModules, nDigis / 2);
89 
90  edm::DetSet<PixelDigi>* detDigis = nullptr;
91  uint32_t detId = 0;
92  for (uint32_t i = 0; i < nDigis; i++) {
93  // check for uninitialized digis
94  // this is set in RawToDigi_kernel in SiPixelRawToClusterGPUKernel.cu
95  if (digis.rawIdArr(i) == 0)
96  continue;
97  // check for noisy/dead pixels (electrons set to 0)
98  if (digis.adc(i) == 0)
99  continue;
100 
101  detId = digis.rawIdArr(i);
102  if (storeDigis_) {
103  detDigis = &collection->find_or_insert(detId);
104  if ((*detDigis).empty())
105  (*detDigis).data.reserve(64); // avoid the first relocations
106  }
107  break;
108  }
109 
110  int32_t nclus = -1;
112 #ifdef EDM_ML_DEBUG
113  auto totClustersFilled = 0;
114 #endif
115 
116  auto fillClusters = [&](uint32_t detId) {
117  if (nclus < 0)
118  return; // this in reality should never happen
120  auto layer = (DetId(detId).subdetId() == 1) ? ttopo.pxbLayer(detId) : 0;
121  auto clusterThreshold = clusterThresholds_.getThresholdForLayerOnCondition(layer == 1);
122  for (int32_t ic = 0; ic < nclus + 1; ++ic) {
123  auto const& acluster = aclusters[ic];
124  // in any case we cannot go out of sync with gpu...
126  edm::LogWarning("SiPixelDigisClustersFromSoA") << "cluster below charge Threshold "
127  << "Layer/DetId/clusId " << layer << '/' << detId << '/' << ic
128  << " size/charge " << acluster.isize << '/' << acluster.charge;
129  // sort by row (x)
130  spc.emplace_back(acluster.isize, acluster.adc, acluster.x, acluster.y, acluster.xmin, acluster.ymin, ic);
131  aclusters[ic].clear();
132 #ifdef EDM_ML_DEBUG
133  ++totClustersFilled;
134  const auto& cluster{spc.back()};
135  LogDebug("SiPixelDigisClustersFromSoA")
136  << "putting in this cluster " << ic << " " << cluster.charge() << " " << cluster.pixelADC().size();
137 #endif
138  std::push_heap(spc.begin(), spc.end(), [](SiPixelCluster const& cl1, SiPixelCluster const& cl2) {
139  return cl1.minPixelRow() < cl2.minPixelRow();
140  });
141  }
142  nclus = -1;
143  // sort by row (x)
144  std::sort_heap(spc.begin(), spc.end(), [](SiPixelCluster const& cl1, SiPixelCluster const& cl2) {
145  return cl1.minPixelRow() < cl2.minPixelRow();
146  });
147  if (spc.empty())
148  spc.abort();
149  };
150 
151  for (uint32_t i = 0; i < nDigis; i++) {
152  // check for uninitialized digis
153  if (digis.rawIdArr(i) == 0)
154  continue;
155  // check for noisy/dead pixels (electrons set to 0)
156  if (digis.adc(i) == 0)
157  continue;
158  if (digis.clus(i) > 9000)
159  continue; // not in cluster; TODO add an assert for the size
160 #ifdef EDM_ML_DEBUG
161  assert(digis.rawIdArr(i) > 109999);
162 #endif
163  if (detId != digis.rawIdArr(i)) {
164  // new module
165  fillClusters(detId);
166 #ifdef EDM_ML_DEBUG
167  assert(nclus == -1);
168 #endif
169  detId = digis.rawIdArr(i);
170  if (storeDigis_) {
171  detDigis = &collection->find_or_insert(detId);
172  if ((*detDigis).empty())
173  (*detDigis).data.reserve(64); // avoid the first relocations
174  else {
175  edm::LogWarning("SiPixelDigisClustersFromSoA")
176  << "Problem det present twice in input! " << (*detDigis).detId();
177  }
178  }
179  }
180  PixelDigi dig(digis.pdigi(i));
181  if (storeDigis_)
182  (*detDigis).data.emplace_back(dig);
183  // fill clusters
184 #ifdef EDM_ML_DEBUG
185  assert(digis.clus(i) >= 0);
187 #endif
188  nclus = std::max(digis.clus(i), nclus);
189  auto row = dig.row();
190  auto col = dig.column();
191  SiPixelCluster::PixelPos pix(row, col);
192  aclusters[digis.clus(i)].add(pix, digis.adc(i));
193  }
194 
195  // fill final clusters
196  if (detId > 0)
197  fillClusters(detId);
198 
199 #ifdef EDM_ML_DEBUG
200  LogDebug("SiPixelDigisClustersFromSoA") << "filled " << totClustersFilled << " clusters";
201 #endif
202 
203  if (produceDigis_)
204  iEvent.put(digiPutToken_, std::move(collection));
205  iEvent.put(clusterPutToken_, std::move(outputClusters));
206 }
207 
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
void produce(edm::StreamID, edm::Event &iEvent, const edm::EventSetup &iSetup) const override
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
~SiPixelDigisClustersFromSoAT() override=default
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
assert(be >=bs)
const edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > topoToken_
constexpr uint16_t numberOfModules
int minPixelRow() const
edm::EDGetTokenT< SiPixelDigisSoA > digiGetToken_
constexpr uint32_t maxNumClustersPerModules
int iEvent
Definition: GenABIO.cc:224
constexpr uint16_t clusterThresholdOtherLayers
const SiPixelClusterThresholds clusterThresholds_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:48
constexpr uint16_t clusterThresholdLayerOne
SiPixelDigisClustersFromSoAT(const edm::ParameterSet &iConfig)
Definition: DetId.h:17
Pixel cluster – collection of neighboring pixels above threshold.
HLT enums.
void emplace_back(Args &&... args)
edm::EDPutTokenT< edm::DetSetVector< PixelDigi > > digiPutToken_
col
Definition: cuy.py:1009
bool add(SiPixelCluster::PixelPos const &p, uint16_t const iadc)
Log< level::Warning, false > LogWarning
edm::EDPutTokenT< SiPixelClusterCollectionNew > clusterPutToken_
def move(src, dest)
Definition: eostools.py:511
#define LogDebug(id)