CMS 3D CMS Logo

VpspScanAlgorithm.cc
Go to the documentation of this file.
6 #include "TProfile.h"
7 #include "TH1.h"
8 #include <iostream>
9 #include <sstream>
10 #include <iomanip>
11 #include <cmath>
12 
13 using namespace sistrip;
14 
15 // -----------------------------------------------------------------------------
16 //
18  : CommissioningAlgorithm(anal), histos_(2, Histo(nullptr, "")) {
19  ;
20 }
21 
22 // ----------------------------------------------------------------------------
23 //
24 void VpspScanAlgorithm::extract(const std::vector<TH1*>& histos) {
25  if (!anal()) {
26  edm::LogWarning(mlCommissioning_) << "[VpspScanAlgorithm::" << __func__ << "]"
27  << " NULL pointer to Analysis object!";
28  return;
29  }
30 
31  // Check number of histograms
32  if (histos.size() != 2) {
34  }
35 
36  // Extract FED key from histo title
37  if (!histos.empty()) {
38  anal()->fedKey(extractFedKey(histos.front()));
39  }
40 
41  // Extract histograms
42  std::vector<TH1*>::const_iterator ihis = histos.begin();
43  for (; ihis != histos.end(); ihis++) {
44  // Check pointer
45  if (!(*ihis)) {
46  continue;
47  }
48 
49  // Check name
50  SiStripHistoTitle title((*ihis)->GetName());
51  if (title.runType() != sistrip::VPSP_SCAN) {
53  continue;
54  }
55 
56  // Extract APV number
57  uint16_t apv = sistrip::invalid_;
58  if (title.extraInfo().find(sistrip::apv_) != std::string::npos) {
59  std::stringstream ss;
60  ss << title.extraInfo().substr(title.extraInfo().find(sistrip::apv_) + (sizeof(sistrip::apv_) - 1), 1);
61  ss >> std::dec >> apv;
62  }
63 
64  if (apv <= 1) {
65  histos_[apv].first = *ihis;
66  histos_[apv].second = (*ihis)->GetName();
67  } else {
69  }
70  }
71 }
72 
73 // -----------------------------------------------------------------------------
74 //
76  if (!anal()) {
77  edm::LogWarning(mlCommissioning_) << "[VpspScanAlgorithm::" << __func__ << "]"
78  << " NULL pointer to base Analysis object!";
79  return;
80  }
81 
83  VpspScanAnalysis* anal = dynamic_cast<VpspScanAnalysis*>(tmp);
84  if (!anal) {
85  edm::LogWarning(mlCommissioning_) << "[VpspScanAlgorithm::" << __func__ << "]"
86  << " NULL pointer to derived Analysis object!";
87  return;
88  }
89 
90  // from deprecated()
91 
92  std::vector<const TProfile*> histos;
93  std::vector<uint16_t> monitorables;
94  for (uint16_t iapv = 0; iapv < 2; iapv++) {
95  monitorables.clear();
96  monitorables.resize(7, sistrip::invalid_);
97 
98  histos.clear();
99  histos.push_back(const_cast<const TProfile*>(dynamic_cast<TProfile*>(histos_[iapv].first)));
100 
101  if (!histos[0]) {
103  continue;
104  }
105 
106  // Find "top" plateau
107  int first = sistrip::invalid_;
108  float top = -1. * sistrip::invalid_;
109  ;
110  for (int k = 5; k < 55; k++) {
111  if (!histos[0]->GetBinEntries(k)) {
112  continue;
113  }
114  if (histos[0]->GetBinContent(k) >= top) {
115  first = k;
116  top = histos[0]->GetBinContent(k);
117  }
118  }
119  if (top < -1. * sistrip::valid_) {
120  top = sistrip::invalid_;
121  } //@@ just want +ve invalid number here
122  if (top > 1. * sistrip::valid_) {
124  continue;
125  }
126  monitorables[5] = static_cast<uint16_t>(top);
127  monitorables[3] = first;
128 
129  // Find "bottom" plateau
130  int last = sistrip::invalid_;
131  float bottom = 1. * sistrip::invalid_;
132  for (int k = 55; k > 5; k--) {
133  if (!histos[0]->GetBinEntries(k)) {
134  continue;
135  }
136  if (histos[0]->GetBinContent(k) <= bottom) {
137  last = k;
138  bottom = histos[0]->GetBinContent(k);
139  }
140  }
141  if (bottom > 1. * sistrip::valid_) {
143  continue;
144  }
145  monitorables[6] = static_cast<uint16_t>(bottom);
146  monitorables[4] = last;
147 
148  // Set optimum baseline level
149  float opt = bottom + (top - bottom) * 1. / 3.;
150  monitorables[1] = static_cast<uint16_t>(opt);
151 
152  // Find optimum VPSP setting
153  uint16_t vpsp = sistrip::invalid_;
154  if (opt < 1. * sistrip::valid_) {
155  uint16_t ivpsp = 5;
156  for (; ivpsp < 55; ivpsp++) {
157  if (histos[0]->GetBinContent(ivpsp) < opt) {
158  break;
159  }
160  }
161  if (ivpsp != 54) {
162  vpsp = ivpsp;
163  monitorables[0] = vpsp;
164  } else {
166  continue;
167  }
168 
169  } else {
171  continue;
172  }
173 
174  // Set analysis results for both APVs
175  anal->vpsp_[iapv] = monitorables[0];
176  anal->adcLevel_[iapv] = monitorables[1];
177  anal->fraction_[iapv] = monitorables[2];
178  anal->topEdge_[iapv] = monitorables[3];
179  anal->bottomEdge_[iapv] = monitorables[4];
180  anal->topLevel_[iapv] = monitorables[5];
181  anal->bottomLevel_[iapv] = monitorables[6];
182  }
183 }
CommissioningAnalysis *const anal() const
static const char unexpectedTask_[]
Utility class that holds histogram title.
std::vector< Histo > histos_
static const char numberOfHistos_[]
static const char unexpectedExtraInfo_[]
static const char noBaselineLevel_[]
static const uint16_t valid_
Definition: Constants.h:17
sistrip classes
uint32_t extractFedKey(const TH1 *const)
void analyse() override
static const char noBottomPlateau_[]
static const char mlCommissioning_[]
virtual void addErrorCode(const std::string &error)
static const char noTopPlateau_[]
const uint32_t & fedKey() const
static const char apv_[]
Histogram-based analysis for VPSP scan.
void extract(const std::vector< TH1 *> &) override
std::pair< TH1 *, std::string > Histo
static const uint16_t invalid_
Definition: Constants.h:16
histos
Definition: combine.py:4
Abstract base for derived classes that provide analysis of commissioning histograms.
static const char noVpspSetting_[]
Log< level::Warning, false > LogWarning
tmp
align.sh
Definition: createJobs.py:716
static const char nullPtr_[]