DABC (Data Acquisition Backbone Core)  2.9.9
dabc_root.cxx
Go to the documentation of this file.
1 // dedicated to convert DABC files to ROOT files
2 // for instance, hierarchy stored in binary file, can be converted to ROOT histograms
3 
4 #include <cstdio>
5 #include <cstring>
6 
7 #include "dabc/version.h"
8 #include "dabc/Buffer.h"
9 #include "dabc/Hierarchy.h"
10 #include "dabc/BinaryFileIO.h"
11 #include "dabc/Iterator.h"
12 
13 #include "TFile.h"
14 #include "TH1.h"
15 #include "TH2.h"
16 #include "TObjArray.h"
17 #include "TAxis.h"
18 
19 void SetAxisLabels(TAxis *ax, std::string lbls, int nbins)
20 {
21  if (lbls.empty() || !ax) return;
22  TObjArray *lst = TString(lbls.c_str()).Tokenize(",");
23 
24  if (lst && (lst->GetSize() > 0))
25  for (int n=0;n<nbins;++n)
26  ax->SetBinLabel(n+1, (n < lst->GetSize()) ? lst->At(n)->GetName() : Form("bin%d",n));
27 
28  delete lst;
29 }
30 
31 int main(int argc, char* argv[]) {
32 
33  printf("dabc_root utility, v %s\n", DABC_RELEASE);
34 
35  const char* inpfile(0), *outfile(0);
36  bool skipzero(false);
37 
38  for (int n=1;n<argc;n++) {
39  if ((strcmp(argv[n],"-h")==0) && (n+1 < argc)) {
40  inpfile = argv[++n];
41  continue;
42  }
43  if ((strcmp(argv[n],"-o")==0) && (n+1 < argc)) {
44  outfile = argv[++n];
45  continue;
46  }
47  if (strcmp(argv[n],"-skip-zero")==0) {
48  skipzero = true;
49  continue;
50  }
51  }
52 
53  if ((inpfile==0) || (outfile==0)) {
54  printf("Some arguments:\n");
55  printf("-h filename - input file with stored hierarchy\n");
56  printf("-o filename.root - output ROOT file with histograms\n");
57  printf("-skip-zero - skip empty histograms\n");
58  return 1;
59  }
60 
61  dabc::BinaryFile inp;
62 
63  if (!inp.OpenReading(inpfile)) {
64  printf("Cannot open file %s for reading\n", inpfile);
65  return 1;
66  }
67 
68  uint64_t size, typ;
69 
70  if (!inp.ReadBufHeader(&size,&typ)) {
71  printf("Cannot read buffer header from file %s\n", inpfile);
72  return 1;
73  }
74 
76  buf.SetTypeId(typ);
77 
78  if (!inp.ReadBufPayload(buf.SegmentPtr(), size)) {
79  printf("Cannot read buffer %lu from file %s\n", (long unsigned) size, inpfile);
80  return 1;
81  }
82 
84  if (!h.ReadFromBuffer(buf)) {
85  printf("Cannot decode hierarchy from buffer\n");
86  return 1;
87  }
88 
89  TFile* f = TFile::Open(outfile,"recreate");
90  if (f==0) {
91  printf("Cannot create %s for output\n", outfile);
92  return 1;
93  }
94 
95  int cnt(0), emptycnt(0), currcnt(0);
96  TString lastdir;
97  dabc::Iterator iter(h);
98  while (iter.next()) {
99  dabc::Hierarchy item = iter.ref();
100  if (item.null()) continue;
101 
102  if (!item.HasField("_dabc_hist")) continue;
103 
104  TString dirname;
105 
106  for (unsigned n=1;n<iter.level();n++) {
107  dabc::Object* prnt = iter.parent(n);
108  if (prnt==0) break;
109 
110  if (n>1) dirname.Append("/");
111  dirname.Append(prnt->GetName());
112  }
113 
114  if (lastdir!=dirname) {
115  if ((lastdir.Length()>0) && (currcnt==0) && skipzero) {
116  if ((gDirectory->GetList()->GetSize()==0) && (gDirectory->GetListOfKeys()->GetSize()==0)) {
117  printf("Deleting empty directory %s\n", lastdir.Data());
118  f->rmdir(lastdir.Data());
119  }
120  }
121 
122  if (!f->GetDirectory(dirname.Data()))
123  f->mkdir(dirname.Data());
124  f->cd(dirname.Data());
125  currcnt = 0;
126  }
127 
128  lastdir = dirname;
129 
130  if (item.Field("_kind").AsStr() == "ROOT.TH1D") {
131  int nbins = item.Field("nbins").AsInt();
132  double* bins = item.GetFieldPtr("bins")->GetDoubleArr();
133 
134  if (skipzero) {
135  bool isany = false;
136  for (int n=0;n<nbins+2;n++) {
137  if (bins[3+n]!=0.) { isany = true; break; }
138  }
139  if (!isany) { emptycnt++; continue; }
140  }
141 
142  TH1D* hist = new TH1D(item.GetName(), item.Field("_title").AsStr().c_str(),
143  nbins, item.Field("left").AsDouble(), item.Field("right").AsDouble());
144 
145  if (item.HasField("xlabels"))
146  SetAxisLabels(hist->GetXaxis(), item.Field("xlabels").AsStr(), nbins);
147 
148  for (int n=0;n<nbins+2;n++) hist->SetBinContent(n,bins[3+n]);
149  hist->ResetStats(); // recalculate statistic
150  hist->Write();
151 
152  cnt++; currcnt++;
153  } else
154  if ((item.Field("_kind").AsStr() == "ROOT.TH2D") || (item.Field("_kind").AsStr() == "ROOT.TH2Poly")) {
155  int nbins1 = item.Field("nbins1").AsInt();
156  int nbins2 = item.Field("nbins2").AsInt();
157  double* bins = item.GetFieldPtr("bins")->GetDoubleArr();
158  if (skipzero) {
159  bool isany = false;
160  for (int n=0;n<(nbins1+2)*(nbins2+2);n++) {
161  if (bins[6+n]!=0.) { isany = true; break; }
162  }
163  if (!isany) { emptycnt++; continue; }
164  }
165  TH2D* hist = new TH2D(item.GetName(), item.Field("_title").AsStr().c_str(),
166  nbins1, item.Field("left1").AsDouble(), item.Field("right1").AsDouble(),
167  nbins2, item.Field("left2").AsDouble(), item.Field("right2").AsDouble());
168 
169  if (item.HasField("xlabels"))
170  SetAxisLabels(hist->GetXaxis(), item.Field("xlabels").AsStr(), nbins1);
171 
172  if (item.HasField("ylabels"))
173  SetAxisLabels(hist->GetYaxis(), item.Field("ylabels").AsStr(), nbins2);
174 
175  for (int n=0;n<(nbins1+2)*(nbins2+2);n++) hist->SetBinContent(n,bins[6+n]);
176  hist->ResetStats(); // recalculate statistic
177  hist->Write();
178 
179  cnt++; currcnt++;
180  }
181 
182  }
183 
184  printf("TOTAL %d histograms converted, skipped %d, storing file %s ...\n", cnt, emptycnt, outfile);
185 
186  delete f;
187 
188  printf("File %s stored\n", outfile);
189 
190  return 0;
191 }
Generic file storage for DABC buffers.
Definition: BinaryFile.h:181
bool ReadBufPayload(void *ptr, uint64_t sz)
Definition: BinaryFile.h:337
bool OpenReading(const char *fname)
Definition: BinaryFile.h:197
bool ReadBufHeader(uint64_t *size, uint64_t *typ=0)
Definition: BinaryFile.h:317
Reference on memory from memory pool.
Definition: Buffer.h:135
static Buffer CreateBuffer(BufferSize_t sz)
This static method create independent buffer for any other memory pools Therefore it can be used in s...
Definition: Buffer.cxx:419
void SetTypeId(unsigned tid)
Definition: Buffer.h:151
void * SegmentPtr(unsigned n=0) const
Returns pointer on the segment, no any boundary checks.
Definition: Buffer.h:171
Represents objects hierarchy of remote (or local) DABC process.
Definition: Hierarchy.h:285
const RecordField & Field(const std::string &name) const
Definition: Hierarchy.h:304
bool ReadFromBuffer(const dabc::Buffer &buf)
Read hierarchy from buffer.
Definition: Hierarchy.cxx:896
Iterator over objects hierarchy
Definition: Iterator.h:36
Object * next(bool goinside=true)
Definition: Iterator.cxx:44
Object * parent(unsigned lvl)
Returns parents of current object.
Definition: Iterator.cxx:91
unsigned level() const
Definition: Iterator.h:51
Reference ref() const
Definition: Iterator.h:50
Base class for most of the DABC classes.
Definition: Object.h:116
const char * GetName() const
Returns name of the object, thread safe
Definition: Object.h:295
double * GetDoubleArr() const
Definition: Record.h:329
std::string AsStr(const std::string &dflt="") const
Definition: Record.cxx:749
int64_t AsInt(int64_t dflt=0) const
Definition: Record.cxx:501
double AsDouble(double dflt=0.) const
Definition: Record.cxx:549
RecordField * GetFieldPtr(const std::string &name) const
Definition: Record.h:513
bool HasField(const std::string &name) const
Definition: Record.h:498
const char * GetName() const
Return name of referenced object, if object not assigned, returns "---".
Definition: Reference.cxx:167
bool null() const
Returns true if reference contains nullptr.
Definition: Reference.h:151
int main(int argc, char *argv[])
Definition: dabc_root.cxx:31
void SetAxisLabels(TAxis *ax, std::string lbls, int nbins)
Definition: dabc_root.cxx:19
#define DABC_RELEASE
Definition: version.h:32