Public Types | Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes | List of all members
om::HistoSet Class Reference

#include "/cvmfs/nova-development.opensciencegrid.org/novasoft/releases/N20-12-03/OnlineMonitoring/producer/HistoSet.h"

Inheritance diagram for om::HistoSet:
om::TickerSubscriber

Public Types

enum  updates_ {
  kRun = 1<<0, kSubrun = 1<<1, k30sec = 1<<2, k1min = 1<<3,
  k5min = 1<<4, k10min = 1<<5, k30min = 1<<6, kHour = 1<<7,
  k24hr = 1<<8, kUTC = 1<<9, kAll = 0xFFFF
}
 

Public Member Functions

 ~HistoSet ()
 
TH1F * GetTH1F (const char *nm)
 
TH2F * GetTH2F (const char *nm)
 
void GetNames (std::list< std::string > &h)
 
void WriteToRootFile (TFile *f)
 
void RunTicker ()
 Complete the TickerSubscriber interface. More...
 
void SubrunTicker ()
 
void ThirtySecTicker ()
 
void OneMinTicker ()
 
void FiveMinTicker ()
 
void TenMinTicker ()
 
void ThirtyMinTicker ()
 
void HourTicker ()
 
void TwentyFourHrTicker ()
 
void DeleteTH1F (TH1F *h)
 
void DeleteTH2F (TH2F *h)
 
TH1 * FindTH1 (const char *nm)
 
TH1F * FindTH1F (const char *nm)
 
TH2F * FindTH2F (const char *nm)
 
void UTCReset (int UTCHour)
 
void UTCResetTH1F (int UTCHour, TH1F *h)
 
void UTCResetTH2F (int UTCHour, TH2F *h)
 

Static Public Member Functions

static HistoSetInstance ()
 

Private Member Functions

 HistoSet ()
 
void MakeTH1FCopies (const char *base, const char *tag, unsigned int lookback)
 
void MakeTH2FCopies (const char *base, const char *tag, unsigned int lookback)
 
void CopyAndResetAll (unsigned int which)
 
void CopyAndResetOne (const std::string &nm, TH1 *h, unsigned int which)
 
void UTCReset ()
 

Private Attributes

std::map< std::string, TH1F * > fTH1F
 The collection of 1D histos. More...
 
std::map< std::string, TH2F * > fTH2F
 The collection of 2D histos. More...
 

Detailed Description

Hold and manage the collection of histograms created by the producer

Definition at line 21 of file HistoSet.h.

Member Enumeration Documentation

Enumerator
kRun 
kSubrun 
k30sec 
k1min 
k5min 
k10min 
k30min 
kHour 
k24hr 
kUTC 
kAll 

Definition at line 15 of file TickerSubscriber.h.

Constructor & Destructor Documentation

HistoSet::~HistoSet ( )

Definition at line 36 of file HistoSet.cxx.

36 { }
HistoSet::HistoSet ( )
private

Definition at line 19 of file HistoSet.cxx.

Referenced by Instance().

Member Function Documentation

void HistoSet::CopyAndResetAll ( unsigned int  which)
private

Definition at line 159 of file HistoSet.cxx.

References CopyAndResetOne(), fTH1F, fTH2F, prev(), om::regex_match(), and string.

Referenced by FiveMinTicker(), HourTicker(), OneMinTicker(), RunTicker(), SubrunTicker(), TenMinTicker(), ThirtyMinTicker(), ThirtySecTicker(), and TwentyFourHrTicker().

160 {
161  const std::string prev("*_prev_*");
162  //
163  // Besure to use copies of the TH1F and TH2F histogram sets since
164  // the original versions will grow if copies get made.
165  //
166  std::map<std::string,TH1F*> fTH1Fcopy(fTH1F);
167  std::map<std::string,TH1F*>::iterator i1 (fTH1Fcopy.begin());
168  std::map<std::string,TH1F*>::iterator i1end(fTH1Fcopy.end());
169  for (; i1!=i1end; ++i1) {
170  //
171  // Don't make copies of copies
172  //
173  if (regex_match(i1->first, prev)!=0) {
174  this->CopyAndResetOne(i1->first, i1->second, which);
175  }
176  }
177 
178  std::map<std::string,TH2F*> fTH2Fcopy(fTH2F);
179  std::map<std::string,TH2F*>::iterator i2 (fTH2Fcopy.begin());
180  std::map<std::string,TH2F*>::iterator i2end(fTH2Fcopy.end());
181  for (; i2!=i2end; ++i2) {
182  //
183  // Don't make copies of copies
184  //
185  if (regex_match(i2->first, prev)!=0) {
186  this->CopyAndResetOne(i2->first, i2->second, which);
187  }
188  }
189 }
void prev()
Definition: show_event.C:91
std::map< std::string, TH1F * > fTH1F
The collection of 1D histos.
Definition: HistoSet.h:69
void CopyAndResetOne(const std::string &nm, TH1 *h, unsigned int which)
Definition: HistoSet.cxx:296
int regex_match(const std::string &s, const std::string &p)
Definition: RegexMatch.cxx:7
std::map< std::string, TH2F * > fTH2F
The collection of 2D histos.
Definition: HistoSet.h:70
enum BeamMode string
void HistoSet::CopyAndResetOne ( const std::string nm,
TH1 *  h,
unsigned int  which 
)
private
Todo:
The code fragments below are carbon copies of one another with 1->2. Is there a more elegant solution to this?

