CMS 3D CMS Logo

List of all members | Public Types | Public Member Functions | Private Member Functions | Private Attributes | Static Private Attributes
CSCHitFromStripOnly Class Reference

#include <CSCHitFromStripOnly.h>

Public Types

typedef std::array< CSCStripData, 100 > PulseHeightMap
 

Public Member Functions

 CSCHitFromStripOnly (const edm::ParameterSet &ps)
 
bool ganged ()
 
std::vector< CSCStripHitrunStrip (const CSCDetId &id, const CSCLayer *layer, const CSCStripDigiCollection::Range &rstripd)
 
void setConditions (const CSCRecoConditions *reco)
 
void setGanged (bool ig)
 
 ~CSCHitFromStripOnly ()
 

Private Member Functions

void fillPulseHeights (const CSCStripDigiCollection::Range &rstripd)
 Store SCA pulseheight information from strips in digis of one layer. More...
 
float findHitOnStripPosition (const std::vector< CSCStripHitData > &data, const int &centerStrip)
 Find position of hit in strip cluster in terms of strip #. More...
 
void findMaxima (const CSCDetId &id)
 Find local maxima. More...
 
bool isDeadStrip (const CSCDetId &id, int centralStrip, int nstrips)
 Is the strip 'bad'? More...
 
bool isNearDeadStrip (const CSCDetId &id, int centralStrip, int nstrips)
 Is either neighbour 'bad'? More...
 
bool isPeakOK (int iStrip, float heightCluster)
 
float makeCluster (int centerStrip)
 Make clusters using local maxima. More...
 
CSCStripHitData makeStripData (int centerStrip, int offset)
 

Private Attributes

CSCPedestalChoicecalcped_
 
int clusterSize
 
float gainWeight [80]
 
bool ganged_
 
CSCDetId id_
 
const CSCLayerlayer_
 
unsigned nstrips_
 
const CSCRecoConditionsrecoConditions_
 
std::vector< float > strips_adc
 
std::vector< float > strips_adcRaw
 
std::vector< int > theClosestMaximum
 
std::vector< int > theConsecutiveStrips
 
std::vector< int > theMaxima
 
PulseHeightMap thePulseHeightMap
 
std::vector< int > theStrips
 
float theThresholdForAPeak
 
float theThresholdForCluster
 
int tmax_cluster
 
bool useCalib
 

Static Private Attributes

static const int theClusterSize = 3
 

Detailed Description

below)

Search for strip with ADC output exceeding theThresholdForAPeak. For each of these strips, build a cluster of strip of size theClusterSize (typically 5 strips). Finally, make a Strip Hit out of these clusters by finding the center-of-mass position of the hit The DetId, strip hit position, and peaking time are stored in a CSCStripHit collection.

Be careful with the ME_1/a chambers: the 48 strips are ganged into 16 channels, Each channel contains signals from three strips, each separated by 16 strips from the next.

Definition at line 34 of file CSCHitFromStripOnly.h.

Member Typedef Documentation

◆ PulseHeightMap

Definition at line 36 of file CSCHitFromStripOnly.h.

Constructor & Destructor Documentation

◆ CSCHitFromStripOnly()

CSCHitFromStripOnly::CSCHitFromStripOnly ( const edm::ParameterSet ps)
explicit

Definition at line 22 of file CSCHitFromStripOnly.cc.

References calcped_, edm::ParameterSet::getParameter(), LogTrace, theThresholdForAPeak, theThresholdForCluster, and useCalib.

23  : recoConditions_(nullptr), calcped_(nullptr), ganged_(false) {
24  useCalib = ps.getParameter<bool>("CSCUseCalibrations");
25  bool useStaticPedestals = ps.getParameter<bool>("CSCUseStaticPedestals");
26  int noOfTimeBinsForDynamicPed = ps.getParameter<int>("CSCNoOfTimeBinsForDynamicPedestal");
27 
28  theThresholdForAPeak = ps.getParameter<double>("CSCStripPeakThreshold");
29  theThresholdForCluster = ps.getParameter<double>("CSCStripClusterChargeCut");
30 
31  LogTrace("CSCRecHit") << "[CSCHitFromStripOnly] CSCUseStaticPedestals = " << useStaticPedestals;
32  if (!useStaticPedestals)
33  LogTrace("CSCRecHit") << "[CSCHitFromStripOnly] CSCNoOfTimeBinsForDynamicPedestal = " << noOfTimeBinsForDynamicPed;
34 
35  if (useStaticPedestals) {
37  } else {
38  if (noOfTimeBinsForDynamicPed == 1) {
40  } else {
41  calcped_ = new CSCDynamicPedestal2(); // NORMAL DEFAULT!
42  }
43  }
44 }
CSCPedestalChoice * calcped_
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
#define LogTrace(id)
const CSCRecoConditions * recoConditions_

◆ ~CSCHitFromStripOnly()

CSCHitFromStripOnly::~CSCHitFromStripOnly ( )

Definition at line 46 of file CSCHitFromStripOnly.cc.

References calcped_.

46 { delete calcped_; }
CSCPedestalChoice * calcped_

Member Function Documentation

◆ fillPulseHeights()

