CMS 3D CMS Logo

CalibrationHistograms.cc
Go to the documentation of this file.
13 #include <iostream>
14 #include <sstream>
15 #include <iomanip>
16 #include "TH1F.h"
17 #include "TFile.h"
18 #include "TMultiGraph.h"
19 #include "TGraph.h"
20 
21 using namespace std;
22 using namespace sistrip;
23 
24 // -----------------------------------------------------------------------------
27  : CommissioningHistograms(pset.getParameter<edm::ParameterSet>("CalibrationParameters"), bei, task) {
28  LogTrace(mlDqmClient_) << "[CalibrationHistograms::" << __func__ << "]"
29  << " Constructing object...";
30 
32  factory_ = auto_ptr<CalibrationScanSummaryFactory>(new CalibrationScanSummaryFactory);
33  else
34  factory_ = auto_ptr<CalibrationSummaryFactory>(new CalibrationSummaryFactory);
35 
37  this->pset().existsAs<double>("targetRiseTime") ? this->pset().getParameter<double>("targetRiseTime") : 50;
39  this->pset().existsAs<double>("targetDecayTime") ? this->pset().getParameter<double>("targetDecayTime") : 125;
41  this->pset().existsAs<bool>("tuneSimultaneously") ? this->pset().getParameter<bool>("tuneSimultaneously") : false;
42 }
43 
44 // -----------------------------------------------------------------------------
47  LogTrace(mlDqmClient_) << "[CalibrationHistograms::" << __func__ << "]"
48  << " Deleting object...";
49 }
50 
51 // -----------------------------------------------------------------------------
54  // Clear map holding analysis objects
55  Analyses::iterator ianal;
56  for (ianal = data().begin(); ianal != data().end(); ianal++) {
57  if (ianal->second) {
58  delete ianal->second;
59  }
60  }
61  data().clear();
62 
63  // Iterate through map containing vectors of profile histograms
64  HistosMap::const_iterator iter = histos().begin();
65 
66  // One entry for each LLD channel --> differnt thousand entries
67  for (; iter != histos().end(); iter++) {
68  if (iter->second.empty()) {
69  edm::LogWarning(mlDqmClient_) << "[CalibrationHistograms::" << __func__ << "]"
70  << " Zero collation histograms found!";
71  continue;
72  }
73 
74  // Retrieve pointers to 1D histos for this FED channel --> all strips in the fiber = 256
75  vector<TH1*> profs;
76  Histos::const_iterator ihis = iter->second.begin();
77  for (; ihis != iter->second.end(); ihis++) {
78  TH1F* prof = ExtractTObject<TH1F>().extract((*ihis)->me_);
79  if (prof) {
80  profs.push_back(prof);
81  }
82  }
83 
84  // Perform histo analysis
85  bool isdeconv = false;
87  isdeconv = true;
88 
90  CalibrationScanAnalysis* anal = new CalibrationScanAnalysis(iter->first, isdeconv);
91  CalibrationScanAlgorithm algo(this->pset(), anal);
92  algo.analysis(profs);
93  data()[iter->first] = anal;
94 
95  // tune the parameters for this a given target
96  for (int iapv = 0; iapv < 2; iapv++) {
99  else
101  algo.fillTunedObservables(iapv);
102  }
103  } else {
104  CalibrationAnalysis* anal = new CalibrationAnalysis(iter->first, isdeconv);
105  CalibrationAlgorithm algo(this->pset(), anal);
106  algo.analysis(profs);
107  data()[iter->first] = anal;
108  }
109  }
110 }
111 
112 // -----------------------------------------------------------------------------
115  Analyses::iterator ianal = data().begin();
116  Analyses::iterator janal = data().end();
117  for (; ianal != janal; ++ianal) {
118  if (ianal->second) {
119  std::stringstream ss;
120  ianal->second->print(ss, 0);
121  ianal->second->print(ss, 1);
122  if (ianal->second->isValid()) {
123  LogTrace(mlDqmClient_) << ss.str();
124  } else {
125  edm::LogWarning(mlDqmClient_) << ss.str();
126  }
127  }
128  }
129 }
130 
131 //-----------------------------------------------------------------------------
134  // Construct path and filename
135  std::stringstream ss;
136  if (!path.empty()) { // create with a specific outputName
137  ss << path;
138  if (ss.str().find(".root") == std::string::npos) {
139  ss << ".root";
140  }
141 
142  } else {
143  // Retrieve SCRATCH directory
144  std::string scratch = "SCRATCH";
145  std::string dir = "";
146  if (std::getenv(scratch.c_str()) != nullptr) {
147  dir = std::getenv(scratch.c_str());
148  }
149 
150  // Add directory path
151  if (!dir.empty()) {
152  ss << dir << "/";
153  } else {
154  ss << "/tmp/";
155  }
156 
157  // Add filename with run number and ".root" extension
158  if (partitionName.empty())
159  ss << sistrip::dqmClientFileName_ << "_" << std::setfill('0') << std::setw(8) << run_number << ".root";
160  else
161  ss << sistrip::dqmClientFileName_ << "_" << partitionName << "_" << std::setfill('0') << std::setw(8)
162  << run_number << ".root";
163  }
164 
165  // Save file with appropriate filename
166  LogTrace(mlDqmClient_) << "[CommissioningHistograms::" << __func__ << "]"
167  << " Saving histograms to root file"
168  << " (This may take some time!)";
169  path = ss.str();
170  bei()->save(path, sistrip::collate_);
171  edm::LogVerbatim(mlDqmClient_) << "[CommissioningHistograms::" << __func__ << "]"
172  << " Saved histograms to root file \"" << ss.str() << "\"!";
173 
174  // In case of calibration-scan, add also the TGraphs
175  // re-open the file
176  TFile* outputFile = TFile::Open(path.c_str(), "UPDATE");
177  outputFile->cd();
178 
179  // get all sub-dirs
180  std::vector<std::string> contents;
181  bei()->getContents(contents);
182 
183  TMultiGraph* graph_isha = new TMultiGraph("riseTime_vs_isha", "");
184  TMultiGraph* graph_vfs = new TMultiGraph("decayTime_vs_vfs", "");
185 
186  bool save_graph_isha = false;
187  bool save_graph_vfs = false;
188 
189  // loop on the analysis objects which are storing all relevant results
190  Analyses::iterator ianal = data().begin();
191  Analyses::iterator janal = data().end();
192  for (; ianal != janal; ++ianal) {
193  if (ianal->second) {
194  CalibrationScanAnalysis* anal = dynamic_cast<CalibrationScanAnalysis*>(ianal->second);
195  SiStripFecKey feckey = anal->fecKey();
196 
197  TString directory;
198  for (auto content : contents) {
199  std::vector<std::string> tokens;
201  std::istringstream tokenStream(content);
202  while (std::getline(tokenStream, token, ':')) {
203  tokens.push_back(token);
204  }
205  directory = Form("%s", tokens.at(0).c_str());
206  if (directory.Contains(Form("FecCrate%d", feckey.fecCrate())) and
207  directory.Contains(Form("FecRing%d", feckey.fecRing())) and
208  directory.Contains(Form("FecSlot%d", feckey.fecSlot())) and
209  directory.Contains(Form("CcuAddr%d", feckey.ccuAddr())) and
210  directory.Contains(Form("CcuChan%d", feckey.ccuChan())))
211  break;
212  }
213 
214  outputFile->cd("DQMData/" + directory);
215 
216  for (size_t igraph = 0; igraph < anal->decayTimeVsVFS().size(); igraph++) {
217  graph_vfs->Add(anal->decayTimeVsVFS()[igraph]);
218  anal->decayTimeVsVFS()[igraph]->Write();
219  save_graph_vfs = true;
220  }
221 
222  for (size_t igraph = 0; igraph < anal->riseTimeVsISHA().size(); igraph++) {
223  graph_isha->Add(anal->riseTimeVsISHA()[igraph]);
224  anal->riseTimeVsISHA()[igraph]->Write();
225  save_graph_isha = true;
226  }
227 
228  for (size_t igraph = 0; igraph < anal->riseTimeVsISHAVsVFS().size(); igraph++)
229  anal->riseTimeVsISHAVsVFS()[igraph]->Write();
230 
231  for (size_t igraph = 0; igraph < anal->decayTimeVsISHAVsVFS().size(); igraph++)
232  anal->decayTimeVsISHAVsVFS()[igraph]->Write();
233 
234  outputFile->cd();
235  }
236  }
237 
238  outputFile->cd();
239  outputFile->cd("DQMData/Collate/SiStrip/ControlView");
240 
241  if (save_graph_isha)
242  graph_isha->Write("riseTime_vs_isha");
243  if (save_graph_vfs)
244  graph_vfs->Write("decayTime_vs_vfs");
245 
246  outputFile->Close();
247 }
void save(std::string &filename, uint32_t run_number=0, std::string partitionName="")
void analysis(const std::vector< TH1 * > &)
const std::vector< TGraph2D * > & decayTimeVsISHAVsVFS()
T getParameter(std::string const &) const
const sistrip::RunType & task() const
bool existsAs(std::string const &parameterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition: ParameterSet.h:160
Analyses & data(bool getMaskedData=false)
const uint16_t & fecRing() const
const edm::ParameterSet & pset() const
const std::vector< TGraph * > & decayTimeVsVFS()
Algorithm for calibration runs.
static const char dqmClientFileName_[]
std::vector< MonitorElement * > getContents(std::string const &path) const
Definition: DQMStore.cc:1520
static const char mlDqmClient_[]
const uint16_t & fecSlot() const
sistrip classes
Utility class that identifies a position within the strip tracker control structure, down to the level of an APV25.
Definition: SiStripFecKey.h:45
Analysis for calibration runs.
CalibrationHistograms(const edm::ParameterSet &pset, DQMStore *, const sistrip::RunType &task=sistrip::CALIBRATION)
Analysis for calibration scans.
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e< void, edm::EventID const &, edm::Timestamp const & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12
Algorithm for calibration runs.
DQMStore *const bei() const
const uint32_t & fecKey() const
#define LogTrace(id)
const uint16_t & fecCrate() const
int extract(std::vector< int > *output, const std::string &dati)
#define debug
Definition: HDRShower.cc:19
std::unique_ptr< Factory > factory_
const uint16_t & ccuAddr() const
#define begin
Definition: vmac.h:32
HLT enums.
void histoAnalysis(bool debug) override
void tuneIndependently(const int &, const float &, const float &)
const uint16_t & ccuChan() const
const std::vector< TGraph2D * > & riseTimeVsISHAVsVFS()
static const char collate_[]
void tuneSimultaneously(const int &, const float &, const float &)
void save(std::string const &filename, std::string const &path="", std::string const &pattern="", std::string const &rewrite="", uint32_t run=0, uint32_t lumi=0, SaveReferenceTag ref=SaveWithReference, int minStatus=dqm::qstatus::STATUS_OK, std::string const &fileupdate="RECREATE")
Definition: DQMStore.cc:2244
const HistosMap & histos() const
const std::vector< TGraph * > & riseTimeVsISHA()