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
 
bool hasAbilityToProduceInBeginLumis () const final
 
bool hasAbilityToProduceInBeginRuns () const final
 
bool hasAbilityToProduceInEndLumis () const final
 
bool hasAbilityToProduceInEndRuns () const final
 

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_
 
std::string pixelCPE_
 
double ptMin_
 
bool verbose
 
edm::EDGetTokenT< reco::VertexCollectionvertices_
 

Additional Inherited Members

- Public Types inherited from edm::stream::EDProducer<>
typedef CacheContexts< T... > CacheTypes
 
typedef CacheTypes::GlobalCache GlobalCache
 
typedef AbilityChecker< T... > HasAbility
 
typedef CacheTypes::LuminosityBlockCache LuminosityBlockCache
 
typedef LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCacheLuminosityBlockContext
 
typedef CacheTypes::LuminosityBlockSummaryCache LuminosityBlockSummaryCache
 
typedef CacheTypes::RunCache RunCache
 
typedef RunContextT< RunCache, GlobalCacheRunContext
 
typedef CacheTypes::RunSummaryCache RunSummaryCache
 

Detailed Description

Definition at line 26 of file JetCoreClusterSplitter.cc.

Constructor & Destructor Documentation

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

Definition at line 64 of file JetCoreClusterSplitter.cc.

65  : verbose(iConfig.getParameter<bool>("verbose")),
66  pixelCPE_(iConfig.getParameter<std::string>("pixelCPE")),
67  ptMin_(iConfig.getParameter<double>("ptMin")),
68  deltaR_(iConfig.getParameter<double>("deltaRmax")),
69  chargeFracMin_(iConfig.getParameter<double>("chargeFractionMin")),
71  consumes<edmNew::DetSetVector<SiPixelCluster>>(iConfig.getParameter<edm::InputTag>("pixelClusters"))),
72  vertices_(consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("vertices"))),
73  cores_(consumes<edm::View<reco::Candidate>>(iConfig.getParameter<edm::InputTag>("cores"))),
74  forceXError_(iConfig.getParameter<double>("forceXError")),
75  forceYError_(iConfig.getParameter<double>("forceYError")),
76  fractionalWidth_(iConfig.getParameter<double>("fractionalWidth")),
77  chargePerUnit_(iConfig.getParameter<double>("chargePerUnit")),
78  centralMIPCharge_(iConfig.getParameter<double>("centralMIPCharge"))
79 
80 {
81  produces<edmNew::DetSetVector<SiPixelCluster>>();
82 }
T getParameter(std::string const &) const
edm::EDGetTokenT< edmNew::DetSetVector< SiPixelCluster > > pixelClusters_
edm::EDGetTokenT< reco::VertexCollection > vertices_
edm::EDGetTokenT< edm::View< reco::Candidate > > cores_
JetCoreClusterSplitter::~JetCoreClusterSplitter ( )
override

Definition at line 84 of file JetCoreClusterSplitter.cc.

84 {}

Member Function Documentation

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

Definition at line 208 of file JetCoreClusterSplitter.cc.

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

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

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

Definition at line 242 of file JetCoreClusterSplitter.cc.

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

242  {
243  std::multimap<float, int> scores;
244  for (unsigned int j = 0; j < distanceMap.size(); j++) {
245  std::pair<float, float> d = closestClusters(distanceMap[j]);
246  scores.insert(std::pair<float, int>(-d.first, j));
247  }
248  return scores;
249 }
d
Definition: ztail.py:151
std::pair< float, float > closestClusters(const std::vector< float > &distanceMap)
std::vector< SiPixelCluster > JetCoreClusterSplitter::fittingSplit ( const SiPixelCluster aCluster,
float  expectedADC,
int  sizeY,
int  sizeX,
float  jetZOverRho,
unsigned int  nSplitted 
)
private

Definition at line 251 of file JetCoreClusterSplitter.cc.

References funct::abs(), ecalMGPA::adc(), centralMIPCharge_, chargePerUnit_, GetRecoTauVFromDQM_MC_cff::cl, gather_cfg::cout, DEFINE_FWK_MODULE, 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(), x, and y.

Referenced by split().

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