void CSCHitFromStripOnly::fillPulseHeights ( const CSCStripDigiCollection::Range rstripd)
private

Store SCA pulseheight information from strips in digis of one layer.

Definition at line 350 of file CSCHitFromStripOnly.cc.

References cms::cuda::assert(), calcped_, filterCSVwithJSON::copy, gainWeight, ganged(), id_, dqmiolumiharvest::j, CSCPedestalChoice::pedestal(), recoConditions_, thePulseHeightMap, tmax, and useCalib.

Referenced by runStrip().

350  {
351  // Loop over strip digis in one CSCLayer and fill PulseHeightMap with pedestal-subtracted
352  // SCA pulse heights.
353 
354  for (auto& ph : thePulseHeightMap)
355  ph.reset();
356 
357  // for storing sca pulseheights once they may no longer be integer (e.g. after ped subtraction)
358  for (CSCStripDigiCollection::const_iterator it = rstripd.first; it != rstripd.second; ++it) {
359  int thisChannel = (*it).getStrip();
360  auto& stripData = thePulseHeightMap[thisChannel - 1];
361  auto& scaRaw = stripData.phRaw_;
362  auto& sca = stripData.ph_;
363 
364  auto const& scaOri = (*it).getADCCounts();
365  assert(scaOri.size() == 8);
366  // Fill sca from scaRaw, implicitly converting to float
367  std::copy(scaOri.begin(), scaOri.end(), scaRaw.begin());
368  std::copy(scaRaw.begin(), scaRaw.end(), sca.begin());
369 
370  //@@ Find bin with largest pulseheight (_before_ ped subtraction - shouldn't matter, right?)
371  int tmax = std::max_element(sca.begin(), sca.end()) - sca.begin(); // counts from 0
372 
373  // get pedestal - calculated as appropriate - for this sca pulse
374  float ped = calcped_->pedestal(sca, recoConditions_, id_, thisChannel);
375 
376  // subtract the pedestal (from BOTH sets of sca pulseheights)
377  std::for_each(sca.begin(), sca.end(), CSCSubtractPedestal(ped));
378  std::for_each(scaRaw.begin(), scaRaw.end(), CSCSubtractPedestal(ped));
379 
380  //@@ Max in first 3 or last time bins is unacceptable, if so set to zero (why?)
381  float phmax = 0.f;
382  if (tmax > 2 && tmax < 7) {
383  phmax = sca[tmax];
384  }
385  stripData.phmax_ = phmax;
386  stripData.tmax_ = tmax;
387 
388  // Fill the map, possibly apply gains from cond data, and unfold ME1A channels
389  // (To apply gains use CSCStripData::op*= which scales only the non-raw sca ph's;
390  // but note that both sca & scaRaw are pedestal-subtracted.)
391 
392  // From StripDigi, thisChannel labels strip channel. Values phmax, tmax, scaRaw, sca belong to thisChannel
393  if (useCalib)
394  stripData *= gainWeight[thisChannel - 1];
395 
396  // for ganged ME1a need to duplicate values on istrip=thisChannel to iStrip+16 and iStrip+32
397  if (ganged()) {
398  for (int j = 1; j < 3; ++j) {
399  thePulseHeightMap[thisChannel - 1 + 16 * j] = stripData;
400  }
401  }
402  }
403 }
CSCPedestalChoice * calcped_
assert(be >=bs)
PulseHeightMap thePulseHeightMap
const CSCRecoConditions * recoConditions_
static const double tmax[3]
std::vector< DigiType >::const_iterator const_iterator
virtual float pedestal(const std::vector< float > &sca, const CSCRecoConditions *cond=nullptr, const CSCDetId id=0, int ichan=0)=0

◆ findHitOnStripPosition()

float CSCHitFromStripOnly::findHitOnStripPosition ( const std::vector< CSCStripHitData > &  data,
const int &  centerStrip 
)
private

Find position of hit in strip cluster in terms of strip #.

Definition at line 546 of file CSCHitFromStripOnly.cc.

References filterCSVwithJSON::copy, data, mps_fire::i, LogTrace, digitizers_cfi::strip, strips_adc, strips_adcRaw, and w().

Referenced by makeCluster().

546  {
547  float strippos = -1.;
548 
549  if (data.empty())
550  return strippos;
551 
552  // biggestStrip is strip with largest pulse height
553  // Use pointer subtraction
554 
555  int biggestStrip = max_element(data.begin(), data.end()) - data.begin();
556  strippos = data[biggestStrip].strip() * 1.;
557 
558  // If more than one strip: use centroid to find center of cluster
559  // but only use time bin == tmax (otherwise, bias centroid).
560  float sum = 0.;
561  float sum_w = 0.;
562 
563  // std::vector<float> w(4);
564  // std::vector<float> wRaw(4);
565 
566  for (size_t i = 0; i != data.size(); ++i) {
567  auto const& w = data[i].ph();
568  auto const& wRaw = data[i].phRaw();
569 
570  // (Require ADC to be > 0.)
571  // No later studies suggest that this only do harm
572  /*
573  for ( size_t j = 0; j != w.size(); ++j ) {
574  if ( w[j] < 0. ) w[j] = 0.001;
575  }
576  */
577 
578  // Fill the data members
579  std::copy(w.begin(), w.end(), std::back_inserter(strips_adc));
580  std::copy(wRaw.begin(), wRaw.end(), std::back_inserter(strips_adcRaw));
581 
582  if (data[i].strip() < 1) {
583  LogTrace("CSCRecHit") << "[CSCHitFromStripOnly::findHitOnStripPosition] problem in indexing of strip, strip= "
584  << data[i].strip();
585  }
586  sum_w += w[1];
587  sum += w[1] * data[i].strip();
588  }
589 
590  if (sum_w > 0.)
591  strippos = sum / sum_w;
592 
593  return strippos;
594 }
T w() const
#define LogTrace(id)
std::vector< float > strips_adcRaw
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:79
std::vector< float > strips_adc