Definition at line 296 of file HistoSet.cxx.

References om::cerr, allTimeWatchdog::endl, FindTH1F(), FindTH2F(), om::HistoData::fLookBack, om::HistoData::fReset, fTH1F, fTH2F, om::HistoData::fType, h1, h2, hd, om::HistoTable::Instance(), om::TickerSubscriber::k10min, om::TickerSubscriber::k1min, om::TickerSubscriber::k24hr, om::TickerSubscriber::k30min, om::TickerSubscriber::k30sec, om::TickerSubscriber::k5min, om::TickerSubscriber::kHour, om::TickerSubscriber::kRun, om::TickerSubscriber::kSubrun, om::kTH1F, om::kTH2F, om::HistoTable::LookUp(), MakeTH1FCopies(), MakeTH2FCopies(), Munits::nm, string, and getGoodRuns4SAM::tag.

Referenced by CopyAndResetAll().

299 {
300  //
301  // Find the specifications for this histogram
302  //
304  const HistoData* hd = 0;
305  hd = ht.LookUp(nm.c_str());
306  if (hd==0) return;
307 
308  //
309  // Check if this histogram is scheduled to be reset
310  //
311  if ((hd->fReset & which)==0) return;
312 
313  //
314  // Reset is scheduled - make a copy and reset the histogram
315  //
316  // Find the histogram in the collections
317  //
318  TH1F* h1 = 0;
319  TH2F* h2 = 0;
320  if (hd->fType == kTH1F) h1 = this->FindTH1F(nm.c_str());
321  if (hd->fType == kTH2F) h2 = this->FindTH2F(nm.c_str());
322  if (h1==0 && h2==0) {
323  std::cerr << __FILE__ << ":" << __LINE__
324  << "Can't find histogram " << nm << std::endl;
325  abort();
326  }
327 
328  //
329  // Prepare a name for the copy
330  //
331  std::string copyname = nm;
332  const char* tag=0;
333  switch (which) {
334  case TickerSubscriber::kRun: tag = "prev_run"; break;
335  case TickerSubscriber::kSubrun: tag = "prev_subrun"; break;
336  case TickerSubscriber::k30sec: tag = "prev_30sec"; break;
337  case TickerSubscriber::k1min: tag = "prev_1min"; break;
338  case TickerSubscriber::k5min: tag = "prev_5min"; break;
339  case TickerSubscriber::k10min: tag = "prev_10min"; break;
340  case TickerSubscriber::k30min: tag = "prev_30min"; break;
341  case TickerSubscriber::kHour: tag = "prev_hour"; break;
342  case TickerSubscriber::k24hr: tag = "prev_24hr"; break;
343  default: abort();
344  }
345  copyname += "_"; copyname += tag; copyname += "_01";
346 
347  this->MakeTH1FCopies(nm.c_str(), tag, hd->fLookBack);
348  this->MakeTH2FCopies(nm.c_str(), tag, hd->fLookBack);
349 
350  //
351  // If the copy histogram does not exist, we'll need to make it. Once
352  // we're sure we have a destination, make the copy and reset the
353  // original histogram.
354  //
355  /// \todo The code fragments below are carbon copies of one
356  /// another with 1->2. Is there a more elegant solution to this?
357  if (h1) {
358  TH1F* h1copy = 0;
359  h1copy = this->FindTH1F(copyname.c_str());
360  if (h1copy==0) {
361  h1copy = new TH1F(*h1);
362  fTH1F[copyname] = h1copy;
363  }
364  h1->Copy(*h1copy);
365  h1copy->SetName(copyname.c_str());
366  h1->Reset();
367  }
368  if (h2) {
369  TH2F* h2copy = 0;
370  h2copy = this->FindTH2F(copyname.c_str());
371  if (h2copy==0) {
372  h2copy = new TH2F(*h2);
373  fTH2F[copyname] = h2copy;
374  }
375  h2->Copy(*h2copy);
376  h2copy->SetName(copyname.c_str());
377  h2->Reset();
378  }
379 }
unsigned int fLookBack
How many copies to save in history.
Definition: HistoData.h:43
const HistoData * LookUp(const char *nm) const
Definition: HistoTable.cxx:293
TH1F * FindTH1F(const char *nm)
Definition: HistoSet.cxx:50
static constexpr Double_t nm
Definition: Munits.h:133
OStream cerr
Definition: OStream.cxx:7
void MakeTH1FCopies(const char *base, const char *tag, unsigned int lookback)
Definition: HistoSet.cxx:200
std::map< std::string, TH1F * > fTH1F
The collection of 1D histos.
Definition: HistoSet.h:69
TH1F * h2
Definition: plot.C:45
TH1F * h1
Histo_t fType
What kind of histogram is this?
Definition: HistoData.h:46
void MakeTH2FCopies(const char *base, const char *tag, unsigned int lookback)
Definition: HistoSet.cxx:248
std::map< std::string, TH2F * > fTH2F
The collection of 2D histos.
Definition: HistoSet.h:70
TH2F * FindTH2F(const char *nm)
Definition: HistoSet.cxx:59
static HistoTable & Instance(const char *f=0, Detector_t d=kALLDET)
Definition: HistoTable.cxx:21
unsigned int fReset
Reset schedule (see TickerSubscriber.h)
Definition: HistoData.h:42
TH1F * hd
Definition: Xdiff_gwt.C:57
enum BeamMode string
void HistoSet::DeleteTH1F ( TH1F *  h)

