BCDataSet.h
1 #ifndef __BCDATASET__H
2 #define __BCDATASET__H
3 
16 /*
17  * Copyright (C) 2007-2018, the BAT core developer team
18  * All rights reserved.
19  *
20  * For the licensing terms see doc/COPYING.
21  * For documentation see http://mpp.mpg.de/bat
22  */
23 
24 // ---------------------------------------------------------
25 
26 #include <string>
27 #include <vector>
28 
29 #include "BCDataPoint.h"
30 #include "BCLog.h"
31 
32 class TGraph;
33 class TGraphErrors;
34 class TGraphAsymmErrors;
35 class TH2;
36 
37 // ---------------------------------------------------------
38 
39 class BCDataSet
40 {
41 public:
42 
49  BCDataSet(unsigned n = 0);
50 
53  virtual ~BCDataSet() {};
54 
61  BCDataPoint& operator[](unsigned index)
62  { return fDataVector[index]; }
63 
66  const BCDataPoint& operator[](unsigned index) const
67  { return fDataVector[index]; }
68 
75  unsigned GetNDataPoints() const
76  { return fDataVector.size(); }
77 
80  unsigned GetNValuesPerPoint() const
81  { return fNValuesPerPoint; }
82 
87  BCDataPoint& GetDataPoint(unsigned index)
88  { return fDataVector.at(index); }
89 
94  const BCDataPoint& GetDataPoint(unsigned index) const
95  { return fDataVector.at(index); }
96 
100  { return fDataVector.back(); }
101 
107  std::vector<double> GetDataComponents(unsigned index) const;
108 
111  bool BoundsExist() const;
112 
116  { return fLowerBounds; }
117 
121  { return fUpperBounds; }
122 
126  { return fUserLowerBounds; }
127 
131  { return fUserUpperBounds; }
132 
137  double GetLowerBound(unsigned index) const;
138 
143  double GetUpperBound(unsigned index) const;
144 
149  double GetRangeWidth(unsigned index) const
150  { return GetUpperBound(index) - GetLowerBound(index); }
151 
156  bool IsFixed(unsigned index) const
157  { return (index < fFixed.size() and fFixed[index]); }
158 
168  void SetNValuesPerPoint(unsigned n);
169 
176  void SetBounds(unsigned index, double lower_bound, double upper_bound, bool fixed = false);
177 
182  void Fix(unsigned i, bool b = true)
183  { if (i < fFixed.size()) fFixed[i] = b; }
184 
197  bool ReadDataFromFile(const std::string& filename, const std::string& treename, const std::string& branchnames, char delim = ',')
198  { return ReadDataFromFileTree(filename, treename, branchnames, delim); };
199 
205  bool ReadDataFromFile(const std::string& filename, int nvariables)
206  { return ReadDataFromFileTxt(filename, nvariables); };
207 
217  bool ReadDataFromFileTree(const std::string& filename, const std::string& treename, const std::string& branchnames, char delim = ',');
218 
225  bool ReadDataFromFileTxt(const std::string& filename, int nvariables);
226 
230  bool AddDataPoint(const BCDataPoint& datapoint);
231 
241  void AdjustBoundForUncertainties(unsigned i, double nSigma, unsigned i_err1, int i_err2 = -1);
242 
245  void Reset()
246  { fDataVector.clear(); SetNValuesPerPoint(0); }
247 
251  void PrintSummary(void (*output)(const std::string&) = BCLog::OutSummary) const;
252 
261  TGraph* GetGraph(unsigned x, unsigned y) const;
262 
274  TGraphErrors* GetGraph(unsigned x, unsigned y, int ex, int ey) const;
275 
289  TGraphAsymmErrors* GetGraph(unsigned x, unsigned y, int ex_below, int ex_above, int ey_below, int ey_above) const;
290 
306  TH2* CreateH2(const char* name, const char* title, unsigned x, unsigned y, unsigned nbins_x = 100, unsigned nbins_y = 100, double x_padding = 0.10, double y_padding = 0.10) const;
307 
310 private:
313  std::vector<BCDataPoint> fDataVector;
314 
316  unsigned fNValuesPerPoint;
317 
319  BCDataPoint fLowerBounds;
320 
322  BCDataPoint fUpperBounds;
323 
325  BCDataPoint fUserLowerBounds;
326 
328  BCDataPoint fUserUpperBounds;
329 
331  std::vector<bool> fFixed;
332 
333 };
334 
335 // ---------------------------------------------------------
336 
337 #endif
const BCDataPoint & GetUpperBounds() const
Definition: BCDataSet.h:120
void SetBounds(unsigned index, double lower_bound, double upper_bound, bool fixed=false)
Set bounds for data values.
Definition: BCDataSet.cxx:267
bool AddDataPoint(const BCDataPoint &datapoint)
Adds a data point to the data set.
Definition: BCDataSet.cxx:209
bool ReadDataFromFileTxt(const std::string &filename, int nvariables)
Reads data from a .txt file.
Definition: BCDataSet.cxx:157
A class representing a data point.
Definition: BCDataPoint.h:34
const BCDataPoint & operator[](unsigned index) const
Raw and fast access.
Definition: BCDataSet.h:66
std::vector< double > GetDataComponents(unsigned index) const
Viewing the data set as a table with one row per point, this method returns a specified column...
Definition: BCDataSet.cxx:37
bool ReadDataFromFileTree(const std::string &filename, const std::string &treename, const std::string &branchnames, char delim= ',')
Reads a TTree from a .root file.
Definition: BCDataSet.cxx:88
A class representing a set of data points.
Definition: BCDataSet.h:39
const BCDataPoint & GetDataPoint(unsigned index) const
Safer, but slower, access to data points.
Definition: BCDataSet.h:94
TH2 * CreateH2(const char *name, const char *title, unsigned x, unsigned y, unsigned nbins_x=100, unsigned nbins_y=100, double x_padding=0.10, double y_padding=0.10) const
Get ROOT TH2 with ranges set to data bounds.
Definition: BCDataSet.cxx:353
unsigned GetNDataPoints() const
Definition: BCDataSet.h:75
void SetNValuesPerPoint(unsigned n)
Set number of values inside each data point.
Definition: BCDataSet.cxx:256
double GetLowerBound(unsigned index) const
Return user-set lower bound on data, if set, otherwise actual lower bound.
Definition: BCDataSet.cxx:64
void Reset()
Resets the content of the data set.
Definition: BCDataSet.h:245
bool ReadDataFromFile(const std::string &filename, int nvariables)
Reads data from a file.
Definition: BCDataSet.h:205
void AdjustBoundForUncertainties(unsigned i, double nSigma, unsigned i_err1, int i_err2=-1)
Recalculate a data axis bound accounting for uncertainties specified by other data axes...
Definition: BCDataSet.cxx:232
BCDataPoint & Back()
Access to last added data point.
Definition: BCDataSet.h:99
double GetRangeWidth(unsigned index) const
Return upper-bound minus lower-bound for data axis, using user-set bounds, if provided, other actual bounds.
Definition: BCDataSet.h:149
bool BoundsExist() const
Definition: BCDataSet.cxx:55
BCDataPoint & GetDataPoint(unsigned index)
Safer, but slower, access to data points.
Definition: BCDataSet.h:87
void Fix(unsigned i, bool b=true)
Set fixed flag of a data axis.
Definition: BCDataSet.h:182
void PrintSummary(void(*output)(const std::string &)=BCLog::OutSummary) const
Print summary to string handler.
Definition: BCDataSet.cxx:283
BCDataPoint & operator[](unsigned index)
Raw and fast access.
Definition: BCDataSet.h:61
BCDataSet(unsigned n=0)
Default constructor.
Definition: BCDataSet.cxx:31
BCDataPoint & GetUserLowerBounds()
Definition: BCDataSet.h:125
bool ReadDataFromFile(const std::string &filename, const std::string &treename, const std::string &branchnames, char delim= ',')
Reads data from a TTree in file.
Definition: BCDataSet.h:197
double GetUpperBound(unsigned index) const
Return user-set upper bound on data, if set, otherwise actual upper bound.
Definition: BCDataSet.cxx:76
TGraph * GetGraph(unsigned x, unsigned y) const
Get data set as ROOT TGraph object,.
Definition: BCDataSet.cxx:295
const BCDataPoint & GetLowerBounds() const
Definition: BCDataSet.h:115
BCDataPoint & GetUserUpperBounds()
Definition: BCDataSet.h:130
bool IsFixed(unsigned index) const
Return wether data axis is fixed.
Definition: BCDataSet.h:156
virtual ~BCDataSet()
Destructor.
Definition: BCDataSet.h:53
unsigned GetNValuesPerPoint() const
Definition: BCDataSet.h:80