◆ findMaxima()

void CSCHitFromStripOnly::findMaxima ( const CSCDetId id)
private

Find local maxima.

Definition at line 411 of file CSCHitFromStripOnly.cc.

References mps_fire::i, isDeadStrip(), isPeakOK(), dqmiolumiharvest::j, dqmdumpme::k, cmsLHEtoEOSManager::l, nstrips_, theClosestMaximum, theConsecutiveStrips, theMaxima, and thePulseHeightMap.

Referenced by runStrip().

411  {
412  theMaxima.clear();
413  theConsecutiveStrips.clear();
414  theClosestMaximum.clear();
415  for (size_t i = 0; i != thePulseHeightMap.size(); ++i) {
416  // sum 3 strips so that hits between strips are not suppressed
417  float heightCluster = 0.;
418 
419  bool maximumFound = false;
420  // Left edge of chamber
421  if (!isDeadStrip(id, i + 1, nstrips_)) { // Is it i or i+1
422  if (i == 0) {
423  heightCluster = thePulseHeightMap[i].phmax() + thePulseHeightMap[i + 1].phmax();
424  // Have found a strip Hit if...
425  if (thePulseHeightMap[i].phmax() >= thePulseHeightMap[i + 1].phmax() && isPeakOK(i, heightCluster)) {
426  maximumFound = true;
427  }
428  // Right edge of chamber
429  } else if (i == thePulseHeightMap.size() - 1) {
430  heightCluster = thePulseHeightMap[i - 1].phmax() + thePulseHeightMap[i].phmax();
431  // Have found a strip Hit if...
432  if (thePulseHeightMap[i].phmax() > thePulseHeightMap[i - 1].phmax() && isPeakOK(i, heightCluster)) {
433  maximumFound = true;
434  }
435  // Any other strips
436  } else {
437  heightCluster =
438  thePulseHeightMap[i - 1].phmax() + thePulseHeightMap[i].phmax() + thePulseHeightMap[i + 1].phmax();
439  // Have found a strip Hit if...
440  if (thePulseHeightMap[i].phmax() > thePulseHeightMap[i - 1].phmax() &&
441  thePulseHeightMap[i].phmax() >= thePulseHeightMap[i + 1].phmax() && isPeakOK(i, heightCluster)) {
442  maximumFound = true;
443  }
444  }
445  }
446  //---- Consecutive strips with charge (real cluster); if too wide - measurement is not accurate
447  if (maximumFound) {
448  int numberOfConsecutiveStrips = 1;
449  float testThreshold = 10.; //---- ADC counts;
450  //---- this is not XTalk corrected so it is correct in first approximation only
451  int j = 0;
452  for (int l = 0; l < 8; ++l) {
453  if (j < 0)
454  edm::LogWarning("FailedStripCountingWrongConsecutiveStripNumber")
455  << "This should never occur!!! Contact CSC expert!";
456  ++j;
457  bool signalPresent = false;
458  for (int k = 0; k < 2; ++k) {
459  j *= -1; //---- check from left and right
460  int anotherConsecutiveStrip = i + j;
461  if (anotherConsecutiveStrip >= 0 && anotherConsecutiveStrip < int(thePulseHeightMap.size())) {
462  if (thePulseHeightMap[anotherConsecutiveStrip].phmax() > testThreshold) {
463  ++numberOfConsecutiveStrips;
464  signalPresent = true;
465  }
466  }
467  }
468  if (!signalPresent) {
469  break;
470  }
471  }
472 
473  bool additional_maxima_found = false;
474  // search for additional maxima if:
475  // - hit is closer than 3 strips from the edge
476  // - enough consecutive strips with signal
477  // - strip charge distribution looks abnormal
478 
479  if (i > 2 && i + 3 < thePulseHeightMap.size() && numberOfConsecutiveStrips > 3) {
480  //try to look for additional maxima at the left side from the main maxima
481 
482  if (((thePulseHeightMap[i + 1].phmax() >= thePulseHeightMap[i - 1].phmax() &&
483  thePulseHeightMap[i + 1].phmax() >= thePulseHeightMap[i - 2].phmax() &&
484  thePulseHeightMap[i + 2].phmax() <= thePulseHeightMap[i - 2].phmax()) ||
485  (thePulseHeightMap[i + 1].phmax() <= thePulseHeightMap[i - 1].phmax() &&
486  thePulseHeightMap[i + 1].phmax() <= thePulseHeightMap[i - 2].phmax())) &&
487  //to avoid close maxima delimitation (this is already present in the code)
488  thePulseHeightMap[i - 1].phmax() >= thePulseHeightMap[i - 2].phmax() &&
489  //no need in a small charge maxima (might need adjustment)
490  thePulseHeightMap[i - 2].phmax() > 20) {
491  additional_maxima_found = true;
492  theMaxima.push_back(i - 2); //insert left maxima first
493  //insert the same number of cosecutive strips, because they belong to both maximas
494  theConsecutiveStrips.push_back(numberOfConsecutiveStrips);
495  theMaxima.push_back(i); //insert main maxima
496  //insert the same number of cosecutive strips, because they belong to both maximas
497  theConsecutiveStrips.push_back(numberOfConsecutiveStrips);
498 
499  } //looking for additional maxima on the left
500 
501  //try to look for additional maxima at the right side from the main maxima
502  if (((thePulseHeightMap[i + 1].phmax() >= thePulseHeightMap[i - 1].phmax() &&
503  thePulseHeightMap[i + 2].phmax() >= thePulseHeightMap[i - 1].phmax()) ||
504  (thePulseHeightMap[i + 1].phmax() <= thePulseHeightMap[i - 1].phmax() &&
505  thePulseHeightMap[i + 2].phmax() <= thePulseHeightMap[i - 1].phmax() &&
506  thePulseHeightMap[i + 2].phmax() >= thePulseHeightMap[i - 2].phmax())) &&
507  //to avoid close maxima delimitation (this is already present in the code)
508  thePulseHeightMap[i + 1].phmax() >= thePulseHeightMap[i + 2].phmax() &&
509  //no need in a small charge maxima (might need adjustment)
510  thePulseHeightMap[i + 2].phmax() > 20) {
511  additional_maxima_found = true;
512  theMaxima.push_back(i); //insert main maxima first
513  //insert the same number of cosecutive strips, because they belong to both maximas
514  theConsecutiveStrips.push_back(numberOfConsecutiveStrips);
515  theMaxima.push_back(i + 2); //insert right maxima
516  //insert the same number of cosecutive strips, because they belong to both maximas
517  theConsecutiveStrips.push_back(numberOfConsecutiveStrips);
518 
519  } //looking for additional maxima on the right
520 
521  //if nothing additional found fill the maxima
522  if (!additional_maxima_found) {
523  theMaxima.push_back(i);
524  theConsecutiveStrips.push_back(numberOfConsecutiveStrips);
525  }
526  } else { //not the case for looking for the additional maxima
527  theMaxima.push_back(i);
528  theConsecutiveStrips.push_back(numberOfConsecutiveStrips);
529  }
530  } //if maximafound
531  } //all pulses
532 } //find maxima procedure
std::vector< int > theConsecutiveStrips
bool isDeadStrip(const CSCDetId &id, int centralStrip, int nstrips)
Is the strip &#39;bad&#39;?
std::vector< int > theClosestMaximum
PulseHeightMap thePulseHeightMap
std::vector< int > theMaxima
bool isPeakOK(int iStrip, float heightCluster)
Log< level::Warning, false > LogWarning

