CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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)
 
void SetFlagsFromDigi (const HcalTopology *topo, HBHERecHit &hbhe, const HBHEDataFrame &digi, const HcalCoder &coder, const HcalCalibrations &calib, int firstSample=3, int samplesToAdd=4)
 
void SetFlagsFromRecHits (const HcalTopology *topo, HBHERecHitCollection &rec)
 
 ~HBHEStatusBitSetter ()
 

Private Attributes

unsigned int firstSample_
 
double hitEnergyMinimum_
 
int hitMultiplicityThreshold_
 
std::vector< int > hpdMultiplicity_
 
HcalLogicalMaplogicalMap_
 
double nominalPedestal_
 
std::vector< std::vector
< double > > 
pulseShapeParameters_
 
unsigned int samplesToAdd_
 

Detailed Description

Definition at line 14 of file HBHEStatusBitSetter.h.

Constructor & Destructor Documentation

HBHEStatusBitSetter::HBHEStatusBitSetter ( )

Definition at line 6 of file HBHEStatusBitSetter.cc.

References hitEnergyMinimum_, hitMultiplicityThreshold_, hpdMultiplicity_, logicalMap_, HcalFrontEndId::maxRmIndex, and nominalPedestal_.

7 {
8  logicalMap_=0;
9 
10  for (int iRm=0;iRm<HcalFrontEndId::maxRmIndex;iRm++) {
11  hpdMultiplicity_.push_back(0);
12  }
13 
14  nominalPedestal_=3.0;
17 }
std::vector< int > hpdMultiplicity_
static const int maxRmIndex
HcalLogicalMap * logicalMap_
HBHEStatusBitSetter::HBHEStatusBitSetter ( double  nominalPedestal,
double  hitEnergyMinimum,
int  hitMultiplicityThreshold,
const std::vector< edm::ParameterSet > &  pulseShapeParameterSets 
)

Definition at line 19 of file HBHEStatusBitSetter.cc.

References edm::ParameterSet::getParameter(), hitEnergyMinimum_, hitMultiplicityThreshold_, hpdMultiplicity_, logicalMap_, HcalFrontEndId::maxRmIndex, nominalPedestal_, and pulseShapeParameters_.

24 {
25  logicalMap_=0;
26  for (int iRm=0;iRm<HcalFrontEndId::maxRmIndex;iRm++) {
27  hpdMultiplicity_.push_back(0);
28  }
29 
30  nominalPedestal_=nominalPedestal;
31  hitEnergyMinimum_=hitEnergyMinimum;
32  hitMultiplicityThreshold_=hitMultiplicityThreshold;
33 
34  for (unsigned int iPSet=0;iPSet<pulseShapeParameterSets.size();iPSet++) {
35  edm::ParameterSet pset=pulseShapeParameterSets.at(iPSet);
36  std::vector<double> params=pset.getParameter<std::vector<double> >("pulseShapeParameters");
37  pulseShapeParameters_.push_back(params);
38  }
39 
40 }
T getParameter(std::string const &) const
std::vector< int > hpdMultiplicity_
static const int maxRmIndex
std::vector< std::vector< double > > pulseShapeParameters_
HcalLogicalMap * logicalMap_
HBHEStatusBitSetter::~HBHEStatusBitSetter ( )

Definition at line 42 of file HBHEStatusBitSetter.cc.

References logicalMap_.

42  {
43  if (logicalMap_!=0) delete logicalMap_;
44 }
HcalLogicalMap * logicalMap_

Member Function Documentation

void HBHEStatusBitSetter::Clear ( )

Definition at line 46 of file HBHEStatusBitSetter.cc.

References hpdMultiplicity_, and i.

Referenced by HcalHitReconstructor::produce().