Definition at line 556 of file HistoSet.cxx.

References fTH1F.

Referenced by om::WatchListManager::EndSubrunCleanUp().

557 {
558  std::map<std::string,TH1F*>::iterator i1 (fTH1F.begin());
559  std::map<std::string,TH1F*>::iterator i1end(fTH1F.end());
560  for (; i1!=i1end; ++i1) {
561  if(i1->second == h) {
562  break;
563  }
564  }
565 
566  if(i1!=i1end) {
567  delete i1->second;
568  i1->second = 0;
569  fTH1F.erase(i1);
570  }
571 }
std::map< std::string, TH1F * > fTH1F
The collection of 1D histos.
Definition: HistoSet.h:69
void HistoSet::DeleteTH2F ( TH2F *  h)

Definition at line 575 of file HistoSet.cxx.

References fTH2F.

Referenced by om::WatchListManager::EndSubrunCleanUp().

576 {
577  std::map<std::string,TH2F*>::iterator i2 (fTH2F.begin());
578  std::map<std::string,TH2F*>::iterator i2end(fTH2F.end());
579  for (; i2!=i2end; ++i2) {
580  if(i2->second == h) {
581  break;
582  }
583  }
584 
585  if(i2!=i2end) {
586  delete i2->second;
587  i2->second = 0;
588  fTH2F.erase(i2);
589  }
590 }
std::map< std::string, TH2F * > fTH2F
The collection of 2D histos.
Definition: HistoSet.h:70
TH1 * HistoSet::FindTH1 ( const char *  nm)

Definition at line 40 of file HistoSet.cxx.

References FindTH1F(), FindTH2F(), and make_syst_table_plots::h.

41 {
42  TH1* h=0;
43  h = this->FindTH1F(nm); if (h!=0) return h;
44  h = this->FindTH2F(nm); if (h!=0) return h;
45  return 0;
46 }
TH1F * FindTH1F(const char *nm)
Definition: HistoSet.cxx:50
static constexpr Double_t nm
Definition: Munits.h:133
TH2F * FindTH2F(const char *nm)
Definition: HistoSet.cxx:59
TH1F * HistoSet::FindTH1F ( const char *  nm)

Definition at line 50 of file HistoSet.cxx.

References fTH1F.

Referenced by CopyAndResetOne(), FindTH1(), om::OnMonProdIPC::FindTH1F(), GetTH1F(), and MakeTH1FCopies().