◆ ganged()

bool CSCHitFromStripOnly::ganged ( )
inline

Definition at line 48 of file CSCHitFromStripOnly.h.

References ganged_.

Referenced by fillPulseHeights(), and runStrip().

48 { return ganged_; }

◆ isDeadStrip()

bool CSCHitFromStripOnly::isDeadStrip ( const CSCDetId id,
int  centralStrip,
int  nstrips 
)
private

Is the strip 'bad'?

Definition at line 602 of file CSCHitFromStripOnly.cc.

References CSCRecoConditions::badStrip(), and recoConditions_.

Referenced by findMaxima(), and runStrip().

602  {
603  return recoConditions_->badStrip(id, centralStrip, nstrips);
604 }
const CSCRecoConditions * recoConditions_
bool badStrip(const CSCDetId &id, int geomStrip, int nstrips) const
Is the strip bad?

◆ isNearDeadStrip()

bool CSCHitFromStripOnly::isNearDeadStrip ( const CSCDetId id,
int  centralStrip,
int  nstrips 
)
private

Is either neighbour 'bad'?

Definition at line 596 of file CSCHitFromStripOnly.cc.

References CSCRecoConditions::nearBadStrip(), and recoConditions_.

596  {
597  //@@ Tim says: not sure I understand this properly... but just moved code to CSCRecoConditions
598  // where it can handle the conversion from strip to channel etc.
599  return recoConditions_->nearBadStrip(id, centralStrip, nstrips);
600 }
bool nearBadStrip(const CSCDetId &id, int geomStrip, int nstrips) const
Is a neighbour bad?
const CSCRecoConditions * recoConditions_

◆ isPeakOK()

bool CSCHitFromStripOnly::isPeakOK ( int  iStrip,
float  heightCluster 
)
private

Definition at line 534 of file CSCHitFromStripOnly.cc.

References mps_fire::i, thePulseHeightMap, theThresholdForAPeak, and theThresholdForCluster.

Referenced by findMaxima().

