CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
JetCoreClusterSplitter Class Reference
Inheritance diagram for JetCoreClusterSplitter:
edm::stream::EDProducer<>

Public Member Functions

 JetCoreClusterSplitter (const edm::ParameterSet &iConfig)
 
void produce (edm::Event &iEvent, const edm::EventSetup &iSetup) override
 
 ~JetCoreClusterSplitter () override
 
- Public Member Functions inherited from edm::stream::EDProducer<>
 EDProducer ()=default
 
 EDProducer (const EDProducer &)=delete
 
bool hasAbilityToProduceInBeginLumis () const final
 
bool hasAbilityToProduceInBeginProcessBlocks () const final
 
bool hasAbilityToProduceInBeginRuns () const final
 
bool hasAbilityToProduceInEndLumis () const final
 
bool hasAbilityToProduceInEndProcessBlocks () const final
 
bool hasAbilityToProduceInEndRuns () const final
 
const EDProduceroperator= (const EDProducer &)=delete
 

Private Member Functions

std::pair< float, float > closestClusters (const std::vector< float > &distanceMap)
 
std::multimap< float, int > distScore (const std::vector< std::vector< float >> &distanceMap)
 
std::vector< SiPixelClusterfittingSplit (const SiPixelCluster &aCluster, float expectedADC, int sizeY, int sizeX, float jetZOverRho, unsigned int nSplitted)
 
std::multimap< float, int > secondDistDiffScore (const std::vector< std::vector< float >> &distanceMap)
 
std::multimap< float, int > secondDistScore (const std::vector< std::vector< float >> &distanceMap)
 
bool split (const SiPixelCluster &aCluster, edmNew::DetSetVector< SiPixelCluster >::FastFiller &filler, float expectedADC, int sizeY, int sizeX, float jetZOverRho)
 

Private Attributes

double centralMIPCharge_
 
double chargeFracMin_
 
double chargePerUnit_
 
edm::EDGetTokenT< edm::View< reco::Candidate > > cores_
 
double deltaR_
 
double forceXError_
 
double forceYError_
 
double fractionalWidth_
 
edm::EDGetTokenT< edmNew::DetSetVector< SiPixelCluster > > pixelClusters_
 
double ptMin_
 
edm::ESGetToken< PixelClusterParameterEstimator, TkPixelCPERecord > const tCPE_
 
edm::ESGetToken< GlobalTrackingGeometry, GlobalTrackingGeometryRecord > const tTrackingGeom_
 
bool verbose
 
edm::EDGetTokenT< reco::VertexCollectionvertices_
 

Additional Inherited Members

- Public Types inherited from edm::stream::EDProducer<>
using CacheTypes = CacheContexts< T... >
 
using GlobalCache = typename CacheTypes::GlobalCache
 
using HasAbility = AbilityChecker< T... >
 
using InputProcessBlockCache = typename CacheTypes::InputProcessBlockCache
 
using LuminosityBlockCache = typename CacheTypes::LuminosityBlockCache
 
using LuminosityBlockContext = LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCache >
 
using LuminosityBlockSummaryCache = typename CacheTypes::LuminosityBlockSummaryCache
 
using RunCache = typename CacheTypes::RunCache
 
using RunContext = RunContextT< RunCache, GlobalCache >
 
using RunSummaryCache = typename CacheTypes::RunSummaryCache
 

Detailed Description

Definition at line 26 of file JetCoreClusterSplitter.cc.

Constructor & Destructor Documentation

◆ JetCoreClusterSplitter()

JetCoreClusterSplitter::JetCoreClusterSplitter ( const edm::ParameterSet iConfig)

Definition at line 67 of file JetCoreClusterSplitter.cc.