Definition at line 88 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, geometry, edm::EventSetup::get(), edm::Event::getByToken(), edmNew::DetSetVector< T >::id(), GlobalTrackingGeometry::idToDet(), GlobalTrackingGeometry::idToDetUnit(), metsig::jet, PixelClusterParameterEstimator::localParametersV(), SiPixelCluster::minPixelRow(), reco::Candidate::momentum(), eostools::move(), convertSQLitetoXML_cfg::output, pixelClusters_, pixelCPE_, reco::Vertex::position(), createTree::pp, edm::ESHandle< T >::product(), ptMin_, edm::Event::put(), MetAnalyzer::pv(), reco::Candidate::px(), reco::Candidate::py(), reco::Candidate::pz(), SiPixelCluster::setSplitClusterErrorX(), SiPixelCluster::setSplitClusterErrorY(), SiPixelCluster::sizeX(), SiPixelCluster::sizeY(), split(), mathSSE::sqrt(), GeomDet::surface(), Surface::toGlobal(), pwdgSkimBPark_cfi::vertices, vertices_, and PV3DBase< T, PVType, FrameType >::z().

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

Definition at line 223 of file JetCoreClusterSplitter.cc.

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

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

Definition at line 233 of file JetCoreClusterSplitter.cc.

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

Referenced by fittingSplit().

233  {
234  std::multimap<float, int> scores;
235  for (unsigned int j = 0; j < distanceMap.size(); j++) {
236  std::pair<float, float> d = closestClusters(distanceMap[j]);
237  scores.insert(std::pair<float, int>(-d.second, j));
238  }
239  return scores;
240 }
d
Definition: ztail.py:151
std::pair< float, float > closestClusters(const std::vector< float > &distanceMap)
bool JetCoreClusterSplitter::split ( const SiPixelCluster aCluster,
edmNew::DetSetVector< SiPixelCluster >::FastFiller &  filler,
float  expectedADC,
int  sizeY,
int  sizeX,
float  jetZOverRho 
)
private

Definition at line 184 of file JetCoreClusterSplitter.cc.

References edmNew::DetSetVector< T >::FastFiller::begin(), SiPixelCluster::charge(), GetRecoTauVFromDQM_MC_cff::cl2, edmNew::DetSetVector< T >::FastFiller::end(), fittingSplit(), mps_fire::i, SiPixelCluster::minPixelRow(), and edmNew::DetSetVector< T >::FastFiller::push_back().

Referenced by produce().

189  {
190  // This function should test several configuration of splitting, and then
191  // return the one with best chi2
192 
193  std::vector<SiPixelCluster> sp = fittingSplit(
194  aCluster, expectedADC, sizeY, sizeX, jetZOverRho, std::floor(aCluster.charge() / expectedADC + 0.5f));
195 
196  // for the config with best chi2
197  for (unsigned int i = 0; i < sp.size(); i++) {
198  filler.push_back(sp[i]);
199  std::push_heap(filler.begin(), filler.end(), [](SiPixelCluster const& cl1, SiPixelCluster const& cl2) {
200  return cl1.minPixelRow() < cl2.minPixelRow();
201  });
202  }
203 
204  return (!sp.empty());
205 }
void push_back(data_type const &d)
int charge() const
int minPixelRow() 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

double JetCoreClusterSplitter::centralMIPCharge_
private

Definition at line 61 of file JetCoreClusterSplitter.cc.

Referenced by fittingSplit(), and produce().

double JetCoreClusterSplitter::chargeFracMin_
private

Definition at line 53 of file JetCoreClusterSplitter.cc.

Referenced by produce().

double JetCoreClusterSplitter::chargePerUnit_
private

Definition at line 60 of file JetCoreClusterSplitter.cc.

Referenced by fittingSplit().

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

Definition at line 56 of file JetCoreClusterSplitter.cc.

Referenced by produce().

double JetCoreClusterSplitter::deltaR_
private

Definition at line 52 of file JetCoreClusterSplitter.cc.

Referenced by produce().

double JetCoreClusterSplitter::forceXError_
private

Definition at line 57 of file JetCoreClusterSplitter.cc.

Referenced by fittingSplit().

double JetCoreClusterSplitter::forceYError_
private

Definition at line 58 of file JetCoreClusterSplitter.cc.

Referenced by fittingSplit().

double JetCoreClusterSplitter::fractionalWidth_
private

Definition at line 59 of file JetCoreClusterSplitter.cc.

Referenced by fittingSplit().

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

Definition at line 54 of file JetCoreClusterSplitter.cc.

Referenced by produce().

std::string JetCoreClusterSplitter::pixelCPE_
private

Definition at line 50 of file JetCoreClusterSplitter.cc.

Referenced by produce().

double JetCoreClusterSplitter::ptMin_
private

Definition at line 51 of file JetCoreClusterSplitter.cc.

Referenced by produce().

bool JetCoreClusterSplitter::verbose
private

Definition at line 49 of file JetCoreClusterSplitter.cc.

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

Definition at line 55 of file JetCoreClusterSplitter.cc.

Referenced by produce().