534  {
535  int i = iStrip;
536  bool peakOK = (thePulseHeightMap[i].phmax() > theThresholdForAPeak && heightCluster > theThresholdForCluster &&
537  // ... and proper peak time; note that the values below are used elsewhere in this file;
538  // they should become parameters or at least constants defined in appropriate place
539  thePulseHeightMap[i].tmax() > 2 && thePulseHeightMap[i].tmax() < 7);
540  return peakOK;
541 }
PulseHeightMap thePulseHeightMap

◆ makeCluster()

float CSCHitFromStripOnly::makeCluster ( int  centerStrip)
private

Make clusters using local maxima.

Definition at line 211 of file CSCHitFromStripOnly.cc.

References clusterSize, data, findHitOnStripPosition(), mps_fire::i, LogTrace, makeStripData(), nstrips_, theClusterSize, and theStrips.

Referenced by runStrip().

211  {
212  float strippos = -1.;
214  std::vector<CSCStripHitData> stripDataV;
215 
216  // We only want to use strip position in terms of strip # for the strip hit. //@@ What other choice is there?
217 
218  // If the cluster size is such that you go beyond the edge of detector, shrink cluster appropriately
219  for (int i = 1; i < theClusterSize / 2 + 1; ++i) {
220  if (centerStrip - i < 1 || centerStrip + i > int(nstrips_)) {
221  // Shrink cluster size, but keep it an odd number of strips.
222  clusterSize = 2 * i - 1;
223  }
224  }
225  for (int i = -clusterSize / 2; i <= clusterSize / 2; ++i) {
226  CSCStripHitData data = makeStripData(centerStrip, i);
227  stripDataV.push_back(data);
228  theStrips.push_back(centerStrip + i);
229  }
230  strippos = findHitOnStripPosition(stripDataV, centerStrip);
231 
232  LogTrace("CSCHitFromStripOnly") << "[CSCHitFromStripOnly::makeCluster] centerStrip= " << centerStrip
233  << " strippos=" << strippos;
234 
235  return strippos;
236 }
std::vector< int > theStrips
#define LogTrace(id)
float findHitOnStripPosition(const std::vector< CSCStripHitData > &data, const int &centerStrip)
Find position of hit in strip cluster in terms of strip #.
static const int theClusterSize
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:79
CSCStripHitData makeStripData(int centerStrip, int offset)

◆ makeStripData()

CSCStripHitData CSCHitFromStripOnly::makeStripData ( int  centerStrip,
int  offset 
)
private

makeStripData

Definition at line 241 of file CSCHitFromStripOnly.cc.

References gpuClustering::adc, clusterSize, filterCSVwithJSON::copy, spr::find(), mps_fire::i, dqmdumpme::k, LogTrace, SiStripPI::min, hltrates_dqm_sourceclient-live_cfg::offset, Validation_hcalonly_cfi::sign, theMaxima, thePulseHeightMap, tmax, tmax_cluster, and validateGeometry_cfg::valid.

Referenced by makeCluster().

