CMS 3D CMS Logo

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

#include <PedsFullNoiseTask.h>

Inheritance diagram for PedsFullNoiseTask:
CommissioningTask

Public Member Functions

 PedsFullNoiseTask (DQMStore *dqm, const FedChannelConnection &conn, const edm::ParameterSet &pset)
 
 ~PedsFullNoiseTask () override
 
- Public Member Functions inherited from CommissioningTask
void bookHistograms ()
 
 CommissioningTask (DQMStore *, const FedChannelConnection &, const std::string &my_name)
 
void eventSetup (const edm::EventSetup *)
 
const uint32_t & fillCntr () const
 
void fillHistograms (const SiStripEventSummary &, const edm::DetSet< SiStripRawDigi > &)
 
void fillHistograms (const SiStripEventSummary &, const uint16_t &fed_id, const std::map< uint16_t, float > &fed_ch)
 
const std::string & myName () const
 
const uint32_t & updateFreq () const
 
void updateFreq (const uint32_t &)
 
void updateHistograms ()
 
virtual ~CommissioningTask ()
 

Private Member Functions

void book () override
 
void fill (const SiStripEventSummary &, const edm::DetSet< SiStripRawDigi > &) override
 
void update () override
 

Private Attributes

bool fillnoiseprofile_
 
TH2S * hist2d_
 
uint16_t nadcnoise_
 
uint16_t nevpeds_
 
CompactHistoSet noisehist_
 
HistoSet noiseprof_
 
uint16_t nskip_
 
uint16_t nstrips_
 
HistoSet pedhist_
 
std::vector< int16_t > peds_
 
bool pedsdone_
 
std::vector< float > pedsfl_
 
bool skipped_
 
bool useavgcm_
 
bool usefloatpeds_
 

Additional Inherited Members

- Protected Member Functions inherited from CommissioningTask
const FedChannelConnectionconnection () const
 
DQMStore *const dqm () const
 
const edm::EventSetup *const eventSetup () const
 
const uint32_t & fecKey () const
 
const uint32_t & fedKey () const
 
void updateHistoSet (HistoSet &, const uint32_t &bin, const float &value)
 
void updateHistoSet (CompactHistoSet &, const uint32_t &bin, const short &value)
 
void updateHistoSet (HistoSet &, const uint32_t &bin)
 
void updateHistoSet (CompactHistoSet &, const uint32_t &bin)
 
void updateHistoSet (HistoSet &, const float &value)
 
void updateHistoSet (CompactHistoSet &)
 
void updateHistoSet (HistoSet &)
 

Detailed Description

Definition at line 20 of file PedsFullNoiseTask.h.

Constructor & Destructor Documentation

PedsFullNoiseTask::PedsFullNoiseTask ( DQMStore dqm,
const FedChannelConnection conn,
const edm::ParameterSet pset 
)

Definition at line 17 of file PedsFullNoiseTask.cc.

References fillnoiseprofile_, edm::ParameterSet::getParameter(), LogTrace, sistrip::mlDqmSource_, nadcnoise_, nevpeds_, nskip_, pedsdone_, skipped_, useavgcm_, and usefloatpeds_.

17  :
18  CommissioningTask(dqm, conn, "PedsFullNoiseTask"),
19  nstrips_(256)
20 {
22  << "[PedsFullNoiseTask::" << __func__ << "]"
23  << " Constructing object...";
24  edm::ParameterSet params = pset.getParameter<edm::ParameterSet>("PedsFullNoiseParameters");
25  nskip_ = params.getParameter<int>("NrEvToSkipAtStart");
26  skipped_ = false;
27  nevpeds_ = params.getParameter<int>("NrEvForPeds");
28  pedsdone_ = false;
29  nadcnoise_ = params.getParameter<int>("NrPosBinsNoiseHist");
30  fillnoiseprofile_ = params.getParameter<bool>("FillNoiseProfile");
31  useavgcm_ = params.getParameter<bool>("UseAverageCommonMode");
32  usefloatpeds_ = params.getParameter<bool>("UseFloatPedestals");
33 }
T getParameter(std::string const &) const
static const char mlDqmSource_[]
#define LogTrace(id)
PedsFullNoiseTask::~PedsFullNoiseTask ( )
override

Definition at line 37 of file PedsFullNoiseTask.cc.

References LogTrace, and sistrip::mlDqmSource_.

38 {
40  << "[PedsFullNoiseTask::" << __func__ << "]"
41  << " Destructing object...";
42 }
static const char mlDqmSource_[]
#define LogTrace(id)

Member Function Documentation

void PedsFullNoiseTask::book ( )
overrideprivatevirtual