69  tCPE_(esConsumes(edm::ESInputTag("", iConfig.getParameter<std::string>("pixelCPE")))),
70  verbose(iConfig.getParameter<bool>("verbose")),
71  ptMin_(iConfig.getParameter<double>("ptMin")),
72  deltaR_(iConfig.getParameter<double>("deltaRmax")),
73  chargeFracMin_(iConfig.getParameter<double>("chargeFractionMin")),
75  consumes<edmNew::DetSetVector<SiPixelCluster>>(iConfig.getParameter<edm::InputTag>("pixelClusters"))),
76  vertices_(consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("vertices"))),
77  cores_(consumes<edm::View<reco::Candidate>>(iConfig.getParameter<edm::InputTag>("cores"))),
78  forceXError_(iConfig.getParameter<double>("forceXError")),
79  forceYError_(iConfig.getParameter<double>("forceYError")),
80  fractionalWidth_(iConfig.getParameter<double>("fractionalWidth")),
81  chargePerUnit_(iConfig.getParameter<double>("chargePerUnit")),
82  centralMIPCharge_(iConfig.getParameter<double>("centralMIPCharge"))
83 
84 {
85  produces<edmNew::DetSetVector<SiPixelCluster>>();
86 }
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
edm::EDGetTokenT< edmNew::DetSetVector< SiPixelCluster > > pixelClusters_
edm::EDGetTokenT< reco::VertexCollection > vertices_
edm::EDGetTokenT< edm::View< reco::Candidate > > cores_
edm::ESGetToken< PixelClusterParameterEstimator, TkPixelCPERecord > const tCPE_
edm::ESGetToken< GlobalTrackingGeometry, GlobalTrackingGeometryRecord > const tTrackingGeom_

◆ ~JetCoreClusterSplitter()

JetCoreClusterSplitter::~JetCoreClusterSplitter ( )
override

Definition at line 88 of file JetCoreClusterSplitter.cc.

88 {}

Member Function Documentation

◆ closestClusters()

std::pair< float, float > JetCoreClusterSplitter::closestClusters ( const std::vector< float > &  distanceMap)
private

Definition at line 207 of file JetCoreClusterSplitter.cc.

References mps_fire::i, and SiStripPI::max.

Referenced by distScore(), secondDistDiffScore(), and secondDistScore().

207  {
208  float minDist = std::numeric_limits<float>::max();
209  float secondMinDist = std::numeric_limits<float>::max();
210  for (unsigned int i = 0; i < distanceMap.size(); i++) {
211  float dist = distanceMap[i];
212  if (dist < minDist) {
213  secondMinDist = minDist;
214  minDist = dist;
215  } else if (dist < secondMinDist) {
216  secondMinDist = dist;
217  }
218  }
219  return std::pair<float, float>(minDist, secondMinDist);
220 }

◆ distScore()

std::multimap< float, int > JetCoreClusterSplitter::distScore ( const std::vector< std::vector< float >> &  distanceMap)
private

Definition at line 241 of file JetCoreClusterSplitter.cc.

References closestClusters(), ztail::d, and dqmiolumiharvest::j.

241  {
242  std::multimap<float, int> scores;
243  for (unsigned int j = 0; j < distanceMap.size(); j++) {
244  std::pair<float, float> d = closestClusters(distanceMap[j]);
245  scores.insert(std::pair<float, int>(-d.first, j));
246  }
247  return scores;
248 }
d
Definition: ztail.py:151
std::pair< float, float > closestClusters(const std::vector< float > &distanceMap)

◆ fittingSplit()

std::vector< SiPixelCluster > JetCoreClusterSplitter::fittingSplit ( const SiPixelCluster aCluster,
float  expectedADC,
int  sizeY,
int  sizeX,
float  jetZOverRho,
unsigned int  nSplitted 
)
private

Definition at line 250 of file JetCoreClusterSplitter.cc.

References funct::abs(), gpuClustering::adc, centralMIPCharge_, chargePerUnit_, haddnano::cl, gather_cfg::cout, MillePedeFileConverter_cfg::e, relativeConstraints::empty, JetChargeProducer_cfi::exp, f, forceXError_, forceYError_, fractionalWidth_, mps_fire::i, createfilelist::int, dqmiolumiharvest::j, dqmdumpme::k, convertSQLitetoXML_cfg::output, AlCaHLTBitMon_ParallelJobs::p, SiPixelCluster::pixels(), edm::second(), secondDistScore(), RecoTauValidation_cfi::sizeX, RecoTauValidation_cfi::sizeY, mathSSE::sqrt(), verbose, x, and y.

Referenced by split().

