CMS 3D CMS Logo

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

#include <ESPedestalClient.h>

Inheritance diagram for ESPedestalClient:
ESClient

Public Member Functions

void endJobAnalyze (DQMStore::IGetter &) override
 
 ESPedestalClient (const edm::ParameterSet &)
 
 ~ESPedestalClient () override
 
- Public Member Functions inherited from ESClient
virtual void endLumiAnalyze (DQMStore::IGetter &)
 
 ESClient (edm::ParameterSet const &)
 
template<typename T >
TgetHisto (MonitorElement *, bool=false, T *=0) const
 
void setup (DQMStore::IBooker &)
 
virtual ~ESClient ()
 

Private Member Functions

void book (DQMStore::IBooker &) override
 

Private Attributes

TF1 * fg_
 
bool fitPedestal_
 
MonitorElementhPed_ [2][2][40][40]
 
MonitorElementhTotN_ [2][2][40][40]
 
std::vector< int > senP_
 
std::vector< int > senX_
 
std::vector< int > senY_
 
std::vector< int > senZ_
 

Additional Inherited Members

- Public Types inherited from ESClient
typedef dqm::legacy::DQMStore DQMStore
 
typedef dqm::legacy::MonitorElement MonitorElement
 
- Protected Attributes inherited from ESClient
bool cloneME_
 
bool debug_
 
bool initialized_
 
std::string prefixME_
 
bool verbose_
 

Detailed Description

Definition at line 14 of file ESPedestalClient.h.

Constructor & Destructor Documentation

◆ ESPedestalClient()

ESPedestalClient::ESPedestalClient ( const edm::ParameterSet ps)

Definition at line 13 of file ESPedestalClient.cc.

References gather_cfg::cout, FCDTask_cfi::fiber, geometryDiff::file, edm::FileInPath::fullPath(), edm::ParameterSet::getUntrackedParameter(), hPed_, hTotN_, mps_fire::i, ALPAKA_ACCELERATOR_NAMESPACE::ecal::reconstruction::internal::endcap::ix(), ALPAKA_ACCELERATOR_NAMESPACE::ecal::reconstruction::internal::endcap::iy(), dqmiolumiharvest::j, dqmdumpme::k, visualization-live-secondInstance_cfg::m, submitDQMOfflineCAF::nLines, senP_, senX_, senY_, senZ_, and AlCaHLTBitMon_QueryRunRegistry::string.

14  : ESClient(ps),
15  fitPedestal_(ps.getUntrackedParameter<bool>("fitPedestal")),
16  fg_(new TF1("fg", "gaus")),
17  senZ_(),
18  senP_(),
19  senX_(),
20  senY_() {
21  for (int i = 0; i < 2; i++)
22  for (int j = 0; j < 2; j++)
23  for (int k = 0; k < 40; k++)
24  for (int m = 0; m < 40; m++) {
25  hPed_[i][j][k][m] = nullptr;
26  hTotN_[i][j][k][m] = nullptr;
27  }
28 
29  std::string lutPath(ps.getUntrackedParameter<edm::FileInPath>("LookupTable").fullPath());
30 
31  // read in look-up table
32  int nLines, iz, ip, ix, iy, fed, kchip, pace, bundle, fiber, optorx;
33  ifstream file(lutPath);
34 
35  if (file.is_open()) {
36  file >> nLines;
37 
38  for (int i = 0; i < nLines; ++i) {
39  file >> iz >> ip >> ix >> iy >> fed >> kchip >> pace >> bundle >> fiber >> optorx;
40 
41  senZ_.push_back(iz);
42  senP_.push_back(ip);
43  senX_.push_back(ix);
44  senY_.push_back(iy);
45  }
46 
47  file.close();
48  } else {
49  cout << "ESPedestalClient : Look up table file can not be found in " << lutPath << endl;
50  }
51 }
std::vector< int > senP_
std::string fullPath() const
Definition: FileInPath.cc:161
MonitorElement * hTotN_[2][2][40][40]
ESClient(edm::ParameterSet const &)
Definition: ESClient.cc:5
std::vector< int > senX_
T getUntrackedParameter(std::string const &, T const &) const
std::vector< int > senZ_
std::vector< int > senY_
ALPAKA_FN_ACC ALPAKA_FN_INLINE uint32_t ix(uint32_t id)
ALPAKA_FN_ACC ALPAKA_FN_INLINE uint32_t iy(uint32_t id)
MonitorElement * hPed_[2][2][40][40]

