Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #include "DQM/SiStripCommissioningClients/interface/PedsOnlyHistograms.h"
0002 #include "CondFormats/SiStripObjects/interface/PedsOnlyAnalysis.h"
0003 #include "DQM/SiStripCommissioningAnalysis/interface/PedsOnlyAlgorithm.h"
0004 #include "DQM/SiStripCommissioningSummary/interface/PedsOnlySummaryFactory.h"
0005 #include "DQM/SiStripCommon/interface/ExtractTObject.h"
0006 #include "DataFormats/SiStripCommon/interface/SiStripConstants.h"
0007 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0008 #include <iostream>
0009 #include <memory>
0010 
0011 #include <sstream>
0012 #include <iomanip>
0013 #include "TProfile.h"
0014 
0015 using namespace std;
0016 using namespace sistrip;
0017 
0018 // -----------------------------------------------------------------------------
0019 /** */
0020 PedsOnlyHistograms::PedsOnlyHistograms(const edm::ParameterSet& pset, DQMStore* bei)
0021     : CommissioningHistograms(pset.getParameter<edm::ParameterSet>("PedsOnlyParameters"), bei, sistrip::PEDS_ONLY) {
0022   factory_ = std::make_unique<PedsOnlySummaryFactory>();
0023   LogTrace(mlDqmClient_) << "[PedsOnlyHistograms::" << __func__ << "]"
0024                          << " Constructing object...";
0025 }
0026 
0027 // -----------------------------------------------------------------------------
0028 /** */
0029 PedsOnlyHistograms::~PedsOnlyHistograms() {
0030   LogTrace(mlDqmClient_) << "[PedsOnlyHistograms::" << __func__ << "]"
0031                          << " Destructing object...";
0032 }
0033 
0034 // -----------------------------------------------------------------------------
0035 /** */
0036 void PedsOnlyHistograms::histoAnalysis(bool debug) {
0037   LogTrace(mlDqmClient_) << "[PedsOnlyHistograms::" << __func__ << "]";
0038 
0039   // Some initialisation
0040   uint16_t valid = 0;
0041   HistosMap::const_iterator iter;
0042   Analyses::iterator ianal;
0043   std::map<std::string, uint16_t> errors;
0044 
0045   // Clear map holding analysis objects
0046   for (ianal = data().begin(); ianal != data().end(); ianal++) {
0047     if (ianal->second) {
0048       delete ianal->second;
0049     }
0050   }
0051   data().clear();
0052 
0053   // Iterate through map containing histograms
0054   for (iter = histos().begin(); iter != histos().end(); iter++) {
0055     // Check vector of histos is not empty
0056     if (iter->second.empty()) {
0057       edm::LogWarning(mlDqmClient_) << "[PedsOnlyHistograms::" << __func__ << "]"
0058                                     << " Zero histograms found!";
0059       continue;
0060     }
0061 
0062     // Retrieve pointers to profile histos
0063     std::vector<TH1*> profs;
0064     Histos::const_iterator ihis = iter->second.begin();
0065     for (; ihis != iter->second.end(); ihis++) {
0066       TProfile* prof = ExtractTObject<TProfile>().extract((*ihis)->me_);
0067       if (prof) {
0068         profs.push_back(prof);
0069       }
0070       //@@ Common mode histos?...
0071       //TH1F* his = ExtractTObject<TH1F>().extract( (*ihis)->me_ );
0072       //if ( his ) { profs.push_back(his); }
0073     }
0074 
0075     // Perform histo analysis
0076     PedsOnlyAnalysis* anal = new PedsOnlyAnalysis(iter->first);
0077     PedsOnlyAlgorithm algo(this->pset(), anal);
0078     algo.analysis(profs);
0079     data()[iter->first] = anal;
0080     if (anal->isValid()) {
0081       valid++;
0082     }
0083     if (!anal->getErrorCodes().empty()) {
0084       errors[anal->getErrorCodes()[0]]++;
0085     }
0086   }
0087 
0088   if (!histos().empty()) {
0089     edm::LogVerbatim(mlDqmClient_) << "[PedsOnlyHistograms::" << __func__ << "]"
0090                                    << " Analyzed histograms for " << histos().size() << " FED channels, of which "
0091                                    << valid << " (" << 100 * valid / histos().size() << "%) are valid.";
0092     if (!errors.empty()) {
0093       uint16_t count = 0;
0094       std::stringstream ss;
0095       ss << std::endl;
0096       std::map<std::string, uint16_t>::const_iterator ii;
0097       for (ii = errors.begin(); ii != errors.end(); ++ii) {
0098         ss << " " << ii->first << ": " << ii->second << std::endl;
0099         count += ii->second;
0100       }
0101       edm::LogWarning(mlDqmClient_) << "[PedsOnlyHistograms::" << __func__ << "]"
0102                                     << " Found " << count << " errors (" << 100 * count / histos().size()
0103                                     << "%): " << ss.str();
0104     }
0105   } else {
0106     edm::LogWarning(mlDqmClient_) << "[PedsOnlyHistograms::" << __func__ << "]"
0107                                   << " No histograms to analyze!";
0108   }
0109 }
0110 
0111 // -----------------------------------------------------------------------------
0112 /** */
0113 void PedsOnlyHistograms::printAnalyses() {
0114   Analyses::iterator ianal = data().begin();
0115   Analyses::iterator janal = data().end();
0116   for (; ianal != janal; ++ianal) {
0117     if (ianal->second) {
0118       std::stringstream ss;
0119       ianal->second->print(ss, 1);
0120       ianal->second->print(ss, 2);
0121       if (ianal->second->isValid()) {
0122         LogTrace(mlDqmClient_) << ss.str();
0123       } else {
0124         edm::LogWarning(mlDqmClient_) << ss.str();
0125       }
0126     }
0127   }
0128 }