Reimplemented from CommissioningTask.

Definition at line 46 of file PedsFullNoiseTask.cc.

References CommissioningTask::connection(), CommissioningTask::dqm(), sistrip::EXPERT_HISTO, CommissioningTask::HistoSet::explicitFill_, CommissioningTask::CompactHistoSet::explicitFill_, sistrip::FED_KEY, CommissioningTask::fedKey(), hist2d_, CommissioningTask::HistoSet::histo(), CommissioningTask::CompactHistoSet::histo(), CommissioningTask::HistoSet::isProfile_, sistrip::LLD_CHAN, LogTrace, sistrip::mlDqmSource_, nadcnoise_, sistrip::extrainfo::noise2D_, noisehist_, noiseprof_, sistrip::extrainfo::noiseProfile_, nstrips_, sistrip::extrainfo::pedestals_, pedhist_, sistrip::PEDS_FULL_NOISE, AlCaHLTBitMon_QueryRunRegistry::string, SiStripHistoTitle::title(), CommissioningTask::HistoSet::vNumOfEntries_, CommissioningTask::CompactHistoSet::vNumOfEntries_, CommissioningTask::HistoSet::vSumOfContents_, and CommissioningTask::HistoSet::vSumOfSquares_.

47 {
48 
49  LogTrace(sistrip::mlDqmSource_) << "[PedsFullNoiseTask::" << __func__ << "]";
50 
51  // pedestal profile histo
52  pedhist_.isProfile_ = true;
53  pedhist_.explicitFill_ = false;
54  if (!pedhist_.explicitFill_) {
58  }
62  fedKey(),
64  connection().lldChannel(),
66  pedhist_.histo( dqm()->bookProfile( titleped, titleped,
67  nstrips_, -0.5, nstrips_*1.-0.5,
68  1025, 0., 1025. ) );
69 
70  // Noise profile
71  noiseprof_.isProfile_ = true;
72  noiseprof_.explicitFill_ = false;
77  }
81  fedKey(),
83  connection().lldChannel(),
85  noiseprof_.histo( dqm()->bookProfile( titlenoise, titlenoise,
86  nstrips_, -0.5, nstrips_*1.-0.5,
87  1025, 0., 1025. ) );
88 
89  // noise 2D compact histo
90  noisehist_.explicitFill_ = false;
92  noisehist_.vNumOfEntries_.resize((nstrips_+2)*2*(nadcnoise_+2), 0);
93  }
97  fedKey(),
99  connection().lldChannel(),
101  noisehist_.histo( dqm()->book2S( titlenoise2d, titlenoise2d,
103  nstrips_, -0.5, nstrips_*1.-0.5 ) );
104  hist2d_ = (TH2S *) noisehist_.histo()->getTH2S();
105 
106 }
std::vector< float > vNumOfEntries_
Utility class that holds histogram title.
const std::string & title() const
static const char mlDqmSource_[]
std::vector< float > vSumOfContents_
static const char noise2D_[]
static const char noiseProfile_[]
#define LogTrace(id)
static const char pedestals_[]
DQMStore *const dqm() const
CompactHistoSet noisehist_
void histo(MonitorElement *)
const uint32_t & fedKey() const
std::vector< double > vSumOfSquares_
const FedChannelConnection & connection() const
void PedsFullNoiseTask::fill ( const SiStripEventSummary summary,
const edm::DetSet< SiStripRawDigi > &  digis 
)
overrideprivatevirtual

Reimplemented from CommissioningTask.

Definition at line 110 of file PedsFullNoiseTask.cc.

References ecalMGPA::adc(), edm::DetSet< T >::data, SiStripEventSummary::event(), fillnoiseprofile_, hist2d_, diffTreeTool::index, LogTrace, sistrip::mlDqmSource_, nadcnoise_, pileupCalc::nbins, nevpeds_, noisehist_, noiseprof_, nskip_, nstrips_, pedhist_, peds_, pedsdone_, pedsfl_, skipped_, CommissioningTask::updateHistoSet(), useavgcm_, usefloatpeds_, CommissioningTask::HistoSet::vNumOfEntries_, and CommissioningTask::HistoSet::vSumOfContents_.