◆ ~ESPedestalClient()

ESPedestalClient::~ESPedestalClient ( )
override

Definition at line 53 of file ESPedestalClient.cc.

References fg_.

53 { delete fg_; }

Member Function Documentation

◆ book()

void ESPedestalClient::book ( DQMStore::IBooker _ibooker)
overrideprivatevirtual

Reimplemented from ESClient.

Definition at line 112 of file ESPedestalClient.cc.

References dqm::implementation::IBooker::book1D(), hPed_, hTotN_, mps_fire::i, ESClient::prefixME_, senP_, senX_, senY_, senZ_, and dqm::implementation::NavigatorBase::setCurrentFolder().

112  {
113  // define histograms
114  _ibooker.setCurrentFolder(prefixME_ + "/ESPedestalClient");
115 
116  char hname[300];
117  for (unsigned i = 0; i < senZ_.size(); ++i) {
118  int iz = (senZ_[i] == 1) ? 0 : 1;
119 
120  sprintf(hname, "Ped Z %d P %d X %d Y %d", senZ_[i], senP_[i], senX_[i], senY_[i]);
121  hPed_[iz][senP_[i] - 1][senX_[i] - 1][senY_[i] - 1] = _ibooker.book1D(hname, hname, 32, 0, 32);
122 
123  sprintf(hname, "Total Noise Z %d P %d X %d Y %d", senZ_[i], senP_[i], senX_[i], senY_[i]);
124  hTotN_[iz][senP_[i] - 1][senX_[i] - 1][senY_[i] - 1] = _ibooker.book1D(hname, hname, 32, 0, 32);
125  }
126 }
std::vector< int > senP_
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:36
MonitorElement * hTotN_[2][2][40][40]
std::vector< int > senX_
std::string prefixME_
Definition: ESClient.h:32
std::vector< int > senZ_
std::vector< int > senY_
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
MonitorElement * hPed_[2][2][40][40]

◆ endJobAnalyze()

void ESPedestalClient::endJobAnalyze ( DQMStore::IGetter _igetter)
overridevirtual

Reimplemented from ESClient.

Definition at line 55 of file ESPedestalClient.cc.

References gather_cfg::cout, ESClient::debug_, BTVHLTOfflineSource_cfi::dirname, fg_, fitPedestal_, dqm::implementation::IGetter::get(), dqm::legacy::MonitorElement::getTH1F(), hPed_, hTotN_, mps_fire::i, ESClient::prefixME_, senP_, senX_, senY_, senZ_, dqm::impl::MonitorElement::setBinContent(), and ESClient::verbose_.

