CMS 3D CMS Logo

PedsFullNoiseTask.cc
Go to the documentation of this file.
1 
3 
11 
12 #include "boost/lexical_cast.hpp"
13 
14 // -----------------------------------------------------------------------------
15 //
17  : CommissioningTask(dqm, conn, "PedsFullNoiseTask"), nstrips_(256) {
18  LogTrace(sistrip::mlDqmSource_) << "[PedsFullNoiseTask::" << __func__ << "]"
19  << " Constructing object...";
20  edm::ParameterSet params = pset.getParameter<edm::ParameterSet>("PedsFullNoiseParameters");
21  nskip_ = params.getParameter<int>("NrEvToSkipAtStart");
22  skipped_ = false;
23  nevpeds_ = params.getParameter<int>("NrEvForPeds");
24  pedsdone_ = false;
25  nadcnoise_ = params.getParameter<int>("NrPosBinsNoiseHist");
26  fillnoiseprofile_ = params.getParameter<bool>("FillNoiseProfile");
27  useavgcm_ = params.getParameter<bool>("UseAverageCommonMode");
28  usefloatpeds_ = params.getParameter<bool>("UseFloatPedestals");
29 }
30 
31 // -----------------------------------------------------------------------------
32 //
34  LogTrace(sistrip::mlDqmSource_) << "[PedsFullNoiseTask::" << __func__ << "]"
35  << " Destructing object...";
36 }
37 
38 // -----------------------------------------------------------------------------
39 //
41  LogTrace(sistrip::mlDqmSource_) << "[PedsFullNoiseTask::" << __func__ << "]";
42 
43  // pedestal profile histo
44  pedhist_.isProfile_ = true;
45  pedhist_.explicitFill_ = false;
46  if (!pedhist_.explicitFill_) {
47  pedhist_.vNumOfEntries_.resize(nstrips_, 0);
49  pedhist_.vSumOfSquares_.resize(nstrips_, 0);
50  }
54  fedKey(),
56  connection().lldChannel(),
58  .title();
59  pedhist_.histo(dqm()->bookProfile(titleped, titleped, nstrips_, -0.5, nstrips_ * 1. - 0.5, 1025, 0., 1025.));
60 
61  // Noise profile
62  noiseprof_.isProfile_ = true;
63  noiseprof_.explicitFill_ = false;
68  }
72  fedKey(),
74  connection().lldChannel(),
76  .title();
77  noiseprof_.histo(dqm()->bookProfile(titlenoise, titlenoise, nstrips_, -0.5, nstrips_ * 1. - 0.5, 1025, 0., 1025.));
78 
79  // noise 2D compact histo
80  noisehist_.explicitFill_ = false;
82  noisehist_.vNumOfEntries_.resize((nstrips_ + 2) * 2 * (nadcnoise_ + 2), 0);
83  }
87  fedKey(),
89  connection().lldChannel(),
91  .title();
92  noisehist_.histo(dqm()->book2S(
93  titlenoise2d, titlenoise2d, 2 * nadcnoise_, -nadcnoise_, nadcnoise_, nstrips_, -0.5, nstrips_ * 1. - 0.5));
94  hist2d_ = (TH2S*)noisehist_.histo()->getTH2S();
95 }
96 
97 // -----------------------------------------------------------------------------
98 //
100  // Check number of digis
101  uint16_t nbins = digis.data.size();
102  if (nbins != nstrips_) {
103  edm::LogWarning(sistrip::mlDqmSource_) << "[PedsFullNoiseTask::" << __func__ << "]"
104  << " " << nstrips_ << " digis expected, but got " << nbins << ". Skipping.";
105  return;
106  }
107 
108  // get the event number of the first event, not necessarily 1 (parallel processing on FUs)
109  static int32_t firstev = summary.event();
110 
111  // skipping events
112  if (!skipped_) {
113  if (static_cast<int32_t>(summary.event()) - firstev < nskip_) {
114  return;
115  } else { // when all events are skipped
116  skipped_ = true;
117  if (nskip_ > 0)
118  LogTrace(sistrip::mlDqmSource_) << "[PedsFullNoiseTask::" << __func__ << "]"
119  << " Done skipping events. Now starting pedestals.";
120  }
121  }
122 
123  // determine pedestals - decoupled from noise determination
124  if (!pedsdone_) {
125  if (static_cast<int32_t>(summary.event()) - firstev < nskip_ + nevpeds_) {
126  // estimate the pedestals
127  for (uint16_t istrip = 0; istrip < nstrips_; ++istrip) {
128  updateHistoSet(pedhist_, istrip, digis.data[istrip].adc());
129  }
130  return;
131  } else { // when pedestals are done
132  pedsdone_ = true;
133  // cache the pedestal values for use in the 2D noise estimation
134  peds_.clear();
135  pedsfl_.clear();
136  for (uint16_t iapv = 0; iapv < 2; ++iapv) {
137  for (uint16_t ibin = 0; ibin < 128; ++ibin) {
138  uint16_t istrip = (iapv * 128) + ibin;
139  if (usefloatpeds_) {
140  pedsfl_.push_back(1. * pedhist_.vSumOfContents_.at(istrip) / pedhist_.vNumOfEntries_.at(istrip));
141  } else {
142  peds_.push_back(
143  static_cast<int16_t>(1. * pedhist_.vSumOfContents_.at(istrip) / pedhist_.vNumOfEntries_.at(istrip)));
144  }
145  }
146  }
147  LogTrace(sistrip::mlDqmSource_) << "[PedsFullNoiseTask::" << __func__ << "]"
148  << " Rough pedestals done. Now starting noise measurements.";
149  }
150  }
151 
152  // fill (or not) the old-style niose profile
153  if (fillnoiseprofile_) {
154  // Calc common mode for both APVs
155  std::vector<int32_t> cm;
156  cm.resize(2, 0);
157  std::vector<uint16_t> adc;
158  for (uint16_t iapv = 0; iapv < 2; iapv++) {
159  adc.clear();
160  adc.reserve(128);
161  for (uint16_t ibin = 0; ibin < 128; ibin++) {
162  if ((iapv * 128) + ibin < nbins) {
163  adc.push_back(digis.data.at((iapv * 128) + ibin).adc());
164  }
165  }
166  sort(adc.begin(), adc.end());
167  // take median as common mode
168  uint16_t index = adc.size() % 2 ? adc.size() / 2 : adc.size() / 2 - 1;
169  cm[iapv] = static_cast<int16_t>(adc[index]);
170  }
171  // 1D noise profile - see also further processing in the update() method
172  for (uint16_t istrip = 0; istrip < nstrips_; ++istrip) {
173  // calculate the noise in the old way, by subtracting the common mode, but without pedestal subtraction
174  int16_t noiseval = static_cast<int16_t>(digis.data.at(istrip).adc()) - cm[istrip / 128];
175  updateHistoSet(noiseprof_, istrip, noiseval);
176  }
177  }
178 
179  // 2D noise histogram
180  std::vector<int16_t> noisevals, noisevalssorted;
181  std::vector<float> noisevalsfl, noisevalssortedfl;
182  for (uint16_t iapv = 0; iapv < 2; ++iapv) {
183  float totadc = 0;
184  noisevals.clear();
185  noisevalsfl.clear();
186  noisevalssorted.clear();
187  noisevalssortedfl.clear();
188  for (uint16_t ibin = 0; ibin < 128; ++ibin) {
189  uint16_t istrip = (iapv * 128) + ibin;
190  // calculate the noise after subtracting the pedestal
191  if (usefloatpeds_) { // if float pedestals -> before FED processing
192  noisevalsfl.push_back(static_cast<float>(digis.data.at(istrip).adc()) - pedsfl_.at(istrip));
193  // now we still have a possible constant shift of the adc values with respect to 0, so we prepare to calculate the median of this shift
194  if (useavgcm_) { // if average CM -> before FED processing
195  totadc += noisevalsfl[ibin];
196  } else { // if median CM -> after FED processing
197  noisevalssortedfl.push_back(noisevalsfl[ibin]);
198  }
199  } else { // if integer pedestals -> after FED processing
200  noisevals.push_back(static_cast<int16_t>(digis.data.at(istrip).adc()) - peds_.at(istrip));
201  // now we still have a possible constant shift of the adc values with respect to 0, so we prepare to calculate the median of this shift
202  if (useavgcm_) { // if average CM -> before FED processing
203  totadc += noisevals[ibin];
204  } else { // if median CM -> after FED processing
205  noisevalssorted.push_back(noisevals[ibin]);
206  }
207  }
208  }
209  // calculate the common mode shift to apply
210  float cmshift = 0;
211  if (useavgcm_) { // if average CM -> before FED processing
212  if (usefloatpeds_) { // if float pedestals -> before FED processing
213  cmshift = totadc / 128;
214  } else { // if integer pedestals -> after FED processing
215  cmshift = static_cast<int16_t>(totadc / 128);
216  }
217  } else { // if median CM -> after FED processing
218  if (usefloatpeds_) { // if float pedestals -> before FED processing
219  // get the median common mode
220  sort(noisevalssortedfl.begin(), noisevalssortedfl.end());
221  uint16_t index = noisevalssortedfl.size() % 2 ? noisevalssortedfl.size() / 2 : noisevalssortedfl.size() / 2 - 1;
222  cmshift = noisevalssortedfl[index];
223  } else { // if integer pedestals -> after FED processing
224  // get the median common mode
225  sort(noisevalssorted.begin(), noisevalssorted.end());
226  uint16_t index = noisevalssorted.size() % 2 ? noisevalssorted.size() / 2 : noisevalssorted.size() / 2 - 1;
227  cmshift = noisevalssorted[index];
228  }
229  }
230  // now loop again to calculate the CM+pedestal subtracted noise values
231  for (uint16_t ibin = 0; ibin < 128; ++ibin) {
232  uint16_t istrip = (iapv * 128) + ibin;
233  // subtract the remaining common mode after subtraction of the rough pedestals
234  float noiseval = (usefloatpeds_ ? noisevalsfl[ibin] : noisevals[ibin]) - cmshift;
235  // retrieve the linear binnr through the histogram
236  uint32_t binnr = hist2d_->GetBin(static_cast<int>(noiseval + nadcnoise_), istrip + 1);
237  // store the noise value in the 2D histo
238  updateHistoSet(noisehist_, binnr); // no value, so weight 1
239  }
240  }
241 }
242 
243 // -----------------------------------------------------------------------------
244 //
246  // pedestals
248 
249  if (fillnoiseprofile_) {
250  // noise profile (does not use HistoSet directly, as want to plot noise as "contents", not "error")
252  for (uint16_t ii = 0; ii < noiseprof_.vNumOfEntries_.size(); ++ii) {
253  float mean = 0.;
254  float spread = 0.;
255  float entries = noiseprof_.vNumOfEntries_[ii];
256  if (entries > 0.) {
257  mean = noiseprof_.vSumOfContents_[ii] / entries;
258  spread = sqrt(noiseprof_.vSumOfSquares_[ii] / entries -
259  mean * mean); // nice way to calculate std dev: Sum (x-<x>)^2 / N
260  }
261  float error = spread / sqrt(entries); // uncertainty on std.dev. when no uncertainty on mean
262  UpdateTProfile::setBinContent(histo, ii + 1, entries, spread, error);
263  }
264  }
265 
266  // noise 2D histo
268 }
269 // -----------------------------------------------------------------------------
T getParameter(std::string const &) const
std::vector< float > vNumOfEntries_
Utility class that holds histogram title.
const std::string & title() const
static const char mlDqmSource_[]
std::vector< float > vSumOfContents_
void updateHistoSet(HistoSet &, const uint32_t &bin, const float &value)
static const char noise2D_[]
DQMStore *const dqm() const
Class containning control, module, detector and connection information, at the level of a FED channel...
const uint32_t & event() const
T sqrt(T t)
Definition: SSEVec.h:19
static void setBinContent(TProfile *const profile, const uint32_t &bin, const double &entries, const double &mean, const double &spread)
void update() override
constexpr int adc(sample_type sample)
get the ADC sample (12 bits)
static const char noiseProfile_[]
#define LogTrace(id)
ii
Definition: cuy.py:590
int extract(std::vector< int > *output, const std::string &dati)
static const char pedestals_[]
CompactHistoSet noisehist_
std::vector< int16_t > peds_
collection_type data
Definition: DetSet.h:81
PedsFullNoiseTask(DQMStore *dqm, const FedChannelConnection &conn, const edm::ParameterSet &pset)
void histo(MonitorElement *)
~PedsFullNoiseTask() override
const uint32_t & fedKey() const
std::vector< double > vSumOfSquares_
void book() override
const FedChannelConnection & connection() const
void fill(const SiStripEventSummary &, const edm::DetSet< SiStripRawDigi > &) override
std::vector< float > pedsfl_