47 {
48  for (unsigned int i=0;i<hpdMultiplicity_.size();i++) hpdMultiplicity_[i]=0;
49 }
int i
Definition: DBlmapReader.cc:9
std::vector< int > hpdMultiplicity_
void HBHEStatusBitSetter::SetFlagsFromDigi ( const HcalTopology topo,
HBHERecHit hbhe,
const HBHEDataFrame digi,
const HcalCoder coder,
const HcalCalibrations calib,
int  firstSample = 3,
int  samplesToAdd = 4 
)

Definition at line 51 of file HBHEStatusBitSetter.cc.

References HcalCoder::adc2fC(), asciidump::at, HcalLogicalMapGenerator::createMap(), CaloRecHit::detid(), CaloRecHit::energy(), castor_dqm_sourceclient_file_cfg::firstSample, firstSample_, relval_steps::gen(), HcalLogicalMap::getHcalFrontEndId(), HcalCaloFlagLabels::HBHEPulseShape, hitEnergyMinimum_, hpdMultiplicity_, getHLTprescales::index, logicalMap_, nominalPedestal_, pulseShapeParameters_, samplesToAdd_, CaloRecHit::setFlagField(), HBHEDataFrame::size(), and findQualityFiles::size.

Referenced by HcalHitReconstructor::produce().

58 {
59  if (logicalMap_==0) {
61  logicalMap_=new HcalLogicalMap(gen.createMap(topo));
62  }
63 
64 
65  // get firstSample, samplesToAdd from database for each hit
67  samplesToAdd_ = samplesToAdd;
68 
69  //increment hit multiplicity
70  if (hbhe.energy()>hitEnergyMinimum_) {
71  int index=logicalMap_->getHcalFrontEndId(hbhe.detid()).rmIndex();
72  hpdMultiplicity_.at(index)++;
73  }
74 
75  //set pulse shape bits
76  // Shuichi's algorithm uses the "correct" charge & pedestals, while Ted's uses "nominal" values.
77  // Perhaps we should correct Ted's algorithm in the future, though that will mean re-tuning thresholds for his cuts. -- Jeff, 28 May 2010
78  //double shuichi_charge_total=0.0;
79  double nominal_charge_total=0.0;
80  double charge_max3=-100.0;
81  double charge_late3=-100.0;
82  unsigned int slice_max3=0;
83  unsigned int size=digi.size();
84 
85  CaloSamples tool;
86  coder.adc2fC(digi,tool);
87 
88  // int capid=-1;
89  for (unsigned int iSlice=0;iSlice<size;iSlice++)
90  {
91  // capid = digi.sample(iSlice).capid();
92  //shuichi_charge_total+=tool[iSlice]-calib.pedestal(capid);
93  nominal_charge_total+=digi[iSlice].nominal_fC()-nominalPedestal_;
94 
95  if (iSlice<2) continue;
96  // digi[i].nominal_fC() could be replaced by tool[iSlice], I think... -- Jeff, 11 April 2011
97  double qsum3=digi[iSlice].nominal_fC() + digi[iSlice-1].nominal_fC() + digi[iSlice-2].nominal_fC() - 3*nominalPedestal_;
98  if (qsum3>charge_max3) {
99  charge_max3=qsum3;
100  slice_max3=iSlice;
101  }
102  }
103 
104  if ((4+slice_max3)>size) return;
105  charge_late3=digi[slice_max3+1].nominal_fC() + digi[slice_max3+2].nominal_fC() + digi[slice_max3+3].nominal_fC() - 3*nominalPedestal_;
106 
107  for (unsigned int iCut=0;iCut<pulseShapeParameters_.size();iCut++) {
108  if (pulseShapeParameters_[iCut].size()!=6) continue;
109  if (nominal_charge_total<pulseShapeParameters_[iCut].at(0) || nominal_charge_total>=pulseShapeParameters_[iCut].at(1)) continue;
110  if ( charge_late3< (pulseShapeParameters_[iCut].at(2)+nominal_charge_total*pulseShapeParameters_[iCut].at(3)) ) continue;
111  if ( charge_late3>=(pulseShapeParameters_[iCut].at(4)+nominal_charge_total*pulseShapeParameters_[iCut].at(5)) ) continue;
113  return;
114  }
115 
116 }
const DetId & detid() const
Definition: CaloRecHit.h:20
int size() const
total number of samples in the digi
Definition: HBHEDataFrame.h:26
std::vector< int > hpdMultiplicity_
void setFlagField(uint32_t value, int base, int width=1)
Definition: CaloRecHit.cc:20
float energy() const
Definition: CaloRecHit.h:17
HcalLogicalMap createMap(const HcalTopology *topo, unsigned int mapIOV=4)
const HcalFrontEndId getHcalFrontEndId(const DetId &)
std::vector< std::vector< double > > pulseShapeParameters_
virtual void adc2fC(const HBHEDataFrame &df, CaloSamples &lf) const =0
HcalLogicalMap * logicalMap_
tuple size
Write out results.
list at
Definition: asciidump.py:428
void HBHEStatusBitSetter::SetFlagsFromRecHits ( const HcalTopology topo,
HBHERecHitCollection rec 
)

