CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
HBHEStatusBitSetter Class Reference

#include <HBHEStatusBitSetter.h>

Public Member Functions

void Clear ()
 
 HBHEStatusBitSetter ()
 
 HBHEStatusBitSetter (double nominalPedestal, double hitEnergyMinimum, int hitMultiplicityThreshold, const std::vector< edm::ParameterSet > &pulseShapeParameterSets)
 
 HBHEStatusBitSetter (const HBHEStatusBitSetter &)=delete
 
HBHEStatusBitSetteroperator= (const HBHEStatusBitSetter &)=delete
 
void rememberHit (const HBHERecHit &hbhe)
 
void SetFlagsFromDigi (HBHERecHit &hbhe, const HBHEDataFrame &digi, const HcalCoder &coder, const HcalCalibrations &calib)
 
void SetFlagsFromRecHits (HBHERecHitCollection &rec)
 
void SetFrontEndMap (const HcalFrontEndMap *m)
 
 ~HBHEStatusBitSetter ()
 

Private Attributes

const HcalFrontEndMapfrontEndMap_
 
double hitEnergyMinimum_
 
int hitMultiplicityThreshold_
 
std::vector< int > hpdMultiplicity_
 
double nominalPedestal_
 
std::vector< std::vector< double > > pulseShapeParameters_
 

Detailed Description

Definition at line 14 of file HBHEStatusBitSetter.h.

Constructor & Destructor Documentation

◆ HBHEStatusBitSetter() [1/3]

HBHEStatusBitSetter::HBHEStatusBitSetter ( )

Definition at line 5 of file HBHEStatusBitSetter.cc.

References hitEnergyMinimum_, hitMultiplicityThreshold_, and nominalPedestal_.

5  : frontEndMap_(nullptr) {
6  nominalPedestal_ = 3.0;
7  hitEnergyMinimum_ = 2.0;
9 }
const HcalFrontEndMap * frontEndMap_

◆ HBHEStatusBitSetter() [2/3]

HBHEStatusBitSetter::HBHEStatusBitSetter ( double  nominalPedestal,
double  hitEnergyMinimum,
int  hitMultiplicityThreshold,
const std::vector< edm::ParameterSet > &  pulseShapeParameterSets 
)

Definition at line 11 of file HBHEStatusBitSetter.cc.

References submitPVValidationJobs::params, muonDTDigis_cfi::pset, pulseShapeParameters_, and hltHbhereco_cfi::pulseShapeParameterSets.

18  frontEndMap_(nullptr) {
19  const unsigned sz = pulseShapeParameterSets.size();
20  pulseShapeParameters_.reserve(sz);
21  for (unsigned iPSet = 0; iPSet < sz; iPSet++) {
23  const std::vector<double>& params = pset.getParameter<std::vector<double> >("pulseShapeParameters");
24  pulseShapeParameters_.push_back(params);
25  }
26 }
const HcalFrontEndMap * frontEndMap_
std::vector< std::vector< double > > pulseShapeParameters_

◆ ~HBHEStatusBitSetter()

HBHEStatusBitSetter::~HBHEStatusBitSetter ( )

Definition at line 28 of file HBHEStatusBitSetter.cc.

28 {}

◆ HBHEStatusBitSetter() [3/3]

HBHEStatusBitSetter::HBHEStatusBitSetter ( const HBHEStatusBitSetter )
delete

Member Function Documentation

◆ Clear()

void HBHEStatusBitSetter::Clear ( )

Definition at line 42 of file HBHEStatusBitSetter.cc.

References hpdMultiplicity_, and mps_fire::i.

42  {
43  const unsigned sz = hpdMultiplicity_.size();
44  for (unsigned i = 0; i < sz; i++)
45  hpdMultiplicity_[i] = 0;
46 }
std::vector< int > hpdMultiplicity_

◆ operator=()

HBHEStatusBitSetter& HBHEStatusBitSetter::operator= ( const HBHEStatusBitSetter )
delete

◆ rememberHit()

void HBHEStatusBitSetter::rememberHit ( const HBHERecHit hbhe)

Definition at line 48 of file HBHEStatusBitSetter.cc.