241  {
242  CSCStripHitData prelimData;
243  int thisStrip = centerStrip + offset;
244 
245  int tmax = thePulseHeightMap[centerStrip - 1].tmax();
246  tmax_cluster = tmax;
247 
248  std::vector<float> adc(4);
249  std::vector<float> adcRaw(4);
250 
251  // Fill adc & adcRaw
252 
253  int istart = tmax - 1;
254  int istop = std::min(tmax + 2, 7); // there are only time bins 0-7
255  adc[3] = 0.1; // in case it isn't filled
256 
257  if (tmax > 2 && tmax < 7) { // for time bins 3-6
258  int ibin = thisStrip - 1;
259  if (thePulseHeightMap[ibin].valid()) {
260  std::copy(
261  thePulseHeightMap[ibin].ph().begin() + istart, thePulseHeightMap[ibin].ph().begin() + istop + 1, adc.begin());
262 
263  std::copy(thePulseHeightMap[ibin].phRaw().begin() + istart,
264  thePulseHeightMap[ibin].phRaw().begin() + istop + 1,
265  adcRaw.begin());
266  }
267  } else {
268  adc[0] = 0.1;
269  adc[1] = 0.1;
270  adc[2] = 0.1;
271  adc[3] = 0.1;
272  adcRaw = adc;
273  LogTrace("CSCRecHit") << "[CSCHitFromStripOnly::makeStripData] Tmax out of range: contact CSC expert!";
274  }
275 
276  if (offset == 0) {
277  prelimData = CSCStripHitData(thisStrip, tmax_cluster, adcRaw, adc);
278  } else {
279  int sign = offset > 0 ? 1 : -1;
280  // If there's another maximum that would like to use part of this cluster,
281  // it gets shared in proportion to the height of the maxima
282  for (int i = 1; i <= clusterSize / 2; ++i) {
283  // Find the direction of the offset
284  int testStrip = thisStrip + sign * i;
285  std::vector<int>::iterator otherMax = find(theMaxima.begin(), theMaxima.end(), testStrip - 1);
286 
287  // No other maxima found, so just store
288  if (otherMax == theMaxima.end()) {
289  prelimData = CSCStripHitData(thisStrip, tmax_cluster, adcRaw, adc);
290  } else {
291  // Another maximum found - share
292  std::vector<float> adc1(4);
293  std::vector<float> adcRaw1(4);
294  std::vector<float> adc2(4);
295  std::vector<float> adcRaw2(4);
296  // In case we only copy (below) into 3 of the 4 bins i.e. when istart=5, istop=7
297  adc1[3] = 0.1;
298  adc2[3] = 0.1;
299  adcRaw1[3] = 0.1;
300  adcRaw2[3] = 0.1;
301 
302  // Fill adcN with content of time bins tmax-1 to tmax+2 (if it exists!)
303  if (tmax > 2 && tmax < 7) { // for time bin tmax from 3-6
304  int ibin = testStrip - 1;
305  int jbin = centerStrip - 1;
306  if (thePulseHeightMap[ibin].valid()) {
307  std::copy(thePulseHeightMap[ibin].ph().begin() + istart,
308  thePulseHeightMap[ibin].ph().begin() + istop + 1,
309  adc1.begin());
310  std::copy(thePulseHeightMap[ibin].phRaw().begin() + istart,
311  thePulseHeightMap[ibin].phRaw().begin() + istop + 1,
312  adcRaw1.begin());
313  }
314 
315  if (thePulseHeightMap[jbin].valid()) {
316  std::copy(thePulseHeightMap[jbin].ph().begin() + istart,
317  thePulseHeightMap[jbin].ph().begin() + istop + 1,
318  adc2.begin());
319 
320  std::copy(thePulseHeightMap[jbin].phRaw().begin() + istart,
321  thePulseHeightMap[jbin].phRaw().begin() + istop + 1,
322  adcRaw2.begin());
323  }
324  } else {
325  adc1.assign(4, 0.1);
326  adcRaw1 = adc1;
327  adc2.assign(4, 0.1);
328  adcRaw2 = adc2;
329  }
330 
331  // Scale shared strip B ('adc') by ratio of peak of ADC counts from central strip A ('adc2')
332  // to sum of A and neighbouring maxima C ('adc1')
333 
334  for (size_t k = 0; k < 4; ++k) {
335  if (adc1[k] > 0 && adc2[k] > 0)
336  adc[k] = adc[k] * adc2[k] / (adc1[k] + adc2[k]);
337  if (adcRaw1[k] > 0 && adcRaw2[k] > 0)
338  adcRaw[k] = adcRaw[k] * adcRaw2[k] / (adcRaw1[k] + adcRaw2[k]);
339  }
340  prelimData = CSCStripHitData(thisStrip, tmax_cluster, adcRaw, adc);
341  }
342  }
343  }
344  return prelimData;
345 }
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
#define LogTrace(id)
PulseHeightMap thePulseHeightMap
static const double tmax[3]
std::vector< int > theMaxima
uint16_t *__restrict__ uint16_t const *__restrict__ adc

◆ runStrip()

std::vector< CSCStripHit > CSCHitFromStripOnly::runStrip ( const CSCDetId id,
const CSCLayer layer,
const CSCStripDigiCollection::Range rstripd 
)

fact (20.10.09);

L1A (Begin looping) Attempt to redefine theStrips, to encode L1A phase bits


L1A (end Looping)

L1A

Print statement to check StripHit content w/ LA1

Definition at line 54 of file CSCHitFromStripOnly.cc.

References funct::abs(), CSCLayer::chamber(), clusterSize, fillPulseHeights(), findMaxima(), gainWeight, ganged(), CSCChamberSpecs::gangedStrips(), triggerObjects_cff::id, id_, createfilelist::int, isDeadStrip(), phase1PixelTopology::layer, layer_, LogTrace, makeCluster(), nstrips_, recoConditions_, CSCDetId::ring(), setGanged(), CSCChamber::specs(), strips_adc, strips_adcRaw, CSCRecoConditions::stripWeights(), theClosestMaximum, theClusterSize, theConsecutiveStrips, theMaxima, theStrips, tmax_cluster, and useCalib.

Referenced by CSCRecHitDBuilder::build().