55  {
56  if (debug_)
57  cout << "ESPedestalClient: endJob" << endl;
58 
59  // Preform pedestal fit
60  char hname[300];
61  int iz = 0;
62  if (fitPedestal_) {
63  if (verbose_)
64  cout << "ESPedestalClient: Fit Pedestal" << endl;
65 
66  for (unsigned i = 0; i < senZ_.size(); ++i) {
67  iz = (senZ_[i] == 1) ? 0 : 1;
68 
69  for (int is = 0; is < 32; ++is) {
70  string dirname = prefixME_ + "/ESPedestalTask/";
71  sprintf(hname, "ADC Z %d P %d X %d Y %d Str %d", senZ_[i], senP_[i], senX_[i], senY_[i], is + 1);
72  MonitorElement *meFit = _igetter.get(dirname + hname);
73 
74  if (meFit == nullptr)
75  continue;
76  TH1F *rootHisto = meFit->getTH1F();
77 
78  rootHisto->Fit(fg_, "Q", "", 500, 1800);
79  rootHisto->Fit(fg_,
80  "RQ",
81  "",
82  fg_->GetParameter(1) - 2. * fg_->GetParameter(2),
83  fg_->GetParameter(1) + 2. * fg_->GetParameter(2));
84  hPed_[iz][senP_[i] - 1][senX_[i] - 1][senY_[i] - 1]->setBinContent(is + 1, (int)(fg_->GetParameter(1) + 0.5));
85  hTotN_[iz][senP_[i] - 1][senX_[i] - 1][senY_[i] - 1]->setBinContent(is + 1, fg_->GetParameter(2));
86  }
87  }
88 
89  } else {
90  if (verbose_)
91  cout << "ESPedestalClient: Use Histogram Mean" << endl;
92 
93  for (unsigned i = 0; i < senZ_.size(); ++i) {
94  iz = (senZ_[i] == 1) ? 0 : 1;
95 
96  for (int is = 0; is < 32; ++is) {
97  string dirname = prefixME_ + "/ESPedestalTask/";
98  sprintf(hname, "ADC Z %d P %d X %d Y %d Str %d", senZ_[i], senP_[i], senX_[i], senY_[i], is + 1);
99  MonitorElement *meMean = _igetter.get(dirname + hname);
100 
101  if (meMean == nullptr)
102  continue;
103  TH1F *rootHisto = meMean->getTH1F();
104 
105  hPed_[iz][senP_[i] - 1][senX_[i] - 1][senY_[i] - 1]->setBinContent(is + 1, (int)(rootHisto->GetMean() + 0.5));
106  hTotN_[iz][senP_[i] - 1][senX_[i] - 1][senY_[i] - 1]->setBinContent(is + 1, rootHisto->GetRMS());
107  }
108  }
109  }
110 }
std::vector< int > senP_
MonitorElement * hTotN_[2][2][40][40]
std::vector< int > senX_
std::string prefixME_
Definition: ESClient.h:32
bool debug_
Definition: ESClient.h:35
std::vector< int > senZ_
std::vector< int > senY_
virtual void setBinContent(int binx, double content)
set content of bin (1-D)
bool verbose_
Definition: ESClient.h:34
virtual TH1F * getTH1F() const
virtual MonitorElement * get(std::string const &fullpath) const
Definition: DQMStore.cc:712
MonitorElement * hPed_[2][2][40][40]

Member Data Documentation

◆ fg_

TF1* ESPedestalClient::fg_
private

Definition at line 28 of file ESPedestalClient.h.

Referenced by endJobAnalyze(), and ~ESPedestalClient().

◆ fitPedestal_

bool ESPedestalClient::fitPedestal_
private

Definition at line 23 of file ESPedestalClient.h.

Referenced by endJobAnalyze().

◆ hPed_

MonitorElement* ESPedestalClient::hPed_[2][2][40][40]
private

Definition at line 25 of file ESPedestalClient.h.

Referenced by book(), endJobAnalyze(), and ESPedestalClient().

◆ hTotN_

MonitorElement* ESPedestalClient::hTotN_[2][2][40][40]
private

Definition at line 26 of file ESPedestalClient.h.

Referenced by book(), endJobAnalyze(), and ESPedestalClient().

◆ senP_

std::vector<int> ESPedestalClient::senP_
private

Definition at line 31 of file ESPedestalClient.h.

Referenced by book(), endJobAnalyze(), and ESPedestalClient().

◆ senX_

std::vector<int> ESPedestalClient::senX_
private

Definition at line 32 of file ESPedestalClient.h.

Referenced by book(), endJobAnalyze(), and ESPedestalClient().

◆ senY_

std::vector<int> ESPedestalClient::senY_
private

Definition at line 33 of file ESPedestalClient.h.

Referenced by book(), endJobAnalyze(), and ESPedestalClient().

◆ senZ_

std::vector<int> ESPedestalClient::senZ_
private

Definition at line 30 of file ESPedestalClient.h.

Referenced by book(), endJobAnalyze(), and ESPedestalClient().