51 {
52  std::map<std::string,TH1F*>::iterator itr(fTH1F.find(nm));
53  if (itr!=fTH1F.end()) return itr->second;
54  return 0;
55 }
static constexpr Double_t nm
Definition: Munits.h:133
std::map< std::string, TH1F * > fTH1F
The collection of 1D histos.
Definition: HistoSet.h:69
TH2F * HistoSet::FindTH2F ( const char *  nm)

Definition at line 59 of file HistoSet.cxx.

References fTH2F.

Referenced by CopyAndResetOne(), FindTH1(), om::OnMonProdIPC::FindTH2F(), GetTH2F(), and MakeTH2FCopies().

60 {
61  std::map<std::string,TH2F*>::iterator itr(fTH2F.find(nm));
62  if (itr!=fTH2F.end()) return itr->second;
63  return 0;
64 }
static constexpr Double_t nm
Definition: Munits.h:133
std::map< std::string, TH2F * > fTH2F
The collection of 2D histos.
Definition: HistoSet.h:70
void HistoSet::FiveMinTicker ( )
virtual

Reimplemented from om::TickerSubscriber.

Definition at line 516 of file HistoSet.cxx.

References CopyAndResetAll(), and om::TickerSubscriber::k5min.

517 {
519 }
void CopyAndResetAll(unsigned int which)
Definition: HistoSet.cxx:159
void HistoSet::GetNames ( std::list< std::string > &  h)

Definition at line 111 of file HistoSet.cxx.

References fTH1F, and fTH2F.

Referenced by om::OnMonProdIPC::HistoList().

112 {
113  std::map<std::string,TH1F*>::iterator i1 (fTH1F.begin());
114  std::map<std::string,TH1F*>::iterator i1end(fTH1F.end());
115  for (; i1!=i1end; ++i1) {
116  if(i1->second != 0) h.push_back(i1->first);
117  }
118  std::map<std::string,TH2F*>::iterator i2 (fTH2F.begin());
119  std::map<std::string,TH2F*>::iterator i2end(fTH2F.end());
120  for (; i2!=i2end; ++i2) {
121  if(i2->second != 0) h.push_back(i2->first);
122  }
123 }
std::map< std::string, TH1F * > fTH1F
The collection of 1D histos.
Definition: HistoSet.h:69
std::map< std::string, TH2F * > fTH2F
The collection of 2D histos.
Definition: HistoSet.h:70
TH1F * HistoSet::GetTH1F ( const char *  nm)

Definition at line 68 of file HistoSet.cxx.

References om::cerr, d, allTimeWatchdog::endl, FindTH1F(), om::HistoData::fNx, fTH1F, om::HistoData::fTitle, om::HistoData::fX1, om::HistoData::fX2, make_syst_table_plots::h, om::HistoTable::Instance(), om::HistoTable::LookUp(), and Munits::nm.

Referenced by om::WatchListManager::Add(), om::FEBRatesByDiblock::FEBRatesByDiblock(), om::MicroErrors::GetMicroSliceSummary(), om::HitMaps::HitMaps(), om::TQPlots::TQPlots(), and om::TriggerPlots::TriggerPlots().

69 {
70  // Check if we've alread booked this histogram
71  TH1F* h = this->FindTH1F(nm);
72  if (h!=0) return h;
73 
74  // Histogram not found. Get the data we need to book it
76  if (d==0) {
77  // Nothing is known about the named histogram?!
78  std::cerr << "\nNothing known about requested histogram " << nm
79  << "\n" << std::endl;
80  return 0;
81  }
82  h = new TH1F(nm,d->fTitle.c_str(),d->fNx,d->fX1,d->fX2);
83  fTH1F[nm] = h;
84  return h;
85 }
double fX1
Low edge of x range.
Definition: HistoData.h:50
const HistoData * LookUp(const char *nm) const
Definition: HistoTable.cxx:293
TH1F * FindTH1F(const char *nm)
Definition: HistoSet.cxx:50
static constexpr Double_t nm
Definition: Munits.h:133
OStream cerr
Definition: OStream.cxx:7
std::string fTitle
Titles for histogram.
Definition: HistoData.h:48
int fNx
Number of bins in x.
Definition: HistoData.h:49
Float_t d
Definition: plot.C:236
double fX2
High edge of x range.
Definition: HistoData.h:51
std::map< std::string, TH1F * > fTH1F
The collection of 1D histos.
Definition: HistoSet.h:69
static HistoTable & Instance(const char *f=0, Detector_t d=kALLDET)
Definition: HistoTable.cxx:21
TH2F * HistoSet::GetTH2F ( const char *  nm)

Definition at line 89 of file HistoSet.cxx.