56  {
57  std::vector<CSCStripHit> hitsInLayer;
58 
59  // cache layer info for ease of access
60  id_ = id;
61  layer_ = layer;
62  nstrips_ = layer->chamber()->specs()->nStrips();
63 
64  setGanged(false);
65  if (id_.ring() == 4 && layer_->chamber()->specs()->gangedStrips())
66  setGanged(true); //@@ ONLY ME1/1A CAN BE GANGED
67 
68  LogTrace("CSCHitFromStripOnly") << "[CSCHitFromStripOnly::runStrip] id= " << id_ << " nstrips= " << nstrips_
69  << " ganged strips? " << ganged();
70 
71  tmax_cluster = 5;
72 
73  // Get gain correction weights for all strips in layer, and cache in gainWeight.
74  // They're used in fillPulseHeights below.
75  // When ME11a is ganged we only need the first 16 values of the 48 filled,
76  // but 17-48 are just duplicates of 1-16 anyway
77 
78  if (useCalib) {
80 
81  // *** START DUMP gainWeight
82  // std::cout << "gainWeight for id= " << id_ << " nstrips= " << nstrips_ << std::endl;
83  // for ( size_t i = 0; i!=10; ++i ) {
84  // for ( size_t j = 0; j!=8; ++j ) {
85  // std::cout << gainWeight[i*8 + j] << " ";
86  // }
87  // std::cout << std::endl;
88  // }
89  // *** END DUMP gainWeight
90  }
91 
92  // Store pulseheights from SCA and find maxima (potential hits)
93  fillPulseHeights(rstripd);
94  findMaxima(id);
95 
96  // Make a Strip Hit out of each strip local maximum
97  for (size_t imax = 0; imax != theMaxima.size(); ++imax) {
98  // Initialize parameters entering the CSCStripHit
100  theStrips.clear();
101  strips_adc.clear();
102  strips_adcRaw.clear();
103 
104  // makeCluster calls findHitOnStripPosition to determine the centroid position
105 
106  // Remember, the array starts at 0, but the stripId starts at 1...
107  float strippos = makeCluster(theMaxima[imax] + 1);
108 
109  //if ( strippos < 0 || tmax_cluster < 3 ){
110  // the strippos (as calculated here) is not used later on in
112  // with the negative charges allowed it can become negative
113  if (tmax_cluster < 3) {
114  theClosestMaximum.push_back(99); // to keep proper vector size
115  continue;
116  }
117  //---- If two maxima are too close the error assigned will be width/sqrt(12) - see CSCXonStrip_MatchGatti.cc
118  int maximum_to_left = 99; //---- If there is one maximum - the distance is set to 99 (strips)
119  int maximum_to_right = 99;
120  if (imax < theMaxima.size() - 1) {
121  maximum_to_right = theMaxima.at(imax + 1) - theMaxima.at(imax);
122  }
123  if (imax > 0) {
124  maximum_to_left = theMaxima.at(imax - 1) - theMaxima.at(imax);
125  }
126  if (std::abs(maximum_to_right) < std::abs(maximum_to_left)) {
127  theClosestMaximum.push_back(maximum_to_right);
128  } else {
129  theClosestMaximum.push_back(maximum_to_left);
130  }
131 
132  //---- Check if a neighbouring strip is a dead strip
133  //bool deadStrip = isNearDeadStrip(id, theMaxima.at(imax));
134  bool deadStripL = isDeadStrip(id, theMaxima.at(imax) - 1, nstrips_);
135  bool deadStripR = isDeadStrip(id, theMaxima.at(imax) + 1, nstrips_);
136  short int aDeadStrip = 0;
137  if (!deadStripL && !deadStripR) {
138  aDeadStrip = 0;
139  } else if (deadStripL && deadStripR) {
140  aDeadStrip = 255;
141  } else {
142  if (deadStripL) {
143  aDeadStrip = theMaxima.at(imax) - 1;
144  } else {
145  aDeadStrip = theMaxima.at(imax) + 1;
146  }
147  }
148 
151  std::vector<int> theL1AStrips;
152  for (int ila = 0; ila < (int)theStrips.size(); ila++) {
153  bool stripMatchCounter = false;
154  for (auto itl1 = rstripd.first; itl1 != rstripd.second; ++itl1) {
155  int stripNproto = (*itl1).getStrip();
156  if (!ganged()) {
157  if (theStrips[ila] == stripNproto) {
158  stripMatchCounter = true;
159  auto sz = (*itl1).getOverlappedSample().size();
160  int L1AbitOnPlace = 0;
161  for (auto iBit = 0UL; iBit < sz; iBit++) {
162  L1AbitOnPlace |= ((*itl1).getL1APhase(iBit) << (15 - iBit));
163  }
164  theL1AStrips.push_back(theStrips[ila] | L1AbitOnPlace);
165  }
166  } else {
167  for (int tripl = 0; tripl < 3; ++tripl) {
168  if (theStrips[ila] == (stripNproto + tripl * 16)) {
169  stripMatchCounter = true;
170  auto sz = (*itl1).getOverlappedSample().size();
171  int L1AbitOnPlace = 0;
172  for (auto iBit = 0UL; iBit < sz; iBit++) {
173  L1AbitOnPlace |= ((*itl1).getL1APhase(iBit) << (15 - iBit));
174  }
175  theL1AStrips.push_back(theStrips[ila] | L1AbitOnPlace);
176  }
177  }
178  }
179  }
180  if (!stripMatchCounter) {
181  theL1AStrips.push_back(theStrips[ila]);
182  }
183  }
185 
186  CSCStripHit striphit(id,
187  strippos,
188  tmax_cluster,
189  theL1AStrips,
190  strips_adc,
191  strips_adcRaw,
192  theConsecutiveStrips.at(imax),
193  theClosestMaximum.at(imax),
194  aDeadStrip);
195  hitsInLayer.push_back(striphit);
196  }
197 
199  /*
200  for(std::vector<CSCStripHit>::const_iterator itSHit=hitsInLayer.begin(); itSHit!=hitsInLayer.end(); ++itSHit){
201  (*itSHit).print();
202  }
203  */
204 
205  return hitsInLayer;
206 }
const CSCLayer * layer_
const CSCChamber * chamber() const
Definition: CSCLayer.h:49
void findMaxima(const CSCDetId &id)
Find local maxima.
std::vector< int > theConsecutiveStrips
bool isDeadStrip(const CSCDetId &id, int centralStrip, int nstrips)
Is the strip &#39;bad&#39;?
std::vector< int > theStrips
std::vector< int > theClosestMaximum
#define LogTrace(id)
constexpr std::array< uint8_t, layerIndexSize > layer
static const int theClusterSize
const CSCRecoConditions * recoConditions_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
float makeCluster(int centerStrip)
Make clusters using local maxima.
std::vector< float > strips_adcRaw
std::vector< int > theMaxima
std::vector< float > strips_adc
void stripWeights(const CSCDetId &id, short int nstrips, float *weights) const
bool gangedStrips() const
int ring() const
Definition: CSCDetId.h:68
void fillPulseHeights(const CSCStripDigiCollection::Range &rstripd)
Store SCA pulseheight information from strips in digis of one layer.
const CSCChamberSpecs * specs() const
Definition: CSCChamber.h:39