255  {
256  std::vector<SiPixelCluster> output;
257 
258  unsigned int meanExp = nSplitted;
259  if (meanExp <= 1) {
260  output.push_back(aCluster);
261  return output;
262  }
263 
264  std::vector<float> clx(meanExp);
265  std::vector<float> cly(meanExp);
266  std::vector<float> cls(meanExp);
267  std::vector<float> oldclx(meanExp);
268  std::vector<float> oldcly(meanExp);
269  std::vector<SiPixelCluster::Pixel> originalpixels = aCluster.pixels();
270  std::vector<std::pair<int, SiPixelCluster::Pixel>> pixels;
271  for (unsigned int j = 0; j < originalpixels.size(); j++) {
272  int sub = originalpixels[j].adc / chargePerUnit_ * expectedADC / centralMIPCharge_;
273  if (sub < 1)
274  sub = 1;
275  int perDiv = originalpixels[j].adc / sub;
276  if (verbose)
277  std::cout << "Splitting " << j << " in [ " << pixels.size() << " , " << pixels.size() + sub
278  << " ], expected numb of clusters: " << meanExp << " original pixel (x,y) " << originalpixels[j].x
279  << " " << originalpixels[j].y << " sub " << sub << std::endl;
280  for (int k = 0; k < sub; k++) {
281  if (k == sub - 1)
282  perDiv = originalpixels[j].adc - perDiv * k;
283  pixels.push_back(std::make_pair(j, SiPixelCluster::Pixel(originalpixels[j].x, originalpixels[j].y, perDiv)));
284  }
285  }
286  std::vector<int> clusterForPixel(pixels.size());
287  // initial values
288  for (unsigned int j = 0; j < meanExp; j++) {
289  oldclx[j] = -999;
290  oldcly[j] = -999;
291  clx[j] = originalpixels[0].x + j;
292  cly[j] = originalpixels[0].y + j;
293  cls[j] = 0;
294  }
295  bool stop = false;
296  int remainingSteps = 100;
297  while (!stop && remainingSteps > 0) {
298  remainingSteps--;
299  // Compute all distances
300  std::vector<std::vector<float>> distanceMapX(originalpixels.size(), std::vector<float>(meanExp));
301  std::vector<std::vector<float>> distanceMapY(originalpixels.size(), std::vector<float>(meanExp));
302  std::vector<std::vector<float>> distanceMap(originalpixels.size(), std::vector<float>(meanExp));
303  for (unsigned int j = 0; j < originalpixels.size(); j++) {
304  if (verbose)
305  std::cout << "Original Pixel pos " << j << " " << pixels[j].second.x << " " << pixels[j].second.y << std::endl;
306  for (unsigned int i = 0; i < meanExp; i++) {
307  distanceMapX[j][i] = 1.f * originalpixels[j].x - clx[i];
308  distanceMapY[j][i] = 1.f * originalpixels[j].y - cly[i];
309  float dist = 0;
310  // float sizeX=2;
311  if (std::abs(distanceMapX[j][i]) > sizeX / 2.f) {
312  dist +=
313  (std::abs(distanceMapX[j][i]) - sizeX / 2.f + 1.f) * (std::abs(distanceMapX[j][i]) - sizeX / 2.f + 1.f);
314  } else {
315  dist += (2.f * distanceMapX[j][i] / sizeX) * (2.f * distanceMapX[j][i] / sizeX);
316  }
317 
318  if (std::abs(distanceMapY[j][i]) > sizeY / 2.f) {
319  dist += 1.f * (std::abs(distanceMapY[j][i]) - sizeY / 2.f + 1.f) *
320  (std::abs(distanceMapY[j][i]) - sizeY / 2.f + 1.f);
321  } else {
322  dist += 1.f * (2.f * distanceMapY[j][i] / sizeY) * (2.f * distanceMapY[j][i] / sizeY);
323  }
324  distanceMap[j][i] = sqrt(dist);
325  if (verbose)
326  std::cout << "Cluster " << i << " Original Pixel " << j << " distances: " << distanceMapX[j][i] << " "
327  << distanceMapY[j][i] << " " << distanceMap[j][i] << std::endl;
328  }
329  }
330  // Compute scores for sequential addition. The first index is the
331  // distance, in whatever metrics we use, while the second is the
332  // pixel index w.r.t which the distance is computed.
333  std::multimap<float, int> scores;
334 
335  // Using different rankings to improve convergence (as Giulio proposed)
336  scores = secondDistScore(distanceMap);
337 
338  // Iterate starting from the ones with furthest second best clusters, i.e.
339  // easy choices
340  std::vector<float> weightOfPixel(pixels.size());
341  for (std::multimap<float, int>::iterator it = scores.begin(); it != scores.end(); it++) {
342  int pixel_index = it->second;
343  if (verbose)
344  std::cout << "Original Pixel " << pixel_index << " with score " << it->first << std::endl;
345  // find cluster that is both close and has some charge still to assign
346  int subpixel_counter = 0;
347  for (auto subpixel = pixels.begin(); subpixel != pixels.end(); ++subpixel, ++subpixel_counter) {
348  if (subpixel->first > pixel_index) {
349  break;
350  } else if (subpixel->first != pixel_index) {
351  continue;
352  } else {
353  float maxEst = 0;
354  int cl = -1;
355  for (unsigned int subcluster_index = 0; subcluster_index < meanExp; subcluster_index++) {
356  float nsig = (cls[subcluster_index] - expectedADC) /
357  (expectedADC * fractionalWidth_); // 20% uncertainty? realistic from Landau?
358  float clQest = 1.f / (1.f + std::exp(nsig)) + 1e-6f; // 1./(1.+exp(x*x-3*3))
359  float clDest = 1.f / (distanceMap[pixel_index][subcluster_index] + 0.05f);
360 
361  if (verbose)
362  std::cout << " Q: " << clQest << " D: " << clDest << " " << distanceMap[pixel_index][subcluster_index]
363  << std::endl;
364  float est = clQest * clDest;
365  if (est > maxEst) {
366  cl = subcluster_index;
367  maxEst = est;
368  }
369  }
370  cls[cl] += subpixel->second.adc;
371  clusterForPixel[subpixel_counter] = cl;
372  weightOfPixel[subpixel_counter] = maxEst;
373  if (verbose)
374  std::cout << "Pixel weight j cl " << weightOfPixel[subpixel_counter] << " " << subpixel_counter << " " << cl
375  << std::endl;
376  }
377  }
378  }
379  // Recompute cluster centers
380  stop = true;
381  for (unsigned int subcluster_index = 0; subcluster_index < meanExp; subcluster_index++) {
382  if (std::abs(clx[subcluster_index] - oldclx[subcluster_index]) > 0.01f)
383  stop = false; // still moving
384  if (std::abs(cly[subcluster_index] - oldcly[subcluster_index]) > 0.01f)
385  stop = false;
386  oldclx[subcluster_index] = clx[subcluster_index];
387  oldcly[subcluster_index] = cly[subcluster_index];
388  clx[subcluster_index] = 0;
389  cly[subcluster_index] = 0;
390  cls[subcluster_index] = 1e-99;
391  }
392  for (unsigned int pixel_index = 0; pixel_index < pixels.size(); pixel_index++) {
393  if (clusterForPixel[pixel_index] < 0)
394  continue;
395  if (verbose)
396  std::cout << "j " << pixel_index << " x " << pixels[pixel_index].second.x << " * y "
397  << pixels[pixel_index].second.y << " * ADC " << pixels[pixel_index].second.adc << " * W "
398  << weightOfPixel[pixel_index] << std::endl;
399  clx[clusterForPixel[pixel_index]] += pixels[pixel_index].second.x * pixels[pixel_index].second.adc;
400  cly[clusterForPixel[pixel_index]] += pixels[pixel_index].second.y * pixels[pixel_index].second.adc;
401  cls[clusterForPixel[pixel_index]] += pixels[pixel_index].second.adc;
402  }
403  for (unsigned int subcluster_index = 0; subcluster_index < meanExp; subcluster_index++) {
404  if (cls[subcluster_index] != 0) {
405  clx[subcluster_index] /= cls[subcluster_index];
406  cly[subcluster_index] /= cls[subcluster_index];
407  }
408  if (verbose)
409  std::cout << "Center for cluster " << subcluster_index << " x,y " << clx[subcluster_index] << " "
410  << cly[subcluster_index] << std::endl;
411  cls[subcluster_index] = 0;
412  }
413  }
414  if (verbose)
415  std::cout << "maxstep " << remainingSteps << std::endl;
416  // accumulate pixel with same cl
417  std::vector<std::vector<SiPixelCluster::Pixel>> pixelsForCl(meanExp);
418  for (int cl = 0; cl < (int)meanExp; cl++) {
419  for (unsigned int j = 0; j < pixels.size(); j++) {
420  if (clusterForPixel[j] == cl and pixels[j].second.adc != 0) { // for each pixel of cluster
421  // cl find the other pixels
422  // with same x,y and
423  // accumulate+reset their adc
424  for (unsigned int k = j + 1; k < pixels.size(); k++) {
425  if (pixels[k].second.adc != 0 and pixels[k].second.x == pixels[j].second.x and
426  pixels[k].second.y == pixels[j].second.y and clusterForPixel[k] == cl) {
427  if (verbose)
428  std::cout << "Resetting all sub-pixel for location " << pixels[k].second.x << ", " << pixels[k].second.y
429  << " at index " << k << " associated to cl " << clusterForPixel[k] << std::endl;
430  pixels[j].second.adc += pixels[k].second.adc;
431  pixels[k].second.adc = 0;
432  }
433  }
434  for (unsigned int p = 0; p < pixels.size(); ++p)
435  if (verbose)
436  std::cout << "index, x, y, ADC: " << p << ", " << pixels[p].second.x << ", " << pixels[p].second.y << ", "
437  << pixels[p].second.adc << " associated to cl " << clusterForPixel[p] << std::endl
438  << "Adding pixel " << pixels[j].second.x << ", " << pixels[j].second.y << " to cluster " << cl
439  << std::endl;
440  pixelsForCl[cl].push_back(pixels[j].second);
441  }
442  }
443  }
444 
445  // std::vector<std::vector<std::vector<SiPixelCluster::PixelPos *> > >
446  //pixelMap(meanExp,std::vector<std::vector<SiPixelCluster::PixelPos *>
447  //>(512,std::vector<SiPixelCluster::Pixel *>(512,0)));
448 
449  for (int cl = 0; cl < (int)meanExp; cl++) {
450  if (verbose)
451  std::cout << "Pixels of cl " << cl << " ";
452  for (unsigned int j = 0; j < pixelsForCl[cl].size(); j++) {
453  SiPixelCluster::PixelPos newpix(pixelsForCl[cl][j].x, pixelsForCl[cl][j].y);
454  if (verbose)
455  std::cout << pixelsForCl[cl][j].x << "," << pixelsForCl[cl][j].y << "|";
456  if (j == 0) {
457  output.emplace_back(newpix, pixelsForCl[cl][j].adc);
458  } else {
459  output.back().add(newpix, pixelsForCl[cl][j].adc);
460  }
461  }
462  if (verbose)
463  std::cout << std::endl;
464  if (!pixelsForCl[cl].empty()) {
465  if (forceXError_ > 0)
466  output.back().setSplitClusterErrorX(forceXError_);
467  if (forceYError_ > 0)
468  output.back().setSplitClusterErrorY(forceYError_);
469  }
470  }
471  // if(verbose) std::cout << "Weights" << std::endl;
472  // if(verbose) print(theWeights,aCluster,1);
473  // if(verbose) std::cout << "Unused charge" << std::endl;
474  // if(verbose) print(theBufferResidual,aCluster);
475 
476  return output;
477 }
std::multimap< float, int > secondDistScore(const std::vector< std::vector< float >> &distanceMap)
U second(std::pair< T, U > const &p)
const std::vector< Pixel > pixels() const
T sqrt(T t)
Definition: SSEVec.h:19
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double f[11][100]
uint16_t *__restrict__ uint16_t const *__restrict__ adc