References om::cerr, d, allTimeWatchdog::endl, FindTH2F(), om::HistoData::fNx, om::HistoData::fNy, fTH2F, om::HistoData::fTitle, om::HistoData::fX1, om::HistoData::fX2, om::HistoData::fY1, om::HistoData::fY2, make_syst_table_plots::h, om::HistoTable::Instance(), om::HistoTable::LookUp(), and Munits::nm.

Referenced by om::WatchListManager::Add(), om::DataBlockErrors::DataBlockErrors(), om::FEBRatesByDiblock::FEBRatesByDiblock(), om::NanoErrors::GetNanoSliceSummary(), om::TQPlots::GetNanoSliceSummary(), om::HitMaps::GetNanoSliceSummary(), om::FEBRatesByDiblock::GetNanoSliceSummary(), om::HitCounts::HitCounts(), om::HitMaps::HitMaps(), om::MicroErrors::MicroErrors(), om::NanoErrors::NanoErrors(), om::RawEventErrors::RawEventErrors(), om::TQPlots::TQPlots(), and om::TriggerPlots::TriggerPlots().

90 {
91  TH2F* h = this->FindTH2F(nm);
92  if (h!=0) return h;
93 
94  // Histogram not found. Get the data we need to book it
96  if (d==0) {
97  // Nothing is known about the named histogram?!
98  std::cerr << "\nNothing known about requested histogram " << nm
99  << "\n" << std::endl;
100  return 0;
101  }
102  h = new TH2F(nm, d->fTitle.c_str(),
103  d->fNx, d->fX1, d->fX2,
104  d->fNy, d->fY1, d->fY2);
105  fTH2F[nm] = h;
106  return h;
107 }
double fX1
Low edge of x range.
Definition: HistoData.h:50
const HistoData * LookUp(const char *nm) const
Definition: HistoTable.cxx:293
static constexpr Double_t nm
Definition: Munits.h:133
OStream cerr
Definition: OStream.cxx:7
std::string fTitle
Titles for histogram.
Definition: HistoData.h:48
int fNx
Number of bins in x.
Definition: HistoData.h:49
double fY2
High edge of y range.
Definition: HistoData.h:54
Float_t d
Definition: plot.C:236
double fX2
High edge of x range.
Definition: HistoData.h:51
std::map< std::string, TH2F * > fTH2F
The collection of 2D histos.
Definition: HistoSet.h:70
double fNy
Number of bins in y.
Definition: HistoData.h:52
TH2F * FindTH2F(const char *nm)
Definition: HistoSet.cxx:59
static HistoTable & Instance(const char *f=0, Detector_t d=kALLDET)
Definition: HistoTable.cxx:21
double fY1
Low edge of y range.
Definition: HistoData.h:53
void HistoSet::HourTicker ( )
virtual

Reimplemented from om::TickerSubscriber.

Definition at line 537 of file HistoSet.cxx.

References CopyAndResetAll(), om::TickerSubscriber::kHour, and UTCReset().

538 {
540 
541  time_t tt = time(0);
542  struct tm t;
543  gmtime_r(&tt, &t);
544  this->UTCReset(t.tm_hour);
545 }
void CopyAndResetAll(unsigned int which)
Definition: HistoSet.cxx:159
::xsd::cxx::tree::time< char, simple_type > time
Definition: Database.h:194
void UTCReset()
Definition: type_traits.h:56
HistoSet & HistoSet::Instance ( )
static
void HistoSet::MakeTH1FCopies ( const char *  base,
const char *  tag,
unsigned int  lookback 
)
private

Manage the collection of copies. If we are supposed to be able to look back through "lookback" number of copies, then search for histograms matching the pattern "base_tag_i". If I find base_tag_7, rename it base_tag_8 and so on. The last histogram in the list gets deleted to make room for the new ones.

Definition at line 200 of file HistoSet.cxx.

References FindTH1F(), fTH1F, MECModelEnuComparisons::i, and calib::j.

Referenced by CopyAndResetOne().

