Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:08:30

0001 #include "DQM/SiStripCommissioningAnalysis/interface/VpspScanAlgorithm.h"
0002 #include "CondFormats/SiStripObjects/interface/VpspScanAnalysis.h"
0003 #include "DataFormats/SiStripCommon/interface/SiStripHistoTitle.h"
0004 #include "DataFormats/SiStripCommon/interface/SiStripEnumsAndStrings.h"
0005 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0006 #include "TProfile.h"
0007 #include "TH1.h"
0008 #include <iostream>
0009 #include <sstream>
0010 #include <iomanip>
0011 #include <cmath>
0012 
0013 using namespace sistrip;
0014 
0015 // -----------------------------------------------------------------------------
0016 //
0017 VpspScanAlgorithm::VpspScanAlgorithm(const edm::ParameterSet& pset, VpspScanAnalysis* const anal)
0018     : CommissioningAlgorithm(anal), histos_(2, Histo(nullptr, "")) {
0019   ;
0020 }
0021 
0022 // ----------------------------------------------------------------------------
0023 //
0024 void VpspScanAlgorithm::extract(const std::vector<TH1*>& histos) {
0025   if (!anal()) {
0026     edm::LogWarning(mlCommissioning_) << "[VpspScanAlgorithm::" << __func__ << "]"
0027                                       << " NULL pointer to Analysis object!";
0028     return;
0029   }
0030 
0031   // Check number of histograms
0032   if (histos.size() != 2) {
0033     anal()->addErrorCode(sistrip::numberOfHistos_);
0034   }
0035 
0036   // Extract FED key from histo title
0037   if (!histos.empty()) {
0038     anal()->fedKey(extractFedKey(histos.front()));
0039   }
0040 
0041   // Extract histograms
0042   std::vector<TH1*>::const_iterator ihis = histos.begin();
0043   for (; ihis != histos.end(); ihis++) {
0044     // Check pointer
0045     if (!(*ihis)) {
0046       continue;
0047     }
0048 
0049     // Check name
0050     SiStripHistoTitle title((*ihis)->GetName());
0051     if (title.runType() != sistrip::VPSP_SCAN) {
0052       anal()->addErrorCode(sistrip::unexpectedTask_);
0053       continue;
0054     }
0055 
0056     // Extract APV number
0057     uint16_t apv = sistrip::invalid_;
0058     if (title.extraInfo().find(sistrip::apv_) != std::string::npos) {
0059       std::stringstream ss;
0060       ss << title.extraInfo().substr(title.extraInfo().find(sistrip::apv_) + (sizeof(sistrip::apv_) - 1), 1);
0061       ss >> std::dec >> apv;
0062     }
0063 
0064     if (apv <= 1) {
0065       histos_[apv].first = *ihis;
0066       histos_[apv].second = (*ihis)->GetName();
0067     } else {
0068       anal()->addErrorCode(sistrip::unexpectedExtraInfo_);
0069     }
0070   }
0071 }
0072 
0073 // -----------------------------------------------------------------------------
0074 //
0075 void VpspScanAlgorithm::analyse() {
0076   if (!anal()) {
0077     edm::LogWarning(mlCommissioning_) << "[VpspScanAlgorithm::" << __func__ << "]"
0078                                       << " NULL pointer to base Analysis object!";
0079     return;
0080   }
0081 
0082   CommissioningAnalysis* tmp = const_cast<CommissioningAnalysis*>(anal());
0083   VpspScanAnalysis* anal = dynamic_cast<VpspScanAnalysis*>(tmp);
0084   if (!anal) {
0085     edm::LogWarning(mlCommissioning_) << "[VpspScanAlgorithm::" << __func__ << "]"
0086                                       << " NULL pointer to derived Analysis object!";
0087     return;
0088   }
0089 
0090   // from deprecated()
0091 
0092   std::vector<const TProfile*> histos;
0093   std::vector<uint16_t> monitorables;
0094   for (uint16_t iapv = 0; iapv < 2; iapv++) {
0095     monitorables.clear();
0096     monitorables.resize(7, sistrip::invalid_);
0097 
0098     histos.clear();
0099     histos.push_back(const_cast<const TProfile*>(dynamic_cast<TProfile*>(histos_[iapv].first)));
0100 
0101     if (!histos[0]) {
0102       anal->addErrorCode(sistrip::nullPtr_);
0103       continue;
0104     }
0105 
0106     // Find "top" plateau
0107     int first = sistrip::invalid_;
0108     float top = -1. * sistrip::invalid_;
0109     ;
0110     for (int k = 5; k < 55; k++) {
0111       if (!histos[0]->GetBinEntries(k)) {
0112         continue;
0113       }
0114       if (histos[0]->GetBinContent(k) >= top) {
0115         first = k;
0116         top = histos[0]->GetBinContent(k);
0117       }
0118     }
0119     if (top < -1. * sistrip::valid_) {
0120       top = sistrip::invalid_;
0121     }  //@@ just want +ve invalid number here
0122     if (top > 1. * sistrip::valid_) {
0123       anal->addErrorCode(sistrip::noTopPlateau_);
0124       continue;
0125     }
0126     monitorables[5] = static_cast<uint16_t>(top);
0127     monitorables[3] = first;
0128 
0129     // Find "bottom" plateau
0130     int last = sistrip::invalid_;
0131     float bottom = 1. * sistrip::invalid_;
0132     for (int k = 55; k > 5; k--) {
0133       if (!histos[0]->GetBinEntries(k)) {
0134         continue;
0135       }
0136       if (histos[0]->GetBinContent(k) <= bottom) {
0137         last = k;
0138         bottom = histos[0]->GetBinContent(k);
0139       }
0140     }
0141     if (bottom > 1. * sistrip::valid_) {
0142       anal->addErrorCode(sistrip::noBottomPlateau_);
0143       continue;
0144     }
0145     monitorables[6] = static_cast<uint16_t>(bottom);
0146     monitorables[4] = last;
0147 
0148     // Set optimum baseline level
0149     float opt = bottom + (top - bottom) * 1. / 3.;
0150     monitorables[1] = static_cast<uint16_t>(opt);
0151 
0152     // Find optimum VPSP setting
0153     uint16_t vpsp = sistrip::invalid_;
0154     if (opt < 1. * sistrip::valid_) {
0155       uint16_t ivpsp = 5;
0156       for (; ivpsp < 55; ivpsp++) {
0157         if (histos[0]->GetBinContent(ivpsp) < opt) {
0158           break;
0159         }
0160       }
0161       if (ivpsp != 54) {
0162         vpsp = ivpsp;
0163         monitorables[0] = vpsp;
0164       } else {
0165         anal->addErrorCode(sistrip::noVpspSetting_);
0166         continue;
0167       }
0168 
0169     } else {
0170       anal->addErrorCode(sistrip::noBaselineLevel_);
0171       continue;
0172     }
0173 
0174     // Set analysis results for both APVs
0175     anal->vpsp_[iapv] = monitorables[0];
0176     anal->adcLevel_[iapv] = monitorables[1];
0177     anal->fraction_[iapv] = monitorables[2];
0178     anal->topEdge_[iapv] = monitorables[3];
0179     anal->bottomEdge_[iapv] = monitorables[4];
0180     anal->topLevel_[iapv] = monitorables[5];
0181     anal->bottomLevel_[iapv] = monitorables[6];
0182   }
0183 }