112 {
113 
114  // Check number of digis
115  uint16_t nbins = digis.data.size();
116  if (nbins != nstrips_) {
118  << "[PedsFullNoiseTask::" << __func__ << "]"
119  << " " << nstrips_ << " digis expected, but got " << nbins << ". Skipping.";
120  return;
121  }
122 
123  // get the event number of the first event, not necessarily 1 (parallel processing on FUs)
124  static int32_t firstev = summary.event();
125 
126  // skipping events
127  if (!skipped_) {
128  if (static_cast<int32_t>(summary.event()) - firstev < nskip_) {
129  return;
130  } else { // when all events are skipped
131  skipped_ = true;
132  if (nskip_ > 0) LogTrace(sistrip::mlDqmSource_) << "[PedsFullNoiseTask::" << __func__ << "]"
133  << " Done skipping events. Now starting pedestals.";
134  }
135  }
136 
137  // determine pedestals - decoupled from noise determination
138  if (!pedsdone_) {
139  if (static_cast<int32_t>(summary.event()) - firstev < nskip_ + nevpeds_) {
140  // estimate the pedestals
141  for ( uint16_t istrip = 0; istrip < nstrips_; ++istrip ) {
142  updateHistoSet( pedhist_, istrip, digis.data[istrip].adc() );
143  }
144  return;
145  } else { // when pedestals are done
146  pedsdone_ = true;
147  // cache the pedestal values for use in the 2D noise estimation
148  peds_.clear(); pedsfl_.clear();
149  for ( uint16_t iapv = 0; iapv < 2; ++iapv ) {
150  for ( uint16_t ibin = 0; ibin < 128; ++ibin ) {
151  uint16_t istrip = (iapv*128)+ibin;
152  if (usefloatpeds_) {
153  pedsfl_.push_back(1.*pedhist_.vSumOfContents_.at(istrip)/pedhist_.vNumOfEntries_.at(istrip));
154  } else {
155  peds_.push_back(static_cast<int16_t>(1.*pedhist_.vSumOfContents_.at(istrip)/pedhist_.vNumOfEntries_.at(istrip)));
156  }
157  }
158  }
159  LogTrace(sistrip::mlDqmSource_) << "[PedsFullNoiseTask::" << __func__ << "]"
160  << " Rough pedestals done. Now starting noise measurements.";
161  }
162  }
163 
164  // fill (or not) the old-style niose profile
165  if (fillnoiseprofile_) {
166  // Calc common mode for both APVs
167  std::vector<int32_t> cm; cm.resize(2, 0);
168  std::vector<uint16_t> adc;
169  for ( uint16_t iapv = 0; iapv < 2; iapv++ ) {
170  adc.clear(); adc.reserve(128);
171  for ( uint16_t ibin = 0; ibin < 128; ibin++ ) {
172  if ( (iapv*128)+ibin < nbins ) {
173  adc.push_back( digis.data.at((iapv*128)+ibin).adc() );
174  }
175  }
176  sort( adc.begin(), adc.end() );
177  // take median as common mode
178  uint16_t index = adc.size()%2 ? adc.size()/2 : adc.size()/2-1;
179  cm[iapv] = static_cast<int16_t>(adc[index]);
180  }
181  // 1D noise profile - see also further processing in the update() method
182  for ( uint16_t istrip = 0; istrip < nstrips_; ++istrip ) {
183  // calculate the noise in the old way, by subtracting the common mode, but without pedestal subtraction
184  int16_t noiseval = static_cast<int16_t>(digis.data.at(istrip).adc()) - cm[istrip/128];
185  updateHistoSet( noiseprof_, istrip, noiseval );
186  }
187  }
188 
189  // 2D noise histogram
190  std::vector<int16_t> noisevals, noisevalssorted;
191  std::vector<float> noisevalsfl, noisevalssortedfl;
192  for ( uint16_t iapv = 0; iapv < 2; ++iapv ) {
193  float totadc = 0;
194  noisevals.clear(); noisevalsfl.clear();
195  noisevalssorted.clear(); noisevalssortedfl.clear();
196  for ( uint16_t ibin = 0; ibin < 128; ++ibin ) {
197  uint16_t istrip = (iapv*128)+ibin;
198  // calculate the noise after subtracting the pedestal
199  if (usefloatpeds_) { // if float pedestals -> before FED processing
200  noisevalsfl.push_back( static_cast<float>(digis.data.at(istrip).adc()) - pedsfl_.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 += noisevalsfl[ibin];
204  } else { // if median CM -> after FED processing
205  noisevalssortedfl.push_back( noisevalsfl[ibin] );
206  }
207  } else { // if integer pedestals -> after FED processing
208  noisevals.push_back( static_cast<int16_t>(digis.data.at(istrip).adc()) - peds_.at(istrip) );
209  // 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
210  if (useavgcm_) { // if average CM -> before FED processing
211  totadc += noisevals[ibin];
212  } else { // if median CM -> after FED processing
213  noisevalssorted.push_back( noisevals[ibin] );
214  }
215  }
216  }
217  // calculate the common mode shift to apply
218  float cmshift = 0;
219  if (useavgcm_) { // if average CM -> before FED processing
220  if (usefloatpeds_) { // if float pedestals -> before FED processing
221  cmshift = totadc/128;
222  } else { // if integer pedestals -> after FED processing
223  cmshift = static_cast<int16_t>(totadc/128);
224  }
225  } else { // if median CM -> after FED processing
226  if (usefloatpeds_) { // if float pedestals -> before FED processing
227  // get the median common mode
228  sort( noisevalssortedfl.begin(), noisevalssortedfl.end() );
229  uint16_t index = noisevalssortedfl.size()%2 ? noisevalssortedfl.size()/2 : noisevalssortedfl.size()/2-1;
230  cmshift = noisevalssortedfl[index];
231  } else { // if integer pedestals -> after FED processing
232  // get the median common mode
233  sort( noisevalssorted.begin(), noisevalssorted.end() );
234  uint16_t index = noisevalssorted.size()%2 ? noisevalssorted.size()/2 : noisevalssorted.size()/2-1;
235  cmshift = noisevalssorted[index];
236  }
237  }
238  // now loop again to calculate the CM+pedestal subtracted noise values
239  for ( uint16_t ibin = 0; ibin < 128; ++ibin ) {
240  uint16_t istrip = (iapv*128)+ibin;
241  // subtract the remaining common mode after subtraction of the rough pedestals
242  float noiseval = (usefloatpeds_ ? noisevalsfl[ibin] : noisevals[ibin]) - cmshift;
243  // retrieve the linear binnr through the histogram
244  uint32_t binnr = hist2d_->GetBin(static_cast<int>(noiseval+nadcnoise_), istrip+1);
245  // store the noise value in the 2D histo
246  updateHistoSet( noisehist_, binnr ); // no value, so weight 1
247  }
248  }
249 
250 }
int adc(sample_type sample)
get the ADC sample (12 bits)
std::vector< float > vNumOfEntries_
static const char mlDqmSource_[]
std::vector< float > vSumOfContents_
void updateHistoSet(HistoSet &, const uint32_t &bin, const float &value)
const uint32_t & event() const
#define LogTrace(id)
CompactHistoSet noisehist_
std::vector< int16_t > peds_
collection_type data
Definition: DetSet.h:78
std::vector< float > pedsfl_
void PedsFullNoiseTask::update ( )
overrideprivatevirtual

Reimplemented from CommissioningTask.

Definition at line 254 of file PedsFullNoiseTask.cc.

References relativeConstraints::error, extract(), fillnoiseprofile_, trackerHits::histo, CommissioningTask::HistoSet::histo(), cuy::ii, SiStripPI::mean, noisehist_, noiseprof_, pedhist_, UpdateTProfile::setBinContent(), mathSSE::sqrt(), CommissioningTask::updateHistoSet(), CommissioningTask::HistoSet::vNumOfEntries_, CommissioningTask::HistoSet::vSumOfContents_, and CommissioningTask::HistoSet::vSumOfSquares_.

Referenced by progressbar.ProgressBar::__next__(), MatrixUtil.Matrix::__setitem__(), MatrixUtil.Steps::__setitem__(), Vispa.Gui.VispaWidget.VispaWidget::autosize(), Vispa.Views.LineDecayView.LineDecayContainer::createObject(), Vispa.Views.LineDecayView.LineDecayContainer::deselectAllObjects(), Vispa.Gui.VispaWidgetOwner.VispaWidgetOwner::deselectAllWidgets(), Vispa.Gui.VispaWidget.VispaWidget::enableAutosizing(), dqm-mbProfile.Profile::finish(), progressbar.ProgressBar::finish(), Vispa.Gui.MenuWidget.MenuWidget::leaveEvent(), Vispa.Gui.VispaWidgetOwner.VispaWidgetOwner::mouseMoveEvent(), Vispa.Gui.MenuWidget.MenuWidget::mouseMoveEvent(), Vispa.Views.LineDecayView.LineDecayContainer::mouseMoveEvent(), Vispa.Gui.VispaWidgetOwner.VispaWidgetOwner::mouseReleaseEvent(), Vispa.Views.LineDecayView.LineDecayContainer::objectMoved(), MatrixUtil.Steps::overwrite(), Vispa.Views.LineDecayView.LineDecayContainer::removeObject(), Vispa.Gui.ConnectableWidget.ConnectableWidget::removePorts(), Vispa.Gui.FindDialog.FindDialog::reset(), Vispa.Gui.PortConnection.PointToPointConnection::select(), Vispa.Gui.VispaWidget.VispaWidget::select(), Vispa.Views.LineDecayView.LineDecayContainer::select(), Vispa.Gui.VispaWidget.VispaWidget::setText(), Vispa.Gui.VispaWidget.VispaWidget::setTitle(), Vispa.Gui.ZoomableWidget.ZoomableWidget::setZoom(), Vispa.Views.LineDecayView.LineDecayContainer::setZoom(), and Vispa.Gui.PortConnection.PointToPointConnection::updateConnection().

255 {
256 
257  // pedestals
259 
260  if (fillnoiseprofile_) {
261  // noise profile (does not use HistoSet directly, as want to plot noise as "contents", not "error")
263  for ( uint16_t ii = 0; ii < noiseprof_.vNumOfEntries_.size(); ++ii ) {
264  float mean = 0.;
265  float spread = 0.;
266  float entries = noiseprof_.vNumOfEntries_[ii];
267  if ( entries > 0. ) {
268  mean = noiseprof_.vSumOfContents_[ii] / entries;
269  spread = sqrt( noiseprof_.vSumOfSquares_[ii] / entries - mean * mean ); // nice way to calculate std dev: Sum (x-<x>)^2 / N
270  }
271  float error = spread / sqrt(entries); // uncertainty on std.dev. when no uncertainty on mean
272  UpdateTProfile::setBinContent( histo, ii+1, entries, spread, error );
273  }
274  }
275 
276  // noise 2D histo
278 
279 }
std::vector< float > vNumOfEntries_
std::vector< float > vSumOfContents_
void updateHistoSet(HistoSet &, const uint32_t &bin, const float &value)
T sqrt(T t)
Definition: SSEVec.h:18
static void setBinContent(TProfile *const profile, const uint32_t &bin, const double &entries, const double &mean, const double &spread)
ii
Definition: cuy.py:588
int extract(std::vector< int > *output, const std::string &dati)
CompactHistoSet noisehist_
void histo(MonitorElement *)
std::vector< double > vSumOfSquares_

Member Data Documentation

bool PedsFullNoiseTask::fillnoiseprofile_
private

Definition at line 53 of file PedsFullNoiseTask.h.

Referenced by fill(), PedsFullNoiseTask(), and update().

TH2S* PedsFullNoiseTask::hist2d_
private

Definition at line 37 of file PedsFullNoiseTask.h.

Referenced by book(), and fill().

uint16_t PedsFullNoiseTask::nadcnoise_
private

Definition at line 49 of file PedsFullNoiseTask.h.

Referenced by book(), fill(), and PedsFullNoiseTask().

uint16_t PedsFullNoiseTask::nevpeds_
private

Definition at line 47 of file PedsFullNoiseTask.h.

Referenced by fill(), and PedsFullNoiseTask().

CompactHistoSet PedsFullNoiseTask::noisehist_
private

Definition at line 36 of file PedsFullNoiseTask.h.

Referenced by book(), fill(), and update().

HistoSet PedsFullNoiseTask::noiseprof_
private

Definition at line 35 of file PedsFullNoiseTask.h.

Referenced by book(), fill(), and update().

uint16_t PedsFullNoiseTask::nskip_
private

Definition at line 43 of file PedsFullNoiseTask.h.

Referenced by fill(), and PedsFullNoiseTask().

uint16_t PedsFullNoiseTask::nstrips_
private

Definition at line 51 of file PedsFullNoiseTask.h.

Referenced by book(), and fill().

HistoSet PedsFullNoiseTask::pedhist_
private

Definition at line 35 of file PedsFullNoiseTask.h.

Referenced by book(), fill(), and update().

std::vector<int16_t> PedsFullNoiseTask::peds_
private

Definition at line 38 of file PedsFullNoiseTask.h.

Referenced by fill().

bool PedsFullNoiseTask::pedsdone_
private

Definition at line 45 of file PedsFullNoiseTask.h.

Referenced by fill(), and PedsFullNoiseTask().

std::vector<float> PedsFullNoiseTask::pedsfl_
private

Definition at line 39 of file PedsFullNoiseTask.h.

Referenced by fill().

bool PedsFullNoiseTask::skipped_
private

Definition at line 41 of file PedsFullNoiseTask.h.

Referenced by fill(), and PedsFullNoiseTask().

bool PedsFullNoiseTask::useavgcm_
private

Definition at line 55 of file PedsFullNoiseTask.h.

Referenced by fill(), and PedsFullNoiseTask().

bool PedsFullNoiseTask::usefloatpeds_
private

Definition at line 57 of file PedsFullNoiseTask.h.

Referenced by fill(), and PedsFullNoiseTask().