Definition at line 118 of file HBHEStatusBitSetter.cc.

References edm::SortedCollection< T, SORT >::begin(), HcalLogicalMapGenerator::createMap(), edm::SortedCollection< T, SORT >::end(), relval_steps::gen(), HcalLogicalMap::getHcalFrontEndId(), HcalCaloFlagLabels::HBHEHpdHitMultiplicity, hitMultiplicityThreshold_, hpdMultiplicity_, getHLTprescales::index, and logicalMap_.

Referenced by HcalHitReconstructor::produce().

118  {
119 
120  if (logicalMap_==0) {
122  logicalMap_=new HcalLogicalMap(gen.createMap(topo));
123  }
124 
125 
126  for (HBHERecHitCollection::iterator iHBHE=rec.begin();iHBHE!=rec.end();++iHBHE) {
127  int index=logicalMap_->getHcalFrontEndId(iHBHE->detid()).rmIndex();
128  if (hpdMultiplicity_.at(index)<hitMultiplicityThreshold_) continue;
129  iHBHE->setFlagField(1,HcalCaloFlagLabels::HBHEHpdHitMultiplicity);
130  }
131 }
std::vector< int > hpdMultiplicity_
HcalLogicalMap createMap(const HcalTopology *topo, unsigned int mapIOV=4)
std::vector< HBHERecHit >::iterator iterator
const_iterator end() const
const HcalFrontEndId getHcalFrontEndId(const DetId &)
HcalLogicalMap * logicalMap_
const_iterator begin() const

Member Data Documentation

unsigned int HBHEStatusBitSetter::firstSample_
private

Definition at line 29 of file HBHEStatusBitSetter.h.

Referenced by SetFlagsFromDigi().

double HBHEStatusBitSetter::hitEnergyMinimum_
private

Definition at line 27 of file HBHEStatusBitSetter.h.

Referenced by HBHEStatusBitSetter(), and SetFlagsFromDigi().

int HBHEStatusBitSetter::hitMultiplicityThreshold_
private

Definition at line 28 of file HBHEStatusBitSetter.h.

Referenced by HBHEStatusBitSetter(), and SetFlagsFromRecHits().

std::vector<int> HBHEStatusBitSetter::hpdMultiplicity_
private
HcalLogicalMap* HBHEStatusBitSetter::logicalMap_
private
double HBHEStatusBitSetter::nominalPedestal_
private

Definition at line 31 of file HBHEStatusBitSetter.h.

Referenced by HBHEStatusBitSetter(), and SetFlagsFromDigi().

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

Definition at line 34 of file HBHEStatusBitSetter.h.

Referenced by HBHEStatusBitSetter(), and SetFlagsFromDigi().

unsigned int HBHEStatusBitSetter::samplesToAdd_
private

Definition at line 30 of file HBHEStatusBitSetter.h.

Referenced by SetFlagsFromDigi().