◆ setConditions()

void CSCHitFromStripOnly::setConditions ( const CSCRecoConditions reco)
inline

◆ setGanged()

void CSCHitFromStripOnly::setGanged ( bool  ig)
inline

Definition at line 49 of file CSCHitFromStripOnly.h.

References ganged_.

Referenced by runStrip().

49 { ganged_ = ig; }

Member Data Documentation

◆ calcped_

CSCPedestalChoice* CSCHitFromStripOnly::calcped_
private

◆ clusterSize

int CSCHitFromStripOnly::clusterSize
private

Definition at line 104 of file CSCHitFromStripOnly.h.

Referenced by makeCluster(), makeStripData(), and runStrip().

◆ gainWeight

float CSCHitFromStripOnly::gainWeight[80]
private

Definition at line 84 of file CSCHitFromStripOnly.h.

Referenced by fillPulseHeights(), and runStrip().

◆ ganged_

bool CSCHitFromStripOnly::ganged_
private

Definition at line 109 of file CSCHitFromStripOnly.h.

Referenced by ganged(), and setGanged().

◆ id_

CSCDetId CSCHitFromStripOnly::id_
private

Definition at line 78 of file CSCHitFromStripOnly.h.

Referenced by fillPulseHeights(), and runStrip().

◆ layer_

const CSCLayer* CSCHitFromStripOnly::layer_
private

Definition at line 79 of file CSCHitFromStripOnly.h.

Referenced by runStrip().

◆ nstrips_

unsigned CSCHitFromStripOnly::nstrips_
private

Definition at line 82 of file CSCHitFromStripOnly.h.

Referenced by findMaxima(), makeCluster(), and runStrip().

◆ recoConditions_

const CSCRecoConditions* CSCHitFromStripOnly::recoConditions_
private

◆ strips_adc

std::vector<float> CSCHitFromStripOnly::strips_adc
private

Definition at line 105 of file CSCHitFromStripOnly.h.

Referenced by findHitOnStripPosition(), and runStrip().

◆ strips_adcRaw

std::vector<float> CSCHitFromStripOnly::strips_adcRaw
private

Definition at line 106 of file CSCHitFromStripOnly.h.

Referenced by findHitOnStripPosition(), and runStrip().

◆ theClosestMaximum

std::vector<int> CSCHitFromStripOnly::theClosestMaximum
private

Definition at line 100 of file CSCHitFromStripOnly.h.

Referenced by findMaxima(), and runStrip().

◆ theClusterSize

const int CSCHitFromStripOnly::theClusterSize = 3
staticprivate

Definition at line 91 of file CSCHitFromStripOnly.h.

Referenced by makeCluster(), and runStrip().

◆ theConsecutiveStrips

std::vector<int> CSCHitFromStripOnly::theConsecutiveStrips
private

Definition at line 99 of file CSCHitFromStripOnly.h.

Referenced by findMaxima(), and runStrip().

◆ theMaxima

std::vector<int> CSCHitFromStripOnly::theMaxima
private

Definition at line 98 of file CSCHitFromStripOnly.h.

Referenced by findMaxima(), makeStripData(), and runStrip().

◆ thePulseHeightMap

PulseHeightMap CSCHitFromStripOnly::thePulseHeightMap
private

Definition at line 96 of file CSCHitFromStripOnly.h.

Referenced by fillPulseHeights(), findMaxima(), isPeakOK(), and makeStripData().

◆ theStrips

std::vector<int> CSCHitFromStripOnly::theStrips
private

Definition at line 107 of file CSCHitFromStripOnly.h.

Referenced by makeCluster(), and runStrip().

◆ theThresholdForAPeak

float CSCHitFromStripOnly::theThresholdForAPeak
private

Definition at line 92 of file CSCHitFromStripOnly.h.

Referenced by CSCHitFromStripOnly(), and isPeakOK().

◆ theThresholdForCluster

float CSCHitFromStripOnly::theThresholdForCluster
private

Definition at line 93 of file CSCHitFromStripOnly.h.

Referenced by CSCHitFromStripOnly(), and isPeakOK().

◆ tmax_cluster

int CSCHitFromStripOnly::tmax_cluster
private

Definition at line 103 of file CSCHitFromStripOnly.h.

Referenced by makeStripData(), and runStrip().

◆ useCalib

bool CSCHitFromStripOnly::useCalib
private

Definition at line 90 of file CSCHitFromStripOnly.h.

Referenced by CSCHitFromStripOnly(), fillPulseHeights(), and runStrip().