CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
SamplingHistograms.cc
Go to the documentation of this file.
9 #include <iostream>
10 #include <memory>
11 
12 #include <sstream>
13 #include <iomanip>
14 #include "TProfile.h"
15 
16 using namespace std;
17 using namespace sistrip;
18 
19 // -----------------------------------------------------------------------------
22  : CommissioningHistograms(pset.getParameter<edm::ParameterSet>("SamplingParameters"), bei, task), sOnCut_(3) {
23  LogTrace(mlDqmClient_) << "[SamplingHistograms::" << __func__ << "]"
24  << " Constructing object...";
25  factory_ = std::make_unique<SamplingSummaryFactory>();
26  // retreive the latency code from the root file
27  std::string dataPath = std::string(sistrip::collate_) + "/" + sistrip::root_ + "/latencyCode";
28  MonitorElement* codeElement = bei->get(dataPath);
29  if (codeElement)
30  latencyCode_ = codeElement->getIntValue();
31  else
32  latencyCode_ = 0;
33 }
34 
35 // -----------------------------------------------------------------------------
38  LogTrace(mlDqmClient_) << "[SamplingHistograms::" << __func__ << "]"
39  << " Deleting object...";
40 }
41 
42 // -----------------------------------------------------------------------------
45  // Clear map holding analysis objects
46  Analyses::iterator ianal;
47  for (ianal = data().begin(); ianal != data().end(); ianal++) {
48  if (ianal->second) {
49  delete ianal->second;
50  }
51  }
52  data().clear();
53 
54  // Iterate through map containing vectors of profile histograms
55  HistosMap::const_iterator iter = histos().begin();
56  for (; iter != histos().end(); iter++) {
57  // Check vector of histos is not empty (should be 1 histo)
58  if (iter->second.empty()) {
59  edm::LogWarning(mlDqmClient_) << "[SamplingHistograms::" << __func__ << "]"
60  << " Zero collation histograms found!";
61  continue;
62  }
63 
64  // Retrieve pointers to profile histos for this FED channel
65  vector<TH1*> profs;
66  Histos::const_iterator ihis = iter->second.begin();
67  for (; ihis != iter->second.end(); ihis++) {
68  TProfile* prof = ExtractTObject<TProfile>().extract((*ihis)->me_);
69  if (prof) {
70  profs.push_back(prof);
71  }
72  }
73 
74  // Perform histo analysis
75  SamplingAnalysis* anal = new SamplingAnalysis(iter->first);
76  anal->setSoNcut(sOnCut_);
77  SamplingAlgorithm algo(this->pset(), anal, latencyCode_);
78  algo.analysis(profs);
79  data()[iter->first] = anal;
80  }
81 }
82 
84  //TODO: should use the parameter set. Why is this crashing ???
85  // sOnCut_ = pset.getParameter<double>("SignalToNoiseCut");
86  sOnCut_ = 3.;
87 }
void analysis(const std::vector< TH1 * > &)
Analysis for latency run.
SamplingHistograms(const edm::ParameterSet &pset, DQMStore *, const sistrip::RunType &task=sistrip::APV_LATENCY)
Algorithm for latency run.
Analyses & data(bool getMaskedData=false)
const edm::ParameterSet & pset() const
void setSoNcut(const float sOnCut)
virtual int64_t getIntValue() const
static const char mlDqmClient_[]
#define LogTrace(id)
void configure(const edm::ParameterSet &, const edm::EventSetup &) override
virtual MonitorElement * get(std::string const &fullpath) const
Definition: DQMStore.cc:673
~SamplingHistograms() override
int extract(std::vector< int > *output, const std::string &dati)
void histoAnalysis(bool debug) override
#define debug
Definition: HDRShower.cc:19
std::unique_ptr< Factory > factory_
static const char root_[]
Log< level::Warning, false > LogWarning
static const char collate_[]
const HistosMap & histos() const