◆ produce()

void JetCoreClusterSplitter::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
override

Definition at line 92 of file JetCoreClusterSplitter.cc.

References funct::abs(), edmNew::DetSetVector< T >::begin(), HltBtagPostValidation_cff::c, centralMIPCharge_, SiPixelCluster::charge(), chargeFracMin_, GetRecoTauVFromDQM_MC_cff::cl2, HLT_FULL_cff::cores, cores_, gather_cfg::cout, PbPb_ZMuSkimMuonDPG_cff::deltaR, deltaR_, edmNew::DetSetVector< T >::end(), f, trigObjTnPSource_cfi::filler, dqmdumpme::first, edm::EventSetup::getData(), edmNew::DetSetVector< T >::id(), iEvent, metsig::jet, SiPixelCluster::minPixelRow(), eostools::move(), convertSQLitetoXML_cfg::output, pixelClusters_, createTree::pp, ptMin_, SiPixelCluster::sizeX(), SiPixelCluster::sizeY(), split(), mathSSE::sqrt(), GeomDet::surface(), tCPE_, Surface::toGlobal(), tTrackingGeom_, verbose, AlignmentTracksFromVertexSelector_cfi::vertices, vertices_, and PV3DBase< T, PVType, FrameType >::z().

92  {
93  using namespace edm;
94  const auto& geometry = &iSetup.getData(tTrackingGeom_);
95 
97  iEvent.getByToken(pixelClusters_, inputPixelClusters);
98 
100  iEvent.getByToken(vertices_, vertices);
101  const reco::Vertex& pv = (*vertices)[0];
102 
104  iEvent.getByToken(cores_, cores);
105 
107  auto output = std::make_unique<edmNew::DetSetVector<SiPixelCluster>>();
108 
109  edmNew::DetSetVector<SiPixelCluster>::const_iterator detIt = inputPixelClusters->begin();
110  for (; detIt != inputPixelClusters->end(); detIt++) {
112  const edmNew::DetSet<SiPixelCluster>& detset = *detIt;
113  const GeomDet* det = geometry->idToDet(detset.id());
114  for (auto cluster = detset.begin(); cluster != detset.end(); cluster++) {
115  const SiPixelCluster& aCluster = *cluster;
116  bool hasBeenSplit = false;
117  bool shouldBeSplit = false;
118  GlobalPoint cPos =
119  det->surface().toGlobal(pp->localParametersV(aCluster, (*geometry->idToDetUnit(detIt->id())))[0].first);
120  GlobalPoint ppv(pv.position().x(), pv.position().y(), pv.position().z());
121  GlobalVector clusterDir = cPos - ppv;
122  for (unsigned int ji = 0; ji < cores->size(); ji++) {
123  if ((*cores)[ji].pt() > ptMin_) {
124  const reco::Candidate& jet = (*cores)[ji];
125  GlobalVector jetDir(jet.px(), jet.py(), jet.pz());
126  if (Geom::deltaR(jetDir, clusterDir) < deltaR_) {
127  // check if the cluster has to be splitted
128 
129  bool isEndCap = (std::abs(cPos.z()) > 30.f); // FIXME: check detID instead!
130  float jetZOverRho = jet.momentum().Z() / jet.momentum().Rho();
131  if (isEndCap)
132  jetZOverRho = jet.momentum().Rho() / jet.momentum().Z();
133  float expSizeY = std::sqrt((1.3f * 1.3f) + (1.9f * 1.9f) * jetZOverRho * jetZOverRho);
134  if (expSizeY < 1.f)
135  expSizeY = 1.f;
136  float expSizeX = 1.5f;
137  if (isEndCap) {
138  expSizeX = expSizeY;
139  expSizeY = 1.5f;
140  } // in endcap col/rows are switched
141  float expCharge = std::sqrt(1.08f + jetZOverRho * jetZOverRho) * centralMIPCharge_;
142 
143  if (aCluster.charge() > expCharge * chargeFracMin_ &&
144  (aCluster.sizeX() > expSizeX + 1 || aCluster.sizeY() > expSizeY + 1)) {
145  shouldBeSplit = true;
146  if (verbose)
147  std::cout << "Trying to split: charge and deltaR " << aCluster.charge() << " "
148  << Geom::deltaR(jetDir, clusterDir) << " size x y " << aCluster.sizeX() << " "
149  << aCluster.sizeY() << " exp. size (x,y) " << expSizeX << " " << expSizeY << " detid "
150  << detIt->id() << std::endl;
151  if (verbose)
152  std::cout << "jetZOverRho=" << jetZOverRho << std::endl;
153 
154  if (split(aCluster, filler, expCharge, expSizeY, expSizeX, jetZOverRho)) {
155  hasBeenSplit = true;
156  }
157  }
158  }
159  }
160  }
161  if (!hasBeenSplit) {
162  SiPixelCluster c = aCluster;
163  if (shouldBeSplit) {
164  // blowup the error if we failed to split a splittable cluster (does
165  // it ever happen)
166  c.setSplitClusterErrorX(c.sizeX() * (100.f / 3.f)); // this is not really blowing up .. TODO: tune
167  c.setSplitClusterErrorY(c.sizeY() * (150.f / 3.f));
168  }
169  filler.push_back(c);
170  std::push_heap(filler.begin(), filler.end(), [](SiPixelCluster const& cl1, SiPixelCluster const& cl2) {
171  return cl1.minPixelRow() < cl2.minPixelRow();
172  });
173  }
174  } // loop over clusters
175  std::sort_heap(filler.begin(), filler.end(), [](SiPixelCluster const& cl1, SiPixelCluster const& cl2) {
176  return cl1.minPixelRow() < cl2.minPixelRow();
177  });
178 
179  } // loop over det
180  iEvent.put(std::move(output));
181 }
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
edm::EDGetTokenT< edmNew::DetSetVector< SiPixelCluster > > pixelClusters_
T z() const
Definition: PV3DBase.h:61
bool split(const SiPixelCluster &aCluster, edmNew::DetSetVector< SiPixelCluster >::FastFiller &filler, float expectedADC, int sizeY, int sizeX, float jetZOverRho)
int sizeY() const
edm::EDGetTokenT< reco::VertexCollection > vertices_
edm::EDGetTokenT< edm::View< reco::Candidate > > cores_
const_iterator end(bool update=false) const
int minPixelRow() const
int charge() const
int iEvent
Definition: GenABIO.cc:224
edm::ESGetToken< PixelClusterParameterEstimator, TkPixelCPERecord > const tCPE_
T sqrt(T t)
Definition: SSEVec.h:19
int sizeX() const
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
edm::ESGetToken< GlobalTrackingGeometry, GlobalTrackingGeometryRecord > const tTrackingGeom_
id_type id(size_t cell) const
double f[11][100]
GlobalPoint toGlobal(const Point2DBase< Scalar, LocalTag > lp) const
Definition: Surface.h:79
boost::transform_iterator< IterHelp, const_IdIter > const_iterator
const_iterator begin(bool update=false) const
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:37
Pixel cluster – collection of neighboring pixels above threshold.
HLT enums.
def move(src, dest)
Definition: eostools.py:511

