CMS 3D CMS Logo

ApvLatencyAlgorithm.cc
Go to the documentation of this file.
6 #include "TProfile.h"
7 #include <iostream>
8 #include <cmath>
9 
10 using namespace sistrip;
11 
12 // ----------------------------------------------------------------------------
13 //
15  : CommissioningAlgorithm(anal), histo_(nullptr, "") {
16  ;
17 }
18 
19 // ----------------------------------------------------------------------------
20 //
21 void ApvLatencyAlgorithm::extract(const std::vector<TH1*>& histos) {
22  if (!anal()) {
23  edm::LogWarning(mlCommissioning_) << "[ApvLatencyAlgorithm::" << __func__ << "]"
24  << " NULL pointer to Analysis object!";
25  return;
26  }
27 
28  // Check
29  if (histos.size() != 1) {
31  }
32 
33  // Extract FED key from histo title
34  if (!histos.empty()) {
35  anal()->fedKey(extractFedKey(histos.front()));
36  }
37 
38  // Extract histograms
39  std::vector<TH1*>::const_iterator ihis = histos.begin();
40  for (; ihis != histos.end(); ihis++) {
41  // Check for NULL pointer
42  if (!(*ihis)) {
43  continue;
44  }
45 
46  // Check name
47  SiStripHistoTitle title((*ihis)->GetName());
48  if (title.runType() != sistrip::APV_LATENCY) {
50  continue;
51  }
52 
53  // Extract timing histo
54  histo_.first = *ihis;
55  histo_.second = (*ihis)->GetName();
56  }
57 }
58 
59 // ----------------------------------------------------------------------------
60 //
62  if (!anal()) {
63  edm::LogWarning(mlCommissioning_) << "[ApvLatencyAlgorithm::" << __func__ << "]"
64  << " NULL pointer to base Analysis object!";
65  return;
66  }
67 
69  ApvLatencyAnalysis* anal = dynamic_cast<ApvLatencyAnalysis*>(tmp);
70  if (!anal) {
71  edm::LogWarning(mlCommissioning_) << "[ApvLatencyAlgorithm::" << __func__ << "]"
72  << " NULL pointer to derived Analysis object!";
73  return;
74  }
75 
76  // was in deprecated()
77 
78  std::vector<const TProfile*> histos;
79  std::vector<unsigned short> monitorables;
80 
81  // was in analysis()
82 
83  histos.clear();
84  histos.push_back(const_cast<const TProfile*>(dynamic_cast<TProfile*>(histo_.first)));
85  if (!histos[0]) {
87  return;
88  }
89 
90  monitorables.clear();
91 
92  //LogDebug("Commissioning|Algorithm") << "[ApvLatencyAlgorithm::analysis]";
93 
94  //extract root histogram
95  //check
96  if (histos.size() != 1) {
97  // edm::LogWarning("Commissioning|Algorithm") << "[ApvLatencyAlgorithm::analysis]: Requires \"const std::vector<const TH1F*>& \" argument to have size 1. Actual size: " << histos.size() << ". Monitorables set to 0.";
98  monitorables.push_back(0);
99  return;
100  }
101  const TProfile* histo = histos[0];
102 
103  //monitorable
104  unsigned short latency;
105 
106  std::vector<unsigned short> binContent;
107  binContent.reserve((unsigned short)histo->GetNbinsX());
108  binContent.resize((unsigned short)histo->GetNbinsX(), 0);
109 
110  for (unsigned short k = 0; k < (unsigned short)histo->GetNbinsX(); k++) { // k is bin number
111 
112  //fill std::vector with histogram contents
113  binContent.push_back((unsigned int)(histo->GetBinContent(k)));
114  }
115 
116  //calculate median
117 
118  sort(binContent.begin(), binContent.end());
119 
120  //calculate mean and mean2 of the readout within cutoffs
121 
122  float meanNoise = 0.; //M.W method
123  float mean2Noise = 0.;
124 
125  for (unsigned short k = (unsigned short)(binContent.size() * .1); k < (unsigned short)(binContent.size() * .9); k++) {
126  meanNoise += binContent[k];
127  mean2Noise += binContent[k] * binContent[k];
128  ;
129  }
130 
131  meanNoise = meanNoise * binContent.size() * 0.8;
132  mean2Noise = mean2Noise * binContent.size() * 0.8;
133  float sigmaNoise = sqrt(fabs(meanNoise * meanNoise - mean2Noise));
134 
135  //loop to look for signal > 5* sigma_noise
136  unsigned short count = 0;
137  unsigned short maxlatency = 0;
138  unsigned int maxhits = 0;
139 
140  for (unsigned short k = 1; k < ((unsigned short)histo->GetNbinsX() + 1); k++) { // k is bin number
141  if (histo->GetBinContent((Int_t)k) > maxhits)
142  maxlatency = k - 1;
143  if ((float)histo->GetBinContent((Int_t)k) > (meanNoise + 5 * sigmaNoise)) {
144  latency = k - 1;
145  count++;
146  }
147  }
148 
149  if (!count) {
150  // LogDebug("Commissioning|Algorithm") << "[ApvLatencyAlgorithm::analysis]: Warning: no signal found > mean + 5*sigma(noise). Returning latency of highest number of recorded hits.";
151  latency = maxlatency;
152  }
153 
154  if (count > 1) {
155  // LogDebug("Commissioning|Algorithm") << "[ApvLatencyAlgorithm::analysis]: Warning: more than one signal found > mean + 5*sigma(noise). Returning latency of highest number of recorded hits.";
156  latency = maxlatency;
157  }
158 
159  //set monitorables
160  monitorables.clear();
161  monitorables.push_back(latency);
162 
163  anal->latency_ = monitorables[0];
164 }
CommissioningAnalysis *const anal() const
static const char unexpectedTask_[]
Utility class that holds histogram title.
static const char numberOfHistos_[]
sistrip classes
uint32_t extractFedKey(const TH1 *const)
static const char mlCommissioning_[]
T sqrt(T t)
Definition: SSEVec.h:23
virtual void addErrorCode(const std::string &error)
const uint32_t & fedKey() const
Analysis for APV latency scan.
void extract(const std::vector< TH1 *> &) override
const Histo & histo() const
histos
Definition: combine.py:4
Abstract base for derived classes that provide analysis of commissioning histograms.
Log< level::Warning, false > LogWarning
tmp
align.sh
Definition: createJobs.py:716
static const char nullPtr_[]