References frontEndMap_, photonIsolationHIProducer_cfi::hbhe, hitEnergyMinimum_, hpdMultiplicity_, and HcalFrontEndMap::lookupRMIndex().

Referenced by SetFlagsFromDigi().

48  {
49  if (frontEndMap_ == nullptr) {
50  edm::LogError("HBHEStatusBitSetter") << "No HcalFrontEndMap in rememberHit";
51  return;
52  }
53  //increment hit multiplicity
54  if (hbhe.energy() > hitEnergyMinimum_) {
55  const int index = frontEndMap_->lookupRMIndex(hbhe.detid());
56  hpdMultiplicity_.at(index)++;
57  }
58 }
std::vector< int > hpdMultiplicity_
Log< level::Error, false > LogError
const HcalFrontEndMap * frontEndMap_
const int lookupRMIndex(DetId fId) const

◆ SetFlagsFromDigi()

void HBHEStatusBitSetter::SetFlagsFromDigi ( HBHERecHit hbhe,
const HBHEDataFrame digi,
const HcalCoder coder,
const HcalCalibrations calib 
)

Definition at line 60 of file HBHEStatusBitSetter.cc.

References HcalCoder::adc2fC(), frontEndMap_, photonIsolationHIProducer_cfi::hbhe, HcalCaloFlagLabels::HBHEPulseShape, nominalPedestal_, pulseShapeParameters_, rememberHit(), HBHEDataFrame::size(), and findQualityFiles::size.

63  {
64  if (frontEndMap_ == nullptr) {
65  edm::LogError("HBHEStatusBitSetter") << "No HcalFrontEndMap in SetFlagsFromDigi";
66  return;
67  }
69 
70  //set pulse shape bits
71  // Shuichi's algorithm uses the "correct" charge & pedestals, while Ted's uses "nominal" values.
72  // Perhaps we should correct Ted's algorithm in the future, though that will mean re-tuning thresholds for his cuts. -- Jeff, 28 May 2010
73  //double shuichi_charge_total=0.0;
74  double nominal_charge_total = 0.0;
75  double charge_max3 = -100.0;
76  double charge_late3 = -100.0;
77  unsigned int slice_max3 = 0;
78  unsigned int size = digi.size();
79 
80  CaloSamples tool;
81  coder.adc2fC(digi, tool);
82 
83  // int capid=-1;
84  for (unsigned int iSlice = 0; iSlice < size; iSlice++) {
85  // capid = digi.sample(iSlice).capid();
86  //shuichi_charge_total+=tool[iSlice]-calib.pedestal(capid);
87  nominal_charge_total += digi[iSlice].nominal_fC() - nominalPedestal_;
88 
89  if (iSlice < 2)
90  continue;
91  // digi[i].nominal_fC() could be replaced by tool[iSlice], I think... -- Jeff, 11 April 2011
92  double qsum3 = digi[iSlice].nominal_fC() + digi[iSlice - 1].nominal_fC() + digi[iSlice - 2].nominal_fC() -
93  3 * nominalPedestal_;
94  if (qsum3 > charge_max3) {
95  charge_max3 = qsum3;
96  slice_max3 = iSlice;
97  }
98  }
99 
100  if ((4 + slice_max3) > size)
101  return;
102  charge_late3 = digi[slice_max3 + 1].nominal_fC() + digi[slice_max3 + 2].nominal_fC() +
103  digi[slice_max3 + 3].nominal_fC() - 3 * nominalPedestal_;
104 
105  for (unsigned int iCut = 0; iCut < pulseShapeParameters_.size(); iCut++) {
106  if (pulseShapeParameters_[iCut].size() != 6)
107  continue;
108  if (nominal_charge_total < pulseShapeParameters_[iCut].at(0) ||
109  nominal_charge_total >= pulseShapeParameters_[iCut].at(1))
110  continue;
111  if (charge_late3 < (pulseShapeParameters_[iCut].at(2) + nominal_charge_total * pulseShapeParameters_[iCut].at(3)))
112  continue;
113  if (charge_late3 >= (pulseShapeParameters_[iCut].at(4) + nominal_charge_total * pulseShapeParameters_[iCut].at(5)))
114  continue;
115  hbhe.setFlagField(1, HcalCaloFlagLabels::HBHEPulseShape);
116  return;
117  }
118 }
size
Write out results.
Log< level::Error, false > LogError
const HcalFrontEndMap * frontEndMap_
void rememberHit(const HBHERecHit &hbhe)
std::vector< std::vector< double > > pulseShapeParameters_
virtual void adc2fC(const HBHEDataFrame &df, CaloSamples &lf) const =0
constexpr int size() const
total number of samples in the digi
Definition: HBHEDataFrame.h:27