◆ secondDistDiffScore()

std::multimap< float, int > JetCoreClusterSplitter::secondDistDiffScore ( const std::vector< std::vector< float >> &  distanceMap)
private

Definition at line 222 of file JetCoreClusterSplitter.cc.

References closestClusters(), ztail::d, and dqmiolumiharvest::j.

223  {
224  std::multimap<float, int> scores;
225  for (unsigned int j = 0; j < distanceMap.size(); j++) {
226  std::pair<float, float> d = closestClusters(distanceMap[j]);
227  scores.insert(std::pair<float, int>(d.second - d.first, j));
228  }
229  return scores;
230 }
d
Definition: ztail.py:151
std::pair< float, float > closestClusters(const std::vector< float > &distanceMap)

◆ secondDistScore()

std::multimap< float, int > JetCoreClusterSplitter::secondDistScore ( const std::vector< std::vector< float >> &  distanceMap)
private

Definition at line 232 of file JetCoreClusterSplitter.cc.

References closestClusters(), ztail::d, and dqmiolumiharvest::j.

Referenced by fittingSplit().

232  {
233  std::multimap<float, int> scores;
234  for (unsigned int j = 0; j < distanceMap.size(); j++) {
235  std::pair<float, float> d = closestClusters(distanceMap[j]);
236  scores.insert(std::pair<float, int>(-d.second, j));
237  }
238  return scores;
239 }
d
Definition: ztail.py:151
std::pair< float, float > closestClusters(const std::vector< float > &distanceMap)