203 {
204  if (lookback==0) return;
205 
206  unsigned int i, j;
207  TH1F* h1i;
208  TH1F* h1j;
209  char ni[256], nj[256];
210  for (j=lookback+1; j>1; --j) {
211  //
212  // Find the histogram to be renamed / copied
213  //
214  i = j-1;
215  sprintf(ni, "%s_%s_%.2d", base, tag, i);
216  h1i = this->FindTH1F(ni);
217  if (h1i==0) continue;
218 
219  //
220  // If we've reached the look back count don't copy the histogram,
221  // delete it
222  //
223  if (i==lookback) {
224  delete h1i;
225  fTH1F[ni] = 0;
226  continue;
227  }
228 
229  //
230  // Look at the target destination. If its not empty, empty it.
231  //
232  sprintf(nj, "%s_%s_%.2d", base, tag, j);
233  h1j = this->FindTH1F(nj);
234  if (h1j!=0) { delete h1j; fTH1F[nj] = 0; }
235 
236  //
237  // Rename the first histogram as the second and move its place in
238  // the table
239  //
240  h1i->SetName(nj);
241  fTH1F[nj] = h1i;
242  fTH1F[ni] = 0;
243  }
244 }
TH1F * FindTH1F(const char *nm)
Definition: HistoSet.cxx:50
const XML_Char int const XML_Char int const XML_Char * base
Definition: expat.h:331
const double j
Definition: BetheBloch.cxx:29
std::map< std::string, TH1F * > fTH1F
The collection of 1D histos.
Definition: HistoSet.h:69
void HistoSet::MakeTH2FCopies ( const char *  base,
const char *  tag,
unsigned int  lookback 
)
private

Definition at line 248 of file HistoSet.cxx.

References FindTH2F(), fTH2F, MECModelEnuComparisons::i, and calib::j.

Referenced by CopyAndResetOne().

251 {
252  if (lookback==0) return;
253 
254  unsigned int i, j;
255  TH2F* h2i;
256  TH2F* h2j;
257  char ni[256], nj[256];
258  for (j=lookback+1; j>1; --j) {
259  //
260  // Find the histogram to be renamed / copied
261  //
262  i = j-1;
263  sprintf(ni, "%s_%s_%.2d", base, tag, i);
264  h2i = this->FindTH2F(ni);
265  if (h2i==0) continue;
266 
267  //
268  // If we've reached the look back count don't copy the histogram,
269  // delete it
270  //
271  if (i==lookback) {
272  delete h2i;
273  fTH2F[ni] = 0;
274  continue;
275  }
276 
277  //
278  // Look at the target destination. If its not empty, empty it.
279  //
280  sprintf(nj, "%s_%s_%.2d", base, tag, j);
281  h2j = this->FindTH2F(nj);
282  if (h2j!=0) { delete h2j; fTH2F[nj] = 0; }
283 
284  //
285  // Rename the first histogram as the second and move its place in
286  // the table
287  //
288  h2i->SetName(nj);
289  fTH2F[nj] = h2i;
290  fTH2F[ni] = 0;
291  }
292 }
const XML_Char int const XML_Char int const XML_Char * base
Definition: expat.h:331
const double j
Definition: BetheBloch.cxx:29
std::map< std::string, TH2F * > fTH2F
The collection of 2D histos.
Definition: HistoSet.h:70
TH2F * FindTH2F(const char *nm)
Definition: HistoSet.cxx:59
void HistoSet::OneMinTicker ( )
virtual

Reimplemented from om::TickerSubscriber.

Definition at line 509 of file HistoSet.cxx.

References CopyAndResetAll(), and om::TickerSubscriber::k1min.

510 {
512 }
void CopyAndResetAll(unsigned int which)
Definition: HistoSet.cxx:159
void HistoSet::RunTicker ( )
virtual

Complete the TickerSubscriber interface.

Reimplemented from om::TickerSubscriber.

Definition at line 488 of file HistoSet.cxx.

References CopyAndResetAll(), and om::TickerSubscriber::kRun.

489 {
491 }
void CopyAndResetAll(unsigned int which)
Definition: HistoSet.cxx:159
void HistoSet::SubrunTicker ( )
virtual

Reimplemented from om::TickerSubscriber.

Definition at line 495 of file HistoSet.cxx.

References CopyAndResetAll(), and om::TickerSubscriber::kSubrun.

496 {
498 }
void CopyAndResetAll(unsigned int which)
Definition: HistoSet.cxx:159
void HistoSet::TenMinTicker ( )
virtual

Reimplemented from om::TickerSubscriber.

Definition at line 523 of file HistoSet.cxx.

References CopyAndResetAll(), and om::TickerSubscriber::k10min.

524 {
526 }
void CopyAndResetAll(unsigned int which)
Definition: HistoSet.cxx:159
void HistoSet::ThirtyMinTicker ( )
virtual

Reimplemented from om::TickerSubscriber.

Definition at line 530 of file HistoSet.cxx.

