CMS 3D CMS Logo

HLTHcalLaserMisfireFilter.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: HLTHcalLaserMisfireFilter
4 // Class: HLTHcalLaserMisfireFilter
5 //
19 // system include files
20 #include <memory>
21 
22 // user include files
26 
33 
34 //
35 // constructors and destructor
36 //
38  inputHBHE_ = config.getParameter<edm::InputTag>("InputHBHE");
39  inputHF_ = config.getParameter<edm::InputTag>("InputHF");
40  minFracDiffHBHELaser_ = config.getParameter<double>("minFracDiffHBHELaser");
41  minFracHFLaser_ = config.getParameter<double>("minFracHFLaser");
42  minADCHBHE_ = config.getParameter<int>("minADCHBHE");
43  minADCHF_ = config.getParameter<int>("minADCHF"),
44  testMode_ = config.getUntrackedParameter<bool>("testMode",false);
45 
46  inputTokenHBHE_ = consumes<HBHEDigiCollection>(inputHBHE_);
47  inputTokenHF_ = consumes<QIE10DigiCollection>(inputHF_);
48 
49  }
50 
52 
55  desc.add<edm::InputTag>("InputHBHE",edm::InputTag("source"));
56  desc.add<edm::InputTag>("InputHF",edm::InputTag("source"));
57  desc.add<int>("minADCHBHE",10);
58  desc.add<int>("minADCHF",10);
59  desc.add<double>("minFracDiffHBHELaser",0.3);
60  desc.add<double>("minFracHFLaser",0.3);
61  desc.addUntracked<bool>("testMode",false);
62  descriptions.add("hltHcalLaserMisfireFilter",desc);
63 }
64 
66  const edm::EventSetup& iSetup) const {
67 
69  iEvent.getByToken(inputTokenHBHE_, hbhe_digi);
70 
72  iEvent.getByToken(inputTokenHF_, hf_digi);
73 
74  // Count digis in good, bad RBXes. ('bad' RBXes see no laser signal)
75  double badrbxfracHBHE(0), goodrbxfracHBHE(0), allrbxfracHF(0);
76  int NbadHBHE = 72*3; // 3 bad RBXes, 72 channels each
77  int NgoodHBHE= 2592*2-NbadHBHE; // remaining HBHE channels are 'good'
78  int NallHF = 864*4;
79 
80  for (auto hbhe = hbhe_digi->begin(); hbhe != hbhe_digi->end(); ++ hbhe){
81  const HBHEDataFrame digi = (const HBHEDataFrame)(*hbhe);
82  HcalDetId myid=(HcalDetId)digi.id();
83  bool isbad(false); // assume channel is not bad
84 
85  bool passCut(false);
86  int maxdigiHB(0);
87  for (int i=0; i<digi.size(); i++)
88  if(digi.sample(i).adc() > maxdigiHB) maxdigiHB = digi.sample(i).adc();
89  if (maxdigiHB > minADCHBHE_) passCut = true;
90 
91  // Three RBX's in HB do not receive any laser light (HBM5, HBM8, HBM9)
92  // They correspond to iphi = 15:18, 27:30, 31:34 respectively and
93  // ieta < 0
94  if ( myid.subdet()==HcalBarrel && myid.ieta()<0 ) {
95  if (myid.iphi()>=15 && myid.iphi()<=18) isbad=true;
96  else if (myid.iphi()>=27 && myid.iphi()<=34) isbad=true;
97  }
98 
99  if (passCut) {
100  if (isbad) {
101  badrbxfracHBHE += 1.;
102  } else goodrbxfracHBHE += 1.;
103  }
104  }
105  goodrbxfracHBHE /= NgoodHBHE;
106  badrbxfracHBHE /= NbadHBHE;
107 
108  for (auto hf = hf_digi->begin(); hf != hf_digi->end(); ++hf) {
109  const QIE10DataFrame digi = (const QIE10DataFrame)(*hf);
110  bool passCut(false);
111  int maxdigiHF(0);
112  for (int i=0; i<digi.samples(); i++)
113  if(digi[i].adc() > maxdigiHF) maxdigiHF = digi[i].adc();
114  if (maxdigiHF > minADCHF_) passCut = true;
115 
116  if (passCut) {
117  allrbxfracHF += 1.;
118  }
119  }
120  allrbxfracHF /= NallHF;
121 
122  if (testMode_)
123  edm::LogVerbatim("Report")
124  << "******************************************************************\n"
125  << "goodrbxfracHBHE: " << goodrbxfracHBHE << " badrbxfracHBHE: "
126  << badrbxfracHBHE << " Size " << hbhe_digi->size() << "\n"
127  << "allrbxfracHF: " << allrbxfracHF << " Size " << hf_digi->size()
128  << "\n******************************************************************";
129 
130  if (((goodrbxfracHBHE-badrbxfracHBHE) < minFracDiffHBHELaser_) ||
131  (allrbxfracHF < minFracHFLaser_)) return false;
132 
133  return true;
134 }
135 
136 // declare this class as a framework plugin
int samples() const
total number of samples in the digi
HLTHcalLaserMisfireFilter(const edm::ParameterSet &)
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
HcalSubdetector subdet() const
get the subdetector
Definition: HcalDetId.h:146
ParameterDescriptionBase * addUntracked(U const &iLabel, T const &value)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
int size() const
total number of samples in the digi
Definition: HBHEDataFrame.h:31
Definition: config.py:1
const_iterator begin() const
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
int ieta() const
get the cell ieta
Definition: HcalDetId.h:159
constexpr int adc(sample_type sample)
get the ADC sample (12 bits)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
constexpr int adc() const
get the ADC sample
Definition: HcalQIESample.h:59
const_iterator end() const
bool filter(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
int iphi() const
get the cell iphi
Definition: HcalDetId.h:161
edm::EDGetTokenT< HBHEDigiCollection > inputTokenHBHE_
HcalQIESample const & sample(int i) const
access a sample
Definition: HBHEDataFrame.h:44
edm::EDGetTokenT< QIE10DigiCollection > inputTokenHF_
void add(std::string const &label, ParameterSetDescription const &psetDescription)
const_iterator end() const
size_type size() const
const HcalDetId & id() const
Definition: HBHEDataFrame.h:27
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
const_iterator begin() const