◆ SetFlagsFromRecHits()

void HBHEStatusBitSetter::SetFlagsFromRecHits ( HBHERecHitCollection rec)

Definition at line 120 of file HBHEStatusBitSetter.cc.

References edm::SortedCollection< T, SORT >::begin(), edm::SortedCollection< T, SORT >::end(), frontEndMap_, HcalCaloFlagLabels::HBHEHpdHitMultiplicity, hitMultiplicityThreshold_, hpdMultiplicity_, and HcalFrontEndMap::lookupRMIndex().

120  {
121  if (frontEndMap_ == nullptr) {
122  edm::LogError("HBHEStatusBitSetter") << "No HcalFrontEndMap in SetFlagsFromRecHits";
123  return;
124  }
125 
126  for (HBHERecHitCollection::iterator iHBHE = rec.begin(); iHBHE != rec.end(); ++iHBHE) {
127  const int index = frontEndMap_->lookupRMIndex(iHBHE->detid());
129  continue;
130  iHBHE->setFlagField(1, HcalCaloFlagLabels::HBHEHpdHitMultiplicity);
131  }
132 }
std::vector< int > hpdMultiplicity_
Log< level::Error, false > LogError
const HcalFrontEndMap * frontEndMap_
const int lookupRMIndex(DetId fId) const
const_iterator begin() const
std::vector< T >::iterator iterator
const_iterator end() const

◆ SetFrontEndMap()

void HBHEStatusBitSetter::SetFrontEndMap ( const HcalFrontEndMap m)

Definition at line 30 of file HBHEStatusBitSetter.cc.

References frontEndMap_, hpdMultiplicity_, visualization-live-secondInstance_cfg::m, and HcalFrontEndMap::maxRMIndex().

30  {
31  frontEndMap_ = m;
32  hpdMultiplicity_.clear();
33  if (frontEndMap_) {
34  const int sz = frontEndMap_->maxRMIndex();
35  hpdMultiplicity_.reserve(sz);
36  for (int iRm = 0; iRm < sz; iRm++) {
37  hpdMultiplicity_.push_back(0);
38  }
39  }
40 }
std::vector< int > hpdMultiplicity_
const int maxRMIndex() const
const HcalFrontEndMap * frontEndMap_

Member Data Documentation

◆ frontEndMap_

const HcalFrontEndMap* HBHEStatusBitSetter::frontEndMap_
private

◆ hitEnergyMinimum_

double HBHEStatusBitSetter::hitEnergyMinimum_
private

Definition at line 36 of file HBHEStatusBitSetter.h.

Referenced by HBHEStatusBitSetter(), and rememberHit().

◆ hitMultiplicityThreshold_

int HBHEStatusBitSetter::hitMultiplicityThreshold_
private

Definition at line 37 of file HBHEStatusBitSetter.h.

Referenced by HBHEStatusBitSetter(), and SetFlagsFromRecHits().

◆ hpdMultiplicity_

std::vector<int> HBHEStatusBitSetter::hpdMultiplicity_
private

Definition at line 40 of file HBHEStatusBitSetter.h.

Referenced by Clear(), rememberHit(), SetFlagsFromRecHits(), and SetFrontEndMap().

◆ nominalPedestal_

double HBHEStatusBitSetter::nominalPedestal_
private

Definition at line 38 of file HBHEStatusBitSetter.h.

Referenced by HBHEStatusBitSetter(), and SetFlagsFromDigi().

◆ pulseShapeParameters_

std::vector<std::vector<double> > HBHEStatusBitSetter::pulseShapeParameters_
private

Definition at line 41 of file HBHEStatusBitSetter.h.

Referenced by HBHEStatusBitSetter(), and SetFlagsFromDigi().