◆ split()

bool JetCoreClusterSplitter::split ( const SiPixelCluster aCluster,
edmNew::DetSetVector< SiPixelCluster >::FastFiller &  filler,
float  expectedADC,
int  sizeY,
int  sizeX,
float  jetZOverRho 
)
private

Definition at line 183 of file JetCoreClusterSplitter.cc.

References SiPixelCluster::charge(), GetRecoTauVFromDQM_MC_cff::cl2, trigObjTnPSource_cfi::filler, fittingSplit(), mps_fire::i, SiPixelCluster::minPixelRow(), RecoTauValidation_cfi::sizeX, and RecoTauValidation_cfi::sizeY.

Referenced by produce().

188  {
189  // This function should test several configuration of splitting, and then
190  // return the one with best chi2
191 
192  std::vector<SiPixelCluster> sp = fittingSplit(
193  aCluster, expectedADC, sizeY, sizeX, jetZOverRho, std::floor(aCluster.charge() / expectedADC + 0.5f));
194 
195  // for the config with best chi2
196  for (unsigned int i = 0; i < sp.size(); i++) {
197  filler.push_back(sp[i]);
198  std::push_heap(filler.begin(), filler.end(), [](SiPixelCluster const& cl1, SiPixelCluster const& cl2) {
199  return cl1.minPixelRow() < cl2.minPixelRow();
200  });
201  }
202 
203  return (!sp.empty());
204 }
int minPixelRow() const
int charge() const
std::vector< SiPixelCluster > fittingSplit(const SiPixelCluster &aCluster, float expectedADC, int sizeY, int sizeX, float jetZOverRho, unsigned int nSplitted)
Pixel cluster – collection of neighboring pixels above threshold.