References CopyAndResetAll(), and om::TickerSubscriber::k30min.

531 {
533 }
void CopyAndResetAll(unsigned int which)
Definition: HistoSet.cxx:159
void HistoSet::ThirtySecTicker ( )
virtual

Reimplemented from om::TickerSubscriber.

Definition at line 502 of file HistoSet.cxx.

References CopyAndResetAll(), and om::TickerSubscriber::k30sec.

503 {
505 }
void CopyAndResetAll(unsigned int which)
Definition: HistoSet.cxx:159
void HistoSet::TwentyFourHrTicker ( )
virtual

Reimplemented from om::TickerSubscriber.

Definition at line 549 of file HistoSet.cxx.

References CopyAndResetAll(), and om::TickerSubscriber::k24hr.

550 {
552 }
void CopyAndResetAll(unsigned int which)
Definition: HistoSet.cxx:159
void HistoSet::UTCReset ( int  UTCHour)

Definition at line 383 of file HistoSet.cxx.

References om::HistoData::fReset, fTH1F, fTH2F, hd, om::HistoTable::Instance(), om::TickerSubscriber::kUTC, om::HistoTable::LookUp(), UTCResetTH1F(), and UTCResetTH2F().

384 {
386  const HistoData* hd = 0;
387 
388  std::map<std::string,TH1F*>::iterator i1 (fTH1F.begin());
389  std::map<std::string,TH1F*>::iterator i1end(fTH1F.end());
390  for (; i1!=i1end; ++i1) {
391  hd = ht.LookUp(i1->first.c_str());
392  if (hd==0) continue;
393  if (hd->fReset == kUTC) {
394  this->UTCResetTH1F(UTCHour, i1->second);
395  }
396  }
397 
398  std::map<std::string,TH2F*>::iterator i2 (fTH2F.begin());
399  std::map<std::string,TH2F*>::iterator i2end(fTH2F.end());
400  for (; i2!=i2end; ++i2) {
401  hd = ht.LookUp(i2->first.c_str());
402  if (hd==0) continue;
403  if (hd->fReset == kUTC) {
404  this->UTCResetTH2F(UTCHour, i2->second);
405  }
406  }
407 }
const HistoData * LookUp(const char *nm) const
Definition: HistoTable.cxx:293
void UTCResetTH1F(int UTCHour, TH1F *h)
Definition: HistoSet.cxx:411
void UTCResetTH2F(int UTCHour, TH2F *h)
Definition: HistoSet.cxx:446
std::map< std::string, TH1F * > fTH1F
The collection of 1D histos.
Definition: HistoSet.h:69
std::map< std::string, TH2F * > fTH2F
The collection of 2D histos.
Definition: HistoSet.h:70
static HistoTable & Instance(const char *f=0, Detector_t d=kALLDET)
Definition: HistoTable.cxx:21
unsigned int fReset
Reset schedule (see TickerSubscriber.h)
Definition: HistoData.h:42
TH1F * hd
Definition: Xdiff_gwt.C:57
void om::HistoSet::UTCReset ( )
private

Referenced by HourTicker().

void HistoSet::UTCResetTH1F ( int  UTCHour,
TH1F *  h 
)

Definition at line 411 of file HistoSet.cxx.

References stan::math::floor(), MECModelEnuComparisons::i, makeTrainCVSamples::int, and nbinsx.

Referenced by UTCReset().

412 {
413  //
414  // Schedule for reset is to zero the plot not for *this* UTC hour,
415  // but for the next UTC hour. So, if we are to start filling the
416  // hour between 8 and 9 we clear the bins between 9 and 10.
417  //
418  int reset_hour = (UTCHour+1)%24;
419 
420  int i;
421  int nbinsx = h->GetNbinsX();
422  double binutc; // UTC hour at center of bin
423  int binutci; // UTC hour at center of bin, floored to int
424 
425  //
426  // Loop over x-axis looking for bins that match the UTC hour we
427  //
428  for (i=1; i<=nbinsx; ++i) {
429  binutc = h->GetXaxis()->GetBinCenter(i);
430  //
431  // If the histogram has been binned to include UTC hours out of
432  // the standard range of [0...24), restore the hour to something
433  // modulo 24 hours.
434  //
435  while (binutc< 0.0) binutc += 24.0;
436  while (binutc>=24.0) binutc -= 24.0;
437  binutci = (int)floor(binutc);
438  if (binutci==reset_hour) {
439  h->SetBinContent(i,0);
440  }
441  }
442 }
Int_t nbinsx
Definition: plot.C:23
fvar< T > floor(const fvar< T > &x)
Definition: floor.hpp:11
void HistoSet::UTCResetTH2F ( int  UTCHour,
TH2F *  h 
)

