CMS 3D CMS Logo

List of all members | Public Member Functions | Private 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)
 
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 Member Functions

 HBHEStatusBitSetter (const HBHEStatusBitSetter &)
 
HBHEStatusBitSetteroperator= (const 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 15 of file HBHEStatusBitSetter.h.

Constructor & Destructor Documentation

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

Definition at line 11 of file HBHEStatusBitSetter.cc.

References edm::ParameterSet::getParameter(), muonDTDigis_cfi::pset, and pulseShapeParameters_.

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 }
T getParameter(std::string const &) const
const HcalFrontEndMap * frontEndMap_
std::vector< std::vector< double > > pulseShapeParameters_
HBHEStatusBitSetter::~HBHEStatusBitSetter ( )

Definition at line 28 of file HBHEStatusBitSetter.cc.

28 { }
HBHEStatusBitSetter::HBHEStatusBitSetter ( const HBHEStatusBitSetter )
private

Member Function Documentation

void HBHEStatusBitSetter::Clear ( )

Definition at line 42 of file HBHEStatusBitSetter.cc.

References hpdMultiplicity_, and mps_fire::i.

Referenced by HcalHitReconstructor::produce().

42  {
43  const unsigned sz = hpdMultiplicity_.size();
44  for (unsigned i=0; i<sz; i++)
45  hpdMultiplicity_[i] = 0;
46 }
std::vector< int > hpdMultiplicity_
HBHEStatusBitSetter& HBHEStatusBitSetter::operator= ( const HBHEStatusBitSetter )
private
void HBHEStatusBitSetter::rememberHit ( const HBHERecHit hbhe)

Definition at line 48 of file HBHEStatusBitSetter.cc.

References CaloRecHit::detid(), CaloRecHit::energy(), frontEndMap_, hitEnergyMinimum_, hpdMultiplicity_, diffTreeTool::index, 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 }
const DetId & detid() const
Definition: CaloRecHit.h:21
std::vector< int > hpdMultiplicity_
const int lookupRMIndex(DetId fId) const
const HcalFrontEndMap * frontEndMap_
float energy() const
Definition: CaloRecHit.h:17
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_, HcalCaloFlagLabels::HBHEPulseShape, nominalPedestal_, pulseShapeParameters_, rememberHit(), CaloRecHit::setFlagField(), HBHEDataFrame::size(), and findQualityFiles::size.

Referenced by HcalHitReconstructor::produce().

63  {
64  if (frontEndMap_==nullptr) {
65  edm::LogError("HBHEStatusBitSetter") << "No HcalFrontEndMap in SetFlagsFromDigi";
66  return;
67  }
68  rememberHit(hbhe);
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) continue;
90  // digi[i].nominal_fC() could be replaced by tool[iSlice], I think... -- Jeff, 11 April 2011
91  double qsum3=digi[iSlice].nominal_fC() + digi[iSlice-1].nominal_fC() + digi[iSlice-2].nominal_fC() - 3*nominalPedestal_;
92  if (qsum3>charge_max3) {
93  charge_max3=qsum3;
94  slice_max3=iSlice;
95  }
96  }
97 
98  if ((4+slice_max3)>size) return;
99  charge_late3=digi[slice_max3+1].nominal_fC() + digi[slice_max3+2].nominal_fC() + digi[slice_max3+3].nominal_fC() - 3*nominalPedestal_;
100 
101  for (unsigned int iCut=0;iCut<pulseShapeParameters_.size();iCut++) {
102  if (pulseShapeParameters_[iCut].size()!=6) continue;
103  if (nominal_charge_total<pulseShapeParameters_[iCut].at(0) || nominal_charge_total>=pulseShapeParameters_[iCut].at(1)) continue;
104  if ( charge_late3< (pulseShapeParameters_[iCut].at(2)+nominal_charge_total*pulseShapeParameters_[iCut].at(3)) ) continue;
105  if ( charge_late3>=(pulseShapeParameters_[iCut].at(4)+nominal_charge_total*pulseShapeParameters_[iCut].at(5)) ) continue;
107  return;
108  }
109 
110 }
size
Write out results.
int size() const
total number of samples in the digi
Definition: HBHEDataFrame.h:26
void setFlagField(uint32_t value, int base, int width=1)
Definition: CaloRecHit.cc:20
const HcalFrontEndMap * frontEndMap_
void rememberHit(const HBHERecHit &hbhe)
virtual void adc2fC(const HBHEDataFrame &df, CaloSamples &lf) const =0
std::vector< std::vector< double > > pulseShapeParameters_
void HBHEStatusBitSetter::SetFlagsFromRecHits ( HBHERecHitCollection rec)

Definition at line 113 of file HBHEStatusBitSetter.cc.

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

Referenced by HcalHitReconstructor::produce().

113  {
114 
115  if (frontEndMap_==nullptr) {
116  edm::LogError("HBHEStatusBitSetter") << "No HcalFrontEndMap in SetFlagsFromRecHits";
117  return;
118  }
119 
120  for (HBHERecHitCollection::iterator iHBHE=rec.begin();iHBHE!=rec.end();++iHBHE) {
121  const int index=frontEndMap_->lookupRMIndex(iHBHE->detid());
122  if (hpdMultiplicity_.at(index)<hitMultiplicityThreshold_) continue;
123  iHBHE->setFlagField(1,HcalCaloFlagLabels::HBHEHpdHitMultiplicity);
124  }
125 }
std::vector< int > hpdMultiplicity_
const int lookupRMIndex(DetId fId) const
const HcalFrontEndMap * frontEndMap_
std::vector< HBHERecHit >::iterator iterator
const_iterator end() const
const_iterator begin() const
void HBHEStatusBitSetter::SetFrontEndMap ( const HcalFrontEndMap m)

Definition at line 30 of file HBHEStatusBitSetter.cc.

References frontEndMap_, hpdMultiplicity_, funct::m, and HcalFrontEndMap::maxRMIndex().

Referenced by HcalHitReconstructor::beginRun().

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 HcalFrontEndMap * frontEndMap_
const int maxRMIndex() const

Member Data Documentation

const HcalFrontEndMap* HBHEStatusBitSetter::frontEndMap_
private
double HBHEStatusBitSetter::hitEnergyMinimum_
private

Definition at line 32 of file HBHEStatusBitSetter.h.

Referenced by HBHEStatusBitSetter(), and rememberHit().

int HBHEStatusBitSetter::hitMultiplicityThreshold_
private

Definition at line 33 of file HBHEStatusBitSetter.h.

Referenced by HBHEStatusBitSetter(), and SetFlagsFromRecHits().

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

Definition at line 36 of file HBHEStatusBitSetter.h.

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

double HBHEStatusBitSetter::nominalPedestal_
private

Definition at line 34 of file HBHEStatusBitSetter.h.

Referenced by HBHEStatusBitSetter(), and SetFlagsFromDigi().

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

Definition at line 37 of file HBHEStatusBitSetter.h.

Referenced by HBHEStatusBitSetter(), and SetFlagsFromDigi().