Member Data Documentation

◆ centralMIPCharge_

double JetCoreClusterSplitter::centralMIPCharge_
private

Definition at line 64 of file JetCoreClusterSplitter.cc.

Referenced by fittingSplit(), and produce().

◆ chargeFracMin_

double JetCoreClusterSplitter::chargeFracMin_
private

Definition at line 56 of file JetCoreClusterSplitter.cc.

Referenced by produce().

◆ chargePerUnit_

double JetCoreClusterSplitter::chargePerUnit_
private

Definition at line 63 of file JetCoreClusterSplitter.cc.

Referenced by fittingSplit().

◆ cores_

edm::EDGetTokenT<edm::View<reco::Candidate> > JetCoreClusterSplitter::cores_
private

Definition at line 59 of file JetCoreClusterSplitter.cc.

Referenced by produce().

◆ deltaR_

double JetCoreClusterSplitter::deltaR_
private

Definition at line 55 of file JetCoreClusterSplitter.cc.

Referenced by produce().

◆ forceXError_

double JetCoreClusterSplitter::forceXError_
private

Definition at line 60 of file JetCoreClusterSplitter.cc.

Referenced by fittingSplit().

◆ forceYError_

double JetCoreClusterSplitter::forceYError_
private

Definition at line 61 of file JetCoreClusterSplitter.cc.