Definition at line 446 of file HistoSet.cxx.

References stan::math::floor(), MECModelEnuComparisons::i, makeTrainCVSamples::int, calib::j, and nbinsx.

Referenced by UTCReset().

447 {
448  //
449  // Schedule for reset is to zero the plot not for *this* UTC hour,
450  // but for the next UTC hour. So, if we are to start filling the
451  // hour between 8 and 9 we clear the bins between 9 and 10.
452  //
453  int reset_hour = (UTCHour+1)%24;
454 
455  int i, j;
456  int nbinsx = h->GetNbinsX();
457  int nbinsy = h->GetNbinsY();
458  double binutc; // UTC hour at center of bin
459  int binutci; // UTC hour at center of bin, floored to int
460 
461  //
462  // Loop over x-axis looking for bins that match the UTC hour we
463  //
464  for (i=1; i<=nbinsx; ++i) {
465  binutc = h->GetXaxis()->GetBinCenter(i);
466  //
467  // If the histogram has been binned to include UTC hours out of
468  // the standard range of [0...24), restore the hour to something
469  // modulo 24 hours.
470  //
471  while (binutc< 0.0) binutc += 24.0;
472  while (binutc>=24.0) binutc -= 24.0;
473  binutci = (int)floor(binutc);
474 
475  //
476  // Reset the columns to zero if the hours match
477  //
478  if (binutci==reset_hour) {
479  for (j=1; j<=nbinsy; ++j) {
480  h->SetBinContent(i,j,0);
481  }
482  }
483  }
484 }
const double j
Definition: BetheBloch.cxx:29
Int_t nbinsx
Definition: plot.C:23
fvar< T > floor(const fvar< T > &x)
Definition: floor.hpp:11
void HistoSet::WriteToRootFile ( TFile *  f)
Todo:
  • Should make folders inside the file to keep the histograms organized

Definition at line 129 of file HistoSet.cxx.

References fTH1F, and fTH2F.

Referenced by om::OnMonProd::NewRunNumber(), and om::OnMonProd::NewSubrunNumber().

130 {
131  f->cd();
132  std::map<std::string,TH1F*>::iterator i1 (fTH1F.begin());
133  std::map<std::string,TH1F*>::iterator i1end(fTH1F.end());
134 
135  // add each histos to file as a new histo (use for debugging)
136  // for (; i1!=i1end; ++i1) i1->second->Write();
137 
138  // overwrite old histos in file with the same name. QUESTION: should
139  // buffer size be zero??? (3rd argument - 0 is the default value)
140  for (; i1!=i1end; ++i1) {
141  i1->second->Write(i1->first.c_str(),TObject::kOverwrite,0);
142  }
143 
144  std::map<std::string,TH2F*>::iterator i2 (fTH2F.begin());
145  std::map<std::string,TH2F*>::iterator i2end(fTH2F.end());
146 
147  // add each histos to file as a new histo (use for debugging)
148  // for (; i2!=i2end; ++i2) i2->second->Write();
149 
150  // overwrite old histos in file with the same name. QUESTION: should
151  // buffer size be zero??? (3rd argument - 0 is the default value)
152  for (; i2!=i2end; ++i2) {
153  i2->second->Write(i2->first.c_str(),TObject::kOverwrite,0);
154  }
155 }
std::map< std::string, TH1F * > fTH1F
The collection of 1D histos.
Definition: HistoSet.h:69
std::map< std::string, TH2F * > fTH2F
The collection of 2D histos.
Definition: HistoSet.h:70

Member Data Documentation

std::map<std::string,TH1F*> om::HistoSet::fTH1F
private

The collection of 1D histos.

Definition at line 69 of file HistoSet.h.

Referenced by CopyAndResetAll(), CopyAndResetOne(), DeleteTH1F(), FindTH1F(), GetNames(), GetTH1F(), MakeTH1FCopies(), UTCReset(), and WriteToRootFile().

std::map<std::string,TH2F*> om::HistoSet::fTH2F
private

The collection of 2D histos.

Definition at line 70 of file HistoSet.h.

Referenced by CopyAndResetAll(), CopyAndResetOne(), DeleteTH2F(), FindTH2F(), GetNames(), GetTH2F(), MakeTH2FCopies(), UTCReset(), and WriteToRootFile().


The documentation for this class was generated from the following files: