Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-03-17 10:56:31

0001 #include "DQM/SiStripCommissioningSources/interface/DaqScopeModeTask.h"
0002 #include "DataFormats/SiStripCommon/interface/SiStripConstants.h"
0003 #include "DataFormats/SiStripCommon/interface/SiStripHistoTitle.h"
0004 #include "DQMServices/Core/interface/DQMStore.h"
0005 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0006 #include "DQM/SiStripCommon/interface/ExtractTObject.h"
0007 #include "DQM/SiStripCommon/interface/UpdateTProfile.h"
0008 
0009 using namespace sistrip;
0010 
0011 // -----------------------------------------------------------------------------
0012 //
0013 DaqScopeModeTask::DaqScopeModeTask(DQMStore* dqm, const FedChannelConnection& conn, const edm::ParameterSet& pset)
0014     : CommissioningTask(dqm, conn, "DaqScopeModeTask"),
0015       scopeFrame_(),
0016       nBins_(256),     //@@ number of strips per FED channel
0017       nBinsSpy_(298),  //@@ in case of spy events includes header and trailing tick
0018       parameters_(pset) {}
0019 
0020 // -----------------------------------------------------------------------------
0021 //
0022 DaqScopeModeTask::~DaqScopeModeTask() {}
0023 
0024 // -----------------------------------------------------------------------------
0025 //
0026 void DaqScopeModeTask::book() {
0027   LogTrace(mlDqmSource_) << "[CommissioningTask::" << __func__ << "]";
0028 
0029   std::string title;
0030   if (not(parameters_.existsAs<bool>("isSpy") and parameters_.getParameter<bool>("isSpy"))) {
0031     //// Scope mode histograms
0032     title = SiStripHistoTitle(sistrip::EXPERT_HISTO,
0033                               sistrip::DAQ_SCOPE_MODE,
0034                               sistrip::FED_KEY,
0035                               fedKey(),
0036                               sistrip::LLD_CHAN,
0037                               connection().lldChannel(),
0038                               sistrip::extrainfo::scopeModeFrame_)
0039                 .title();
0040 
0041     scopeFrame_.histo(dqm()->book1D(title, title, nBins_, -0.5, nBins_ - 0.5));
0042 
0043     scopeFrame_.vNumOfEntries_.resize(nBins_, 0);
0044     scopeFrame_.vSumOfContents_.resize(nBins_, 0);
0045     scopeFrame_.vSumOfSquares_.resize(nBins_, 0);
0046     scopeFrame_.isProfile_ = false;
0047   } else {  // use spy run to measure pedestal and tick height
0048 
0049     std::string extra_info;
0050     peds_.resize(2);
0051     extra_info = sistrip::extrainfo::pedestals_;
0052 
0053     peds_[0].isProfile_ = true;
0054     title = SiStripHistoTitle(sistrip::EXPERT_HISTO,
0055                               sistrip::DAQ_SCOPE_MODE,
0056                               sistrip::FED_KEY,
0057                               fedKey(),
0058                               sistrip::LLD_CHAN,
0059                               connection().lldChannel(),
0060                               extra_info)
0061                 .title();
0062 
0063     peds_[0].histo(dqm()->bookProfile(title, title, nBins_, -0.5, nBins_ * 1. - 0.5, 1025, 0., 1025.));
0064 
0065     peds_[0].vNumOfEntries_.resize(nBins_, 0);
0066     peds_[0].vSumOfContents_.resize(nBins_, 0);
0067     peds_[0].vSumOfSquares_.resize(nBins_, 0);
0068 
0069     // Noise histogram
0070     extra_info = sistrip::extrainfo::noise_;
0071     peds_[1].isProfile_ = true;
0072 
0073     title = SiStripHistoTitle(sistrip::EXPERT_HISTO,
0074                               sistrip::DAQ_SCOPE_MODE,
0075                               sistrip::FED_KEY,
0076                               fedKey(),
0077                               sistrip::LLD_CHAN,
0078                               connection().lldChannel(),
0079                               extra_info)
0080                 .title();
0081 
0082     peds_[1].histo(dqm()->bookProfile(title, title, nBins_, -0.5, nBins_ * 1. - 0.5, 1025, 0., 1025.));
0083 
0084     peds_[1].vNumOfEntries_.resize(nBins_, 0);
0085     peds_[1].vSumOfContents_.resize(nBins_, 0);
0086     peds_[1].vSumOfSquares_.resize(nBins_, 0);
0087 
0088     // Common mode histograms
0089     cm_.resize(2);
0090     int nbins = 1024;
0091 
0092     for (uint16_t iapv = 0; iapv < 2; iapv++) {
0093       title = SiStripHistoTitle(sistrip::EXPERT_HISTO,
0094                                 sistrip::DAQ_SCOPE_MODE,
0095                                 sistrip::FED_KEY,
0096                                 fedKey(),
0097                                 sistrip::APV,
0098                                 connection().i2cAddr(iapv),
0099                                 sistrip::extrainfo::commonMode_)
0100                   .title();
0101 
0102       cm_[iapv].histo(dqm()->book1D(title, title, nbins, -0.5, nbins * 1. - 0.5));
0103       cm_[iapv].isProfile_ = false;
0104 
0105       cm_[iapv].vNumOfEntries_.resize(nbins, 0);
0106       cm_[iapv].vNumOfEntries_.resize(nbins, 0);
0107     }
0108 
0109     // high and low header histograms
0110     title = SiStripHistoTitle(sistrip::EXPERT_HISTO,
0111                               sistrip::DAQ_SCOPE_MODE,
0112                               sistrip::FED_KEY,
0113                               fedKey(),
0114                               sistrip::LLD_CHAN,
0115                               connection().lldChannel(),
0116                               sistrip::extrainfo::scopeModeHeaderLow_)
0117                 .title();
0118 
0119     lowHeader_.histo(dqm()->book1D(title, title, nbins, -0.5, 1024 * 1. - 0.5));
0120     lowHeader_.isProfile_ = false;
0121     lowHeader_.vNumOfEntries_.resize(nbins, 0);
0122 
0123     title = SiStripHistoTitle(sistrip::EXPERT_HISTO,
0124                               sistrip::DAQ_SCOPE_MODE,
0125                               sistrip::FED_KEY,
0126                               fedKey(),
0127                               sistrip::LLD_CHAN,
0128                               connection().lldChannel(),
0129                               sistrip::extrainfo::scopeModeHeaderHigh_)
0130                 .title();
0131 
0132     highHeader_.histo(dqm()->book1D(title, title, nbins, -0.5, 1024 * 1. - 0.5));
0133     highHeader_.isProfile_ = false;
0134     highHeader_.vNumOfEntries_.resize(nbins, 0);
0135 
0136     //// Scope mode histograms
0137     title = SiStripHistoTitle(sistrip::EXPERT_HISTO,
0138                               sistrip::DAQ_SCOPE_MODE,
0139                               sistrip::FED_KEY,
0140                               fedKey(),
0141                               sistrip::LLD_CHAN,
0142                               connection().lldChannel(),
0143                               sistrip::extrainfo::scopeModeFrame_)
0144                 .title();
0145 
0146     scopeFrame_.histo(dqm()->bookProfile(title, title, nBinsSpy_, -0.5, nBinsSpy_ * 1. - 0.5, 1025, 0., 1025.));
0147 
0148     scopeFrame_.vNumOfEntries_.resize(nBinsSpy_, 0);
0149     scopeFrame_.vSumOfContents_.resize(nBinsSpy_, 0);
0150     scopeFrame_.vSumOfSquares_.resize(nBinsSpy_, 0);
0151     scopeFrame_.isProfile_ = true;
0152   }
0153 }
0154 
0155 // -----------------------------------------------------------------------------
0156 //
0157 void DaqScopeModeTask::fill(const SiStripEventSummary& summary, const edm::DetSet<SiStripRawDigi>& digis) {
0158   // Only fill every 'N' events
0159   if (not(parameters_.existsAs<bool>("isSpy") and parameters_.getParameter<bool>("isSpy"))) {
0160     if (!updateFreq() || fillCntr() % updateFreq()) {
0161       return;
0162     }
0163   }
0164 
0165   if (digis.data.size() != nBins_) {  //@@ check scope mode length?
0166     edm::LogWarning(mlDqmSource_) << "[DaqScopeModeTask::" << __func__ << "]"
0167                                   << " Unexpected number of digis (" << digis.data.size()
0168                                   << ") wrt number of histogram bins (" << nBins_ << ")!";
0169   }
0170 
0171   if (not(parameters_.existsAs<bool>("isSpy") and parameters_.getParameter<bool>("isSpy"))) {
0172     uint16_t bins = digis.data.size() < nBins_ ? digis.data.size() : nBins_;
0173     for (uint16_t ibin = 0; ibin < bins; ibin++) {
0174       updateHistoSet(scopeFrame_, ibin, digis.data[ibin].adc());
0175     }
0176   } else {
0177     // fill the pedestal histograms as done in the pedestal task
0178     // Check number of digis
0179     uint16_t nbins = peds_[0].vNumOfEntries_.size();
0180     if (digis.data.size() < nbins) {
0181       nbins = digis.data.size();
0182     }
0183 
0184     uint16_t napvs = nbins / 128;
0185     std::vector<uint32_t> cm;
0186     cm.resize(napvs, 0);
0187 
0188     // Calc common mode for both APVs
0189     std::vector<uint16_t> adc;
0190     for (uint16_t iapv = 0; iapv < napvs; iapv++) {
0191       adc.clear();
0192       adc.reserve(128);
0193 
0194       for (uint16_t ibin = 0; ibin < 128; ibin++) {
0195         if ((iapv * 128) + ibin < nbins) {
0196           adc.push_back(digis.data[(iapv * 128) + ibin].adc());
0197         }
0198       }
0199 
0200       sort(adc.begin(), adc.end());
0201       uint16_t index = adc.size() % 2 ? adc.size() / 2 : adc.size() / 2 - 1;
0202       if (!adc.empty()) {
0203         cm[iapv] = static_cast<uint32_t>(adc[index]);
0204       }
0205     }
0206 
0207     for (uint16_t ibin = 0; ibin < nbins; ibin++) {
0208       float digiVal = digis.data[ibin].adc();
0209       updateHistoSet(peds_[0], ibin, digiVal);  // peds and raw noise
0210       float diff = digiVal - static_cast<float>(cm[ibin / 128]);
0211       updateHistoSet(peds_[1], ibin, diff);  // residuals and real noise
0212     }
0213 
0214     if (cm.size() < cm_.size()) {
0215       edm::LogWarning(mlDqmSource_) << "[PedestalsTask::" << __func__ << "]"
0216                                     << " Fewer CM values than expected: " << cm.size();
0217     }
0218 
0219     updateHistoSet(cm_[0], cm[0]);
0220     updateHistoSet(cm_[1], cm[1]);
0221   }
0222 }
0223 
0224 // -----------------------------------------------------------------------------
0225 //
0226 void DaqScopeModeTask::fill(const SiStripEventSummary& summary,
0227                             const edm::DetSet<SiStripRawDigi>& digis,
0228                             const edm::DetSet<SiStripRawDigi>& digisAlt) {
0229   // Only fill every 'N' events
0230   if (not(parameters_.existsAs<bool>("isSpy") and parameters_.getParameter<bool>("isSpy"))) {
0231     if (!updateFreq() || fillCntr() % updateFreq()) {
0232       return;
0233     }
0234   }
0235 
0236   if (digis.data.size() != nBins_) {  //@@ check scope mode length?
0237     edm::LogWarning(mlDqmSource_) << "[DaqScopeModeTask::" << __func__ << "]"
0238                                   << " Unexpected number of digis (" << digis.data.size()
0239                                   << ") wrt number of histogram bins (" << nBins_ << ")!";
0240   }
0241 
0242   if (not(parameters_.existsAs<bool>("isSpy") and parameters_.getParameter<bool>("isSpy"))) {
0243     uint16_t bins = digis.data.size() < nBins_ ? digis.data.size() : nBins_;
0244     for (uint16_t ibin = 0; ibin < bins; ibin++) {
0245       updateHistoSet(scopeFrame_, ibin, digis.data[ibin].adc());
0246     }
0247   } else {
0248     // fill the pedestal histograms as done in the pedestal task
0249     // Check number of digis
0250     uint16_t nbins = peds_[0].vNumOfEntries_.size();
0251     if (digis.data.size() < nbins) {
0252       nbins = digis.data.size();
0253     }
0254 
0255     uint16_t napvs = nbins / 128;
0256     std::vector<uint32_t> cm;
0257     cm.resize(napvs, 0);
0258 
0259     // Calc common mode for both APVs
0260     std::vector<uint16_t> adc;
0261     for (uint16_t iapv = 0; iapv < napvs; iapv++) {
0262       adc.clear();
0263       adc.reserve(128);
0264       for (uint16_t ibin = 0; ibin < 128; ibin++) {
0265         if ((iapv * 128) + ibin < nbins) {
0266           adc.push_back(digis.data[(iapv * 128) + ibin].adc());
0267         }
0268       }
0269       sort(adc.begin(), adc.end());
0270       uint16_t index = adc.size() % 2 ? adc.size() / 2 : adc.size() / 2 - 1;
0271       if (!adc.empty()) {
0272         cm[iapv] = static_cast<uint32_t>(adc[index]);
0273       }
0274     }
0275 
0276     /// Calculate pedestal
0277     for (uint16_t ibin = 0; ibin < nbins; ibin++) {
0278       float digiVal = digis.data[ibin].adc();
0279       updateHistoSet(peds_[0], ibin, digiVal);  // peds and raw noise
0280       float diff = digiVal - static_cast<float>(cm[ibin / 128]);
0281       updateHistoSet(peds_[1], ibin, diff);  // residuals and real noise
0282     }
0283 
0284     if (cm.size() < cm_.size()) {
0285       edm::LogWarning(mlDqmSource_) << "[PedestalsTask::" << __func__ << "]"
0286                                     << " Fewer CM values than expected: " << cm.size();
0287     }
0288 
0289     updateHistoSet(cm_[0], cm[0]);
0290     updateHistoSet(cm_[1], cm[1]);
0291 
0292     uint16_t bins = digisAlt.data.size() < nBinsSpy_ ? digisAlt.data.size() : nBinsSpy_;
0293     for (uint16_t ibin = 0; ibin < bins; ibin++) {
0294       updateHistoSet(scopeFrame_, ibin, digisAlt.data[ibin].adc());
0295     }
0296 
0297     // Header low and high for both APVs
0298     std::vector<uint32_t> adcHeader_high;
0299     std::vector<uint32_t> adcHeader_low;
0300 
0301     float threshold_high = (digisAlt.data[286].adc() + digisAlt.data[287].adc()) / 4;
0302     float threshold_low = 100;
0303     int minNumberForHeader = 4;
0304     bool goodHeaderFound = false;
0305     int nConsecutiveHigh = 0;
0306     adcHeader_high.clear();
0307     adcHeader_high.reserve(30);
0308     adcHeader_low.clear();
0309     adcHeader_low.reserve(30);
0310 
0311     for (uint16_t ibin = 6; ibin < 11; ibin++) {
0312       if (digisAlt.data[ibin].adc() > threshold_high) {
0313         nConsecutiveHigh++;
0314       }
0315     }
0316 
0317     if (nConsecutiveHigh > minNumberForHeader)
0318       goodHeaderFound = true;  // if nConsecutiveHigh > 4 --> good header found
0319     if (goodHeaderFound == false)
0320       return;
0321     for (uint16_t ibin = 0; ibin < 30; ibin++) {
0322       if (digisAlt.data[ibin].adc() > threshold_high &&
0323           goodHeaderFound) {  // save of samples above avg(trailing ticks)/4
0324         adcHeader_high.push_back(digisAlt.data[ibin].adc());
0325       }
0326       if (digisAlt.data[ibin].adc() < threshold_low && goodHeaderFound) {
0327         adcHeader_low.push_back(digisAlt.data[ibin].adc());
0328       }
0329     }
0330     if (adcHeader_low.empty() || adcHeader_high.empty()) {
0331       return;
0332     }
0333     for (uint16_t i = 0; i < adcHeader_low.size(); i++) {
0334       updateHistoSet(lowHeader_, adcHeader_low[i]);
0335     }
0336     for (uint16_t i = 0; i < adcHeader_high.size(); i++) {
0337       updateHistoSet(highHeader_, adcHeader_high[i]);
0338     }
0339   }
0340 }
0341 
0342 // -----------------------------------------------------------------------------
0343 //
0344 void DaqScopeModeTask::fill(const SiStripEventSummary& summary,
0345                             const edm::DetSet<SiStripRawDigi>& digis,
0346                             const edm::DetSet<SiStripRawDigi>& digisAlt,
0347                             const std::vector<uint16_t>& stripOnCluster) {
0348   // Only fill every 'N' events
0349   if (not(parameters_.existsAs<bool>("isSpy") and parameters_.getParameter<bool>("isSpy"))) {
0350     if (!updateFreq() || fillCntr() % updateFreq()) {
0351       return;
0352     }
0353   }
0354 
0355   if (digis.data.size() != nBins_) {  //@@ check scope mode length?
0356     edm::LogWarning(mlDqmSource_) << "[DaqScopeModeTask::" << __func__ << "]"
0357                                   << " Unexpected number of digis (" << digis.data.size()
0358                                   << ") wrt number of histogram bins (" << nBins_ << ")!";
0359   }
0360 
0361   if (not(parameters_.existsAs<bool>("isSpy") and parameters_.getParameter<bool>("isSpy"))) {
0362     uint16_t bins = digis.data.size() < nBins_ ? digis.data.size() : nBins_;
0363     for (uint16_t ibin = 0; ibin < bins; ibin++) {
0364       updateHistoSet(scopeFrame_, ibin, digis.data[ibin].adc());
0365     }
0366   } else {
0367     // fill the pedestal histograms as done in the pedestal task
0368     uint16_t nbins = peds_[0].vNumOfEntries_.size();
0369     if (digis.data.size() < nbins) {
0370       nbins = digis.data.size();
0371     }
0372     uint16_t napvs = nbins / 128;
0373     std::vector<uint32_t> cm;
0374     cm.resize(napvs, 0);
0375     // Calc common mode for both APVs
0376     std::vector<uint16_t> adc;
0377     for (uint16_t iapv = 0; iapv < napvs; iapv++) {
0378       adc.clear();
0379       adc.reserve(128);
0380       for (uint16_t ibin = 0; ibin < 128; ibin++) {
0381         if ((iapv * 128) + ibin < nbins) {
0382           if (std::find(stripOnCluster.begin(), stripOnCluster.end(), (iapv * 128) + ibin) ==
0383               stripOnCluster.end())  // if not found, strip is good
0384             adc.push_back(digis.data[(iapv * 128) + ibin].adc());
0385         }
0386       }
0387       sort(adc.begin(), adc.end());
0388       uint16_t index = adc.size() % 2 ? adc.size() / 2 : adc.size() / 2 - 1;
0389       if (!adc.empty()) {
0390         cm[iapv] = static_cast<uint32_t>(adc[index]);
0391       }
0392     }
0393 
0394     /// Calculate pedestal
0395     for (uint16_t ibin = 0; ibin < nbins; ibin++) {
0396       if (std::find(stripOnCluster.begin(), stripOnCluster.end(), ibin) != stripOnCluster.end()) {
0397         continue;
0398       }
0399       float digiVal = digis.data[ibin].adc();
0400       updateHistoSet(peds_[0], ibin, digiVal);  // peds and raw noise
0401       float diff = digiVal - static_cast<float>(cm[ibin / 128]);
0402       updateHistoSet(peds_[1], ibin, diff);  // residuals and real noise
0403     }
0404 
0405     if (cm.size() < cm_.size()) {
0406       edm::LogWarning(mlDqmSource_) << "[PedestalsTask::" << __func__ << "]"
0407                                     << " Fewer CM values than expected: " << cm.size();
0408     }
0409 
0410     updateHistoSet(cm_[0], cm[0]);
0411     updateHistoSet(cm_[1], cm[1]);
0412 
0413     uint16_t bins = digisAlt.data.size() < nBinsSpy_ ? digisAlt.data.size() : nBinsSpy_;
0414     for (uint16_t ibin = 0; ibin < bins; ibin++) {
0415       updateHistoSet(scopeFrame_, ibin, digisAlt.data[ibin].adc());
0416     }
0417     // Header low and high for both APVs
0418     std::vector<uint32_t> adcHeader_high;
0419     std::vector<uint32_t> adcHeader_low;
0420 
0421     float threshold_high = (digisAlt.data[286].adc() + digisAlt.data[287].adc()) / 4;
0422     float threshold_low = 120;
0423     int minNumberForHeader = 4;
0424     bool goodHeaderFound = false;
0425     int nConsecutiveHigh = 0;
0426     adcHeader_high.clear();
0427     adcHeader_high.reserve(30);
0428     adcHeader_low.clear();
0429     adcHeader_low.reserve(30);
0430 
0431     for (uint16_t ibin = 6; ibin < 11; ibin++) {
0432       if (digisAlt.data[ibin].adc() > threshold_high) {
0433         nConsecutiveHigh++;
0434       }
0435     }
0436 
0437     if (nConsecutiveHigh > minNumberForHeader)
0438       goodHeaderFound = true;  // if nConsecutiveHigh > 4 --> good header found
0439     if (goodHeaderFound == false)
0440       return;
0441     for (uint16_t ibin = 0; ibin < 30; ibin++) {
0442       if (digisAlt.data[ibin].adc() > threshold_high &&
0443           goodHeaderFound) {  // save of samples above avg(trailing ticks)/4
0444         adcHeader_high.push_back(digisAlt.data[ibin].adc());
0445       }
0446       if (digisAlt.data[ibin].adc() < threshold_low && goodHeaderFound) {
0447         adcHeader_low.push_back(digisAlt.data[ibin].adc());
0448       }
0449     }
0450     if (adcHeader_low.empty() || adcHeader_high.empty()) {
0451       return;
0452     }
0453     for (uint16_t i = 0; i < adcHeader_low.size(); i++) {
0454       updateHistoSet(lowHeader_, adcHeader_low[i]);
0455     }
0456     for (uint16_t i = 0; i < adcHeader_high.size(); i++) {
0457       updateHistoSet(highHeader_, adcHeader_high[i]);
0458     }
0459   }
0460 }
0461 
0462 // -----------------------------------------------------------------------------
0463 //
0464 void DaqScopeModeTask::update() {
0465   if (not(parameters_.existsAs<bool>("isSpy") and parameters_.getParameter<bool>("isSpy")))
0466     updateHistoSet(scopeFrame_);
0467   else {
0468     updateHistoSet(peds_[0]);
0469     TProfile* histo = ExtractTObject<TProfile>().extract(peds_[1].histo());
0470     for (uint16_t ii = 0; ii < peds_[1].vNumOfEntries_.size(); ++ii) {
0471       float mean = 0.;
0472       float spread = 0.;
0473       float entries = peds_[1].vNumOfEntries_[ii];
0474       if (entries > 0.) {
0475         mean = peds_[1].vSumOfContents_[ii] / entries;
0476         spread = sqrt(peds_[1].vSumOfSquares_[ii] / entries - mean * mean);
0477       }
0478 
0479       float noise = spread;
0480       float error = 0;  // sqrt(entries) / entries;
0481       UpdateTProfile::setBinContent(histo, ii + 1, entries, noise, error);
0482     }
0483 
0484     updateHistoSet(cm_[0]);
0485     updateHistoSet(cm_[1]);
0486 
0487     updateHistoSet(scopeFrame_);
0488     updateHistoSet(lowHeader_);
0489     updateHistoSet(highHeader_);
0490   }
0491 }