Referenced by fittingSplit().

◆ fractionalWidth_

double JetCoreClusterSplitter::fractionalWidth_
private

Definition at line 62 of file JetCoreClusterSplitter.cc.

Referenced by fittingSplit().

◆ pixelClusters_

edm::EDGetTokenT<edmNew::DetSetVector<SiPixelCluster> > JetCoreClusterSplitter::pixelClusters_
private

Definition at line 57 of file JetCoreClusterSplitter.cc.

Referenced by produce().

◆ ptMin_

double JetCoreClusterSplitter::ptMin_
private

Definition at line 54 of file JetCoreClusterSplitter.cc.

Referenced by produce().

◆ tCPE_

edm::ESGetToken<PixelClusterParameterEstimator, TkPixelCPERecord> const JetCoreClusterSplitter::tCPE_
private

Definition at line 51 of file JetCoreClusterSplitter.cc.

Referenced by produce().

◆ tTrackingGeom_

edm::ESGetToken<GlobalTrackingGeometry, GlobalTrackingGeometryRecord> const JetCoreClusterSplitter::tTrackingGeom_
private

Definition at line 50 of file JetCoreClusterSplitter.cc.

Referenced by produce().

◆ verbose

bool JetCoreClusterSplitter::verbose
private

Definition at line 53 of file JetCoreClusterSplitter.cc.

Referenced by fittingSplit(), and produce().

◆ vertices_

edm::EDGetTokenT<reco::VertexCollection> JetCoreClusterSplitter::vertices_
private

Definition at line 58 of file JetCoreClusterSplitter.cc.

Referenced by produce().