MINT2
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
HyperHistogram Class Reference

#include <HyperHistogram.h>

Inheritance diagram for HyperHistogram:
HistogramBase HyperFunction

Public Member Functions

 HyperHistogram (const HyperCuboid &binningRange, const HyperPointSet &points, HyperBinningAlgorithms::Alg alg=HyperBinningAlgorithms::MINT, AlgOption opt0=AlgOption::Empty(), AlgOption opt1=AlgOption::Empty(), AlgOption opt2=AlgOption::Empty(), AlgOption opt3=AlgOption::Empty(), AlgOption opt4=AlgOption::Empty(), AlgOption opt5=AlgOption::Empty(), AlgOption opt6=AlgOption::Empty(), AlgOption opt7=AlgOption::Empty(), AlgOption opt8=AlgOption::Empty(), AlgOption opt9=AlgOption::Empty())
 
 HyperHistogram (const BinningBase &binning)
 
 HyperHistogram (TString filename, TString option="MEMRES READ")
 
 HyperHistogram (std::vector< TString > filename)
 
 HyperHistogram (TString targetFilename, std::vector< TString > filename)
 
 HyperHistogram (const HyperHistogram &other)
 
HyperHistogramoperator= (const HyperHistogram &other)
 
int getDimension () const
 
void setNames (HyperName names)
 
HyperName getNames () const
 
int fill (const HyperPoint &coords, double weight)
 
int fill (const HyperPoint &coords)
 
void fill (const HyperPointSet &points)
 
virtual void merge (const HistogramBase &other)
 
void merge (TString filenameother)
 
int estimateCapacity (std::vector< TString > filename, TString binningType)
 
void project (TH1D *histogram, const HyperCuboid &cuboid, double content, int dimension) const
 
void project (TH1D *histogram, const HyperVolume &hyperVolume, double content, int dimension) const
 
TH1D project (int dim=0, int bins=100, TString name="projection") const
 
void drawProjection (TString path, int dim=0, int bins=100) const
 
void drawAllProjections (TString path, int bins) const
 
void compareProjection (TString path, int dim, const HyperHistogram &other, int bins=100) const
 
void compareAllProjections (TString path, const HyperHistogram &other, int bins=100) const
 
HyperHistogram slice (std::vector< int > sliceDims, const HyperPoint &slicePoint) const
 
std::vector< HyperHistogramslice (std::vector< int > sliceDims, const HyperPointSet &slicePoints) const
 
void draw2DSlice (TString path, int sliceDimX, int sliceDimY, const HyperPoint &slicePoint, TString options="") const
 
void draw2DSliceSet (TString path, int sliceDimX, int sliceDimY, int sliceSetDim, int nSlices, const HyperPoint &slicePoint, TString options="") const
 
void draw2DSliceSet (TString path, int sliceDimX, int sliceDimY, int nSlices, const HyperPoint &slicePoint, TString options="") const
 
void draw2DSliceSet (TString path, int nSlices, const HyperPoint &slicePoint, TString options="") const
 
void drawRandom2DSlice (TString path, TRandom *random=gRandom, TString options="") const
 
HyperCuboid getLimits () const
 
const BinningBasegetBinning () const
 
virtual double getVal (const HyperPoint &point) const
 
std::vector< double > getVal (const HyperPointSet &points) const
 
virtual double getBinVolume (int bin) const
 
void save (TString filename)
 
TString getBinningType (TString filename)
 
void load (TString filename, TString option="MEMRES READ")
 
void loadEmpty (TString filename, TString option="MEMRES READ", TString binningType="HyperBinning")
 
void setContentsFromFunc (const HyperFunction &func)
 
void printFull () const
 
void saveToTxtFile (TString filename, bool incError=true) const
 
void saveToTxtFileNoLinks (TString filename, bool incError=true) const
 
void draw (TString path, TString options="")
 
void drawDensity (TString path, TString options="")
 
void mergeBinsWithSameContent ()
 
virtual ~HyperHistogram ()
 
- Public Member Functions inherited from HistogramBase
 HistogramBase (int nBins)
 
virtual ~HistogramBase ()
 
int checkBinNumber (int bin) const
 
void resetBinContents (int nBins)
 
void clear ()
 
void fillBase (int binNum, double weight)
 
void setBinContent (int bin, double val)
 
void setBinError (int bin, double val)
 
double getBinContent (int bin) const
 
double getBinError (int bin) const
 
int getNBins () const
 
void divide (const HistogramBase &other)
 
void multiply (const HistogramBase &other)
 
void add (const HistogramBase &other)
 
void minus (const HistogramBase &other)
 
void pulls (const HistogramBase &other)
 
void pulls (const HistogramBase &other1, const HistogramBase &other2)
 
void asymmetry (const HistogramBase &other)
 
void asymmetry (const HistogramBase &other1, const HistogramBase &other2)
 
void drawPullHistogram (const HistogramBase &other, TString name, int nBins=50, double pmLimits=3.5) const
 
double chi2 (const HistogramBase &other) const
 
double pvalue (const HistogramBase &other, int ndof=-1) const
 
double chi2sig (const HistogramBase &other, int ndof=-1) const
 
double integral () const
 
double integralError () const
 
void randomiseWithinErrors (int seed)
 
double getMin () const
 
double getMax () const
 
void setMin (double min)
 
void setMax (double max)
 
double getMinDensity () const
 
double getMaxDensity () const
 
void setMinDensity (double min)
 
void setMaxDensity (double max)
 
void saveBase ()
 
void saveBase (TString filename)
 
void loadBase (TString filename)
 
void normalise (double area=1.0)
 
double getFrequencyDensity (int bin) const
 
void reserveCapacity (int nElements)
 
void makeFrequencyDensity ()
 
void print ()
 
- Public Member Functions inherited from HyperFunction
 HyperFunction ()
 
 HyperFunction (const HyperCuboid &limits)
 
void reweightDataset (HyperPointSet &points)
 
void setFuncLimits (const HyperCuboid &limits)
 
const HyperCuboidgetFuncLimits () const
 
TH2D make2DFuncSlice (TString name, int sliceDimX, int sliceDimY, const HyperPoint &slicePoint, int nbins=100) const
 
void draw2DFuncSlice (TString path, int sliceDimX, int sliceDimY, const HyperPoint &slicePoint, int nbins=100) const
 
void draw2DFuncSliceSet (TString path, int sliceDimX, int sliceDimY, int sliceSetDim, int nSlices, const HyperPoint &slicePoint, int nbins=100) const
 
void draw2DFuncSliceSet (TString path, int sliceDimX, int sliceDimY, int nSlices, const HyperPoint &slicePoint, int nbins=100) const
 
void draw2DFuncSliceSet (TString path, int nSlices, const HyperPoint &slicePoint, int nbins=100) const
 
double getDifference (const HyperFunction &other, const HyperPoint &point)
 
void fillCorrelations (TH2D &hist, const HyperFunction &other, const HyperPointSet &points)
 
virtual ~HyperFunction ()
 

Protected Member Functions

 HyperHistogram ()
 

Protected Attributes

BinningBase_binning
 
- Protected Attributes inherited from HistogramBase
int _nBins
 
std::vector< double > _binContents
 
std::vector< double > _sumW2
 
double _min
 
double _max
 
double _minDensity
 
double _maxDensity
 

Detailed Description

HyperPlot, Author: Sam Harnew, sam.h.nosp@m.arne.nosp@m.w@gma.nosp@m.il.c.nosp@m.om , Date: Dec 2015

Binning Algorithms:

Binning Algorithm Options:

Definition at line 53 of file HyperHistogram.h.

Constructor & Destructor Documentation

◆ HyperHistogram() [1/7]

HyperHistogram::HyperHistogram ( )
protected

Private constructor

Definition at line 238 of file HyperHistogram.cpp.

238  :
239  HistogramBase(0),
240  _binning(0)
241 {
242  WELCOME_LOG << "Good day from the HyperHistogram() Constructor";
243 }
HistogramBase(int nBins)
BinningBase * _binning
#define WELCOME_LOG

◆ HyperHistogram() [2/7]

HyperHistogram::HyperHistogram ( const HyperCuboid binningRange,
const HyperPointSet points,
HyperBinningAlgorithms::Alg  alg = HyperBinningAlgorithms::MINT,
AlgOption  opt0 = AlgOption::Empty(),
AlgOption  opt1 = AlgOption::Empty(),
AlgOption  opt2 = AlgOption::Empty(),
AlgOption  opt3 = AlgOption::Empty(),
AlgOption  opt4 = AlgOption::Empty(),
AlgOption  opt5 = AlgOption::Empty(),
AlgOption  opt6 = AlgOption::Empty(),
AlgOption  opt7 = AlgOption::Empty(),
AlgOption  opt8 = AlgOption::Empty(),
AlgOption  opt9 = AlgOption::Empty() 
)

Constuctor that adaptively bins the HyperPointSet provided, within the limits provided, and using the specified binning algorithm. Additionally, binning options can be selected which are passed to the binning algorithm.

Binning Algorithms:

Binning Algorithm Options:

Definition at line 48 of file HyperHistogram.cpp.

62  :
63  HistogramBase(0),
64  HyperFunction(binningRange),
65  _binning(0)
66 {
67 
68  HyperBinningAlgorithms algSetup(alg);
69  algSetup.addAlgOption(opt0);
70  algSetup.addAlgOption(opt1);
71  algSetup.addAlgOption(opt2);
72  algSetup.addAlgOption(opt3);
73  algSetup.addAlgOption(opt4);
74  algSetup.addAlgOption(opt5);
75  algSetup.addAlgOption(opt6);
76  algSetup.addAlgOption(opt7);
77  algSetup.addAlgOption(opt8);
78  algSetup.addAlgOption(opt9);
79 
80  HyperBinningMaker* binnningMaker = algSetup.getHyperBinningMaker(binningRange, points);
81  binnningMaker->makeBinning();
82 
83  HyperHistogram* hist = binnningMaker->getHyperBinningHistogram();
84 
85  *this = *hist;
86  delete hist;
87  delete binnningMaker;
88 
89  //This inherets from a HyperFunction. Although non-essential, it's useful for
90  //the function to have some limits for it's domain.
92 
93 }
HistogramBase(int nBins)
HyperHistogram * getHyperBinningHistogram() const
BinningBase * _binning
virtual void makeBinning()=0
HyperCuboid getLimits() const
void setFuncLimits(const HyperCuboid &limits)

◆ HyperHistogram() [3/7]

HyperHistogram::HyperHistogram ( const BinningBase binning)

The most basic constructor - just pass anything that is derrived from BinningBase

Definition at line 9 of file HyperHistogram.cpp.

9  :
10  HistogramBase(binning.getNumBins()),
11  _binning(binning.clone())
12 {
13  WELCOME_LOG << "Good day from the HyperHistogram() Constructor";
15 }
virtual BinningBase * clone() const =0
HistogramBase(int nBins)
virtual int getNumBins() const =0
BinningBase * _binning
HyperCuboid getLimits() const
#define WELCOME_LOG
void setFuncLimits(const HyperCuboid &limits)

◆ HyperHistogram() [4/7]

HyperHistogram::HyperHistogram ( TString  filename,
TString  option = "MEMRES READ" 
)

Load a HyperHistogram from filename. If

Definition at line 109 of file HyperHistogram.cpp.

109  :
110  HistogramBase(0),
111  _binning(0)
112 {
113  WELCOME_LOG << "Good day from the HyperHistogram() Constructor";
114 
115  if (option.Contains("Empty")){
116  loadEmpty(filename, option);
117  }
118  else{
119  load(filename, option);
120  }
121 
122  //This inherets from a HyperFunction. Although non-essential, it's useful for
123  //the function to have some limits for it's domain.
125 }
void loadEmpty(TString filename, TString option="MEMRES READ", TString binningType="HyperBinning")
HistogramBase(int nBins)
BinningBase * _binning
HyperCuboid getLimits() const
#define WELCOME_LOG
void setFuncLimits(const HyperCuboid &limits)
void load(TString filename, TString option="MEMRES READ")

◆ HyperHistogram() [5/7]

HyperHistogram::HyperHistogram ( std::vector< TString >  filename)

Load an array of HyperHistograms from different files and merge them into a memory resident HyperBinning

Definition at line 131 of file HyperHistogram.cpp.

131  :
132  HistogramBase(0),
133  _binning(0)
134 {
135  WELCOME_LOG << "Good day from the HyperHistogram() Constructor";
136 
137  int nFiles = filename.size();
138  if (nFiles == 0){
139  ERROR_LOG << "The list of filenames you provided to HyperHistogram is empty" << std::endl;
140  }
141 
142  INFO_LOG << "Loading HyperHistogram at: " << filename.at(0) << std::endl;
143  load(filename.at(0), "MEMRES");
144 
145  for (int i = 1; i < nFiles; i++){
146  INFO_LOG << "Loading and merging HyperHistogram at: " << filename.at(i) << std::endl;
147  merge(filename.at(i));
148  }
149 
150  //This inherets from a HyperFunction. Although non-essential, it's useful for
151  //the function to have some limits for it's domain.
153 
154 }
#define INFO_LOG
virtual void merge(const HistogramBase &other)
#define ERROR_LOG
HistogramBase(int nBins)
BinningBase * _binning
HyperCuboid getLimits() const
#define WELCOME_LOG
void setFuncLimits(const HyperCuboid &limits)
void load(TString filename, TString option="MEMRES READ")

◆ HyperHistogram() [6/7]

HyperHistogram::HyperHistogram ( TString  targetFilename,
std::vector< TString >  filename 
)

Load an array of HyperHistograms from different files and merge them into a disk resident HyperBinning stored at 'targetFilename'

Definition at line 160 of file HyperHistogram.cpp.

160  :
161  HistogramBase(0),
162  _binning(0)
163 {
164 
165  WELCOME_LOG << "Good day from the HyperHistogram() Constructor";
166 
167  int nFiles = filename.size();
168  if (nFiles == 0){
169  ERROR_LOG << "The list of filenames you provided to HyperHistogram is empty" << std::endl;
170  }
171 
172  //Getting the binning type from the first file
173  TString binningType = getBinningType( filename.at(0) );
174 
175  INFO_LOG << "Creating HyperHistogram at: " << targetFilename << " with binning type " << binningType <<std::endl;
176  loadEmpty(targetFilename, "DISKRES", binningType );
177 
178  int nResBins = estimateCapacity(filename, binningType);
179  INFO_LOG << "I estimate there will be " << nResBins << " in total - resizing the Histogram accordingly" << std::endl;
180 
181  for (int i = 0; i < nFiles; i++){
182  INFO_LOG << "Loading and merging HyperHistogram at: " << filename.at(i) << std::endl;
183  merge(filename.at(i));
184  }
185 
186  //This inherets from a HyperFunction. Although non-essential, it's useful for
187  //the function to have some limits for it's domain.
188  setFuncLimits( getLimits() );
189 }
void loadEmpty(TString filename, TString option="MEMRES READ", TString binningType="HyperBinning")
#define INFO_LOG
virtual void merge(const HistogramBase &other)
#define ERROR_LOG
HistogramBase(int nBins)
BinningBase * _binning
HyperCuboid getLimits() const
int estimateCapacity(std::vector< TString > filename, TString binningType)
#define WELCOME_LOG
void setFuncLimits(const HyperCuboid &limits)
TString getBinningType(TString filename)

◆ HyperHistogram() [7/7]

HyperHistogram::HyperHistogram ( const HyperHistogram other)

Definition at line 226 of file HyperHistogram.cpp.

226  :
227  HistogramBase(other),
228  HyperFunction(other),
229  _binning( other._binning->clone() )
230 {
231 
232 }
virtual BinningBase * clone() const =0
HistogramBase(int nBins)
BinningBase * _binning

◆ ~HyperHistogram()

HyperHistogram::~HyperHistogram ( )
virtual

Destructor

Definition at line 1269 of file HyperHistogram.cpp.

1269  {
1270 
1271  bool diskResidentBinning = 0;
1272  TString filename = "";
1273 
1274  if (_binning != 0){
1275  diskResidentBinning = _binning->isDiskResident();
1276  filename = _binning->filename();
1277  delete _binning;
1278  _binning = 0;
1279  }
1280 
1281  if (diskResidentBinning){
1282 
1283  TFile *file = new TFile(filename, "update");
1284  std::string object_to_remove="HistogramBase";
1285  gDirectory->Delete(object_to_remove.c_str());
1286 
1287  saveBase();
1288 
1289  file->Close();
1290  }
1291 
1292  GOODBYE_LOG << "Goodbye from the HyperHistogram() Constructor";
1293 }
virtual bool isDiskResident() const
Definition: BinningBase.cpp:54
virtual TString filename() const
Definition: BinningBase.cpp:57
BinningBase * _binning
#define GOODBYE_LOG

Member Function Documentation

◆ compareAllProjections()

void HyperHistogram::compareAllProjections ( TString  path,
const HyperHistogram other,
int  bins = 100 
) const
Todo:
remember how this works

Definition at line 1048 of file HyperHistogram.cpp.

1048  {
1049  for(int i = 0; i < _binning->getDimension(); i++){
1050  TString thisPath = path + "_"; thisPath += i;
1051  compareProjection(thisPath, i, other, bins);
1052  }
1053 }
BinningBase * _binning
const int & getDimension() const
Definition: BinningBase.cpp:30
void compareProjection(TString path, int dim, const HyperHistogram &other, int bins=100) const

◆ compareProjection()

void HyperHistogram::compareProjection ( TString  path,
int  dim,
const HyperHistogram other,
int  bins = 100 
) const
Todo:
remember how this works

Definition at line 1036 of file HyperHistogram.cpp.

1036  {
1037  TH1D projection = project(dim, bins);
1038  TH1D projectionOther = other.project(dim, bins, "projection2");
1039  RootPlotter1D plotter(&projection, 300, 300);
1040  plotter.add(&projectionOther);
1041  plotter.setMin(0.0);
1042  plotter.plotWithRatio(path);
1043 }
void project(TH1D *histogram, const HyperCuboid &cuboid, double content, int dimension) const

◆ draw()

void HyperHistogram::draw ( TString  path,
TString  options = "" 
)

Draw the HyperHistogram - the drawing class used depends on the dimensionality of the data. This just plots the raw bin contents, not the frequency density.

Definition at line 594 of file HyperHistogram.cpp.

594  {
595 
596  if (_binning->getDimension() == 1){
597  HyperBinningPainter1D painter(this);
598  painter.draw(path, options);
599  }
600  else if (_binning->getDimension() == 2){
601 
602  HyperBinningPainter2D painter(this);
603  painter.draw(path, options);
604  }
605  else{
606  HyperBinningPainter painter(this);
607  painter.draw(path, options);
608  }
609 
610 }
BinningBase * _binning
const int & getDimension() const
Definition: BinningBase.cpp:30

◆ draw2DSlice()

void HyperHistogram::draw2DSlice ( TString  path,
int  sliceDimX,
int  sliceDimY,
const HyperPoint slicePoint,
TString  options = "" 
) const

Definition at line 825 of file HyperHistogram.cpp.

825  {
826 
827  std::vector<int > sliceDims;
828 
829  for (int i = 0; i < slicePoint.getDimension(); i++){
830  if (i == sliceDimX) continue;
831  if (i == sliceDimY) continue;
832  sliceDims.push_back( i );
833  }
834 
835  HyperHistogram sliceHist = slice( sliceDims, slicePoint );
836  sliceHist.draw(path, options);
837 
838 }
void draw(TString path, TString options="")
int getDimension() const
Definition: HyperPoint.h:99
HyperHistogram slice(std::vector< int > sliceDims, const HyperPoint &slicePoint) const

◆ draw2DSliceSet() [1/3]

void HyperHistogram::draw2DSliceSet ( TString  path,
int  sliceDimX,
int  sliceDimY,
int  sliceSetDim,
int  nSlices,
const HyperPoint slicePoint,
TString  options = "" 
) const

Definition at line 861 of file HyperHistogram.cpp.

861  {
862 
863  std::vector<int > _sliceDims;
864 
865  for (int i = 0; i < slicePoint.getDimension(); i++){
866  if (i == sliceDimX) continue;
867  if (i == sliceDimY) continue;
868  _sliceDims.push_back( i );
869  }
870 
871  HyperPointSet slicePoints(getDimension());
872 
873 
874  HyperPoint slicePointCp(slicePoint);
875 
876  double min = _binning->getMin(sliceSetDim);
877  double max = _binning->getMax(sliceSetDim);
878  double width = (max - min)/double(nSlices);
879 
880  std::vector<TString> paths;
881 
882  for (int i = 0; i < nSlices; i++){
883  double val = min + width*(i + 0.5);
884  slicePointCp.at(sliceSetDim) = val;
885 
886  slicePoints.push_back(slicePointCp);
887 
888  TString uniquePath = path;
889  uniquePath += "_sliceNum";
890  uniquePath += i;
891 
892  paths.push_back(uniquePath);
893 
894  }
895 
896  std::vector<HyperHistogram> hists = slice(_sliceDims, slicePoints);
897 
898  for (unsigned i = 0; i < hists.size(); i++){
899  hists.at(i).draw(paths.at(i), options);
900  }
901 
902 }
double getMax(int dimension) const
Definition: BinningBase.cpp:46
int getDimension() const
double getMin(int dimension) const
Definition: BinningBase.cpp:42
BinningBase * _binning
int getDimension() const
Definition: HyperPoint.h:99
HyperHistogram slice(std::vector< int > sliceDims, const HyperPoint &slicePoint) const

◆ draw2DSliceSet() [2/3]

void HyperHistogram::draw2DSliceSet ( TString  path,
int  sliceDimX,
int  sliceDimY,
int  nSlices,
const HyperPoint slicePoint,
TString  options = "" 
) const

Definition at line 904 of file HyperHistogram.cpp.

904  {
905 
906  //Get the slice dimesnions from the given dimensions
907  std::vector<int > sliceDims;
908 
909  for (int i = 0; i < slicePoint.getDimension(); i++){
910  if (i == sliceDimX) continue;
911  if (i == sliceDimY) continue;
912  sliceDims.push_back( i );
913  }
914 
915  //Fill a HyperPointSet with slice points, and create a name for each in paths
916  HyperPointSet slicePoints(getDimension());
917  std::vector<TString> paths;
918 
919  //Loop over each slice dimesion - will take nSlices along each of these
920  for (unsigned j = 0; j < sliceDims.size(); j++){
921 
922  HyperPoint slicePointCp(slicePoint);
923 
924  int sliceSetDim = sliceDims.at(j);
925 
926  double min = _binning->getMin(sliceSetDim);
927  double max = _binning->getMax(sliceSetDim);
928  double width = (max - min)/double(nSlices);
929 
930  TString pathj = path;
931  pathj += "_scanDim";
932  pathj += sliceSetDim;
933 
934  for (int i = 0; i < nSlices; i++){
935  double val = min + width*(i + 0.5);
936  slicePointCp.at(sliceSetDim) = val;
937 
938  slicePoints.push_back(slicePointCp);
939 
940  TString uniquePath = pathj;
941  uniquePath += "_sliceNum";
942  uniquePath += i;
943 
944  paths.push_back(uniquePath);
945 
946  }
947 
948  }
949 
950 
951 
952  std::vector<HyperHistogram> hists = slice(sliceDims, slicePoints);
953 
954  for (unsigned i = 0; i < hists.size(); i++){
955  hists.at(i).draw(paths.at(i), options);
956  }
957 
958 
959 
960 
961 }
double getMax(int dimension) const
Definition: BinningBase.cpp:46
int getDimension() const
double getMin(int dimension) const
Definition: BinningBase.cpp:42
BinningBase * _binning
int getDimension() const
Definition: HyperPoint.h:99
HyperHistogram slice(std::vector< int > sliceDims, const HyperPoint &slicePoint) const

◆ draw2DSliceSet() [3/3]

void HyperHistogram::draw2DSliceSet ( TString  path,
int  nSlices,
const HyperPoint slicePoint,
TString  options = "" 
) const

Definition at line 963 of file HyperHistogram.cpp.

963  {
964 
965 
966 
967  for (int i = 0; i < slicePoint.getDimension(); i++){
968  for (int j = 0; j < slicePoint.getDimension(); j++){
969 
970  if (i >= j) continue;
971 
972  TString thsPath = path;
973  thsPath += "_";
974  thsPath += i;
975  thsPath += "vs";
976  thsPath += j;
977 
978  draw2DSliceSet(thsPath, i, j, nSlices, slicePoint, options);
979  }
980  }
981 
982 
983 }
void draw2DSliceSet(TString path, int sliceDimX, int sliceDimY, int sliceSetDim, int nSlices, const HyperPoint &slicePoint, TString options="") const
int getDimension() const
Definition: HyperPoint.h:99

◆ drawAllProjections()

void HyperHistogram::drawAllProjections ( TString  path,
int  bins 
) const
Todo:
remember how this works

Definition at line 1024 of file HyperHistogram.cpp.

1024  {
1025 
1026  for(int i = 0; i < _binning->getDimension(); i++){
1027  TString thisPath = path + "_"; thisPath += i;
1028  drawProjection(thisPath, i, bins);
1029  }
1030 
1031 }
void drawProjection(TString path, int dim=0, int bins=100) const
BinningBase * _binning
const int & getDimension() const
Definition: BinningBase.cpp:30

◆ drawDensity()

void HyperHistogram::drawDensity ( TString  path,
TString  options = "" 
)

Draw the frequency density of the HyperHistogram

  • the drawing class used depends on the dimensionality of the data.

Definition at line 616 of file HyperHistogram.cpp.

616  {
617 
618  if (_binning->getDimension() == 1){
619  HyperBinningPainter1D painter(this);
620  painter.useDensity(true);
621  painter.draw(path, options);
622  }
623  else if (_binning->getDimension() == 2){
624  HyperBinningPainter2D painter(this);
625  painter.useDensity(true);
626  painter.draw(path, options);
627  }
628  else{
629  HyperBinningPainter painter(this);
630  painter.useDensity(true);
631  painter.draw(path, options);
632  }
633 
634 }
BinningBase * _binning
const int & getDimension() const
Definition: BinningBase.cpp:30

◆ drawProjection()

void HyperHistogram::drawProjection ( TString  path,
int  dim = 0,
int  bins = 100 
) const
Todo:
remember how this works

Definition at line 1012 of file HyperHistogram.cpp.

1012  {
1013 
1014  TH1D projection = project(dim, bins);
1015  RootPlotter1D plotter(&projection, 300, 300);
1016  plotter.setMin(0.0);
1017  plotter.plot(path);
1018 
1019 }
void project(TH1D *histogram, const HyperCuboid &cuboid, double content, int dimension) const

◆ drawRandom2DSlice()

void HyperHistogram::drawRandom2DSlice ( TString  path,
TRandom *  random = gRandom,
TString  options = "" 
) const

Definition at line 840 of file HyperHistogram.cpp.

840  {
841 
842  int dim = getDimension();
843 
844  if (dim < 3){
845  ERROR_LOG << "Why would you take a 2D slice of something with less than 3 dim." << std::endl;
846  return;
847  }
848 
849  int slicedimx = random->Integer(dim);
850  int slicedimy = random->Integer(dim);
851  while( slicedimx == slicedimy ){
852  slicedimy = random->Integer(dim);
853  }
854 
855  HyperPoint slicepoint = getLimits().getRandomPoint(random);
856 
857  draw2DSlice(path, slicedimx, slicedimy, slicepoint, options);
858 
859 }
int getDimension() const
HyperPoint getRandomPoint(TRandom *random=gRandom) const
#define ERROR_LOG
void draw2DSlice(TString path, int sliceDimX, int sliceDimY, const HyperPoint &slicePoint, TString options="") const
HyperCuboid getLimits() const

◆ estimateCapacity()

int HyperHistogram::estimateCapacity ( std::vector< TString >  filename,
TString  binningType 
)

Definition at line 192 of file HyperHistogram.cpp.

192  {
193 
194  int nbins = 0;
195  int nvols = 0;
196 
197  for (unsigned i = 0; i < filename.size(); i++){
198  TFile* file = new TFile(filename.at(i), "READ");
199  if (file == 0){
200  ERROR_LOG << "HyperHistogram::estimateCapacity - " << filename.at(i) << " does not exist" << std::endl;
201  return 0;
202  }
203  TTree* hist = dynamic_cast<TTree*>( file->Get(binningType) );
204  TTree* base = dynamic_cast<TTree*>( file->Get("HistogramBase") );
205  if (hist == 0){
206  ERROR_LOG << "HyperHistogram::estimateCapacity - " << filename.at(i) << " does not contain tree " << binningType << std::endl;
207  return 0;
208  }
209  if (base == 0){
210  ERROR_LOG << "HyperHistogram::estimateCapacity - " << filename.at(i) << " does not contain tree HistogramBase" << std::endl;
211  return 0;
212  }
213  nbins += base->GetEntries();
214  nvols += hist->GetEntries();
215  file->Close();
216  }
217 
218  reserveCapacity(nbins);
219  _binning->reserveCapacity(nvols);
220 
221  return nbins;
222 
223 }
#define ERROR_LOG
BinningBase * _binning
void reserveCapacity(int nElements)
virtual void reserveCapacity(int nElements)
Definition: BinningBase.cpp:61

◆ fill() [1/3]

int HyperHistogram::fill ( const HyperPoint coords,
double  weight 
)

Fill the HyperHistogram with a HyperPoint and aspecified weight

Definition at line 262 of file HyperHistogram.cpp.

262  {
263 
264  int binNumber = _binning->getBinNum(coords);
265  this->fillBase(binNumber, weight);
266  return binNumber;
267 }
void fillBase(int binNum, double weight)
virtual int getBinNum(const HyperPoint &coords) const =0
BinningBase * _binning

◆ fill() [2/3]

int HyperHistogram::fill ( const HyperPoint coords)

Fill the HyperHistogram with a HyperPoint. If the HyperPoint has a weight, use it.

Definition at line 273 of file HyperHistogram.cpp.

273  {
274 
275  int binNumber = _binning->getBinNum(coords);
276  this->fillBase(binNumber, coords.getWeight(0));
277  return binNumber;
278 }
void fillBase(int binNum, double weight)
virtual int getBinNum(const HyperPoint &coords) const =0
BinningBase * _binning
double getWeight(int i=0) const
Definition: Weights.cpp:40

◆ fill() [3/3]

void HyperHistogram::fill ( const HyperPointSet points)

Add a HyperPointSet to the HyperHistogram - if any of the HyperPoints are weighted, they will be used.

Definition at line 315 of file HyperHistogram.cpp.

315  {
316 
317  for(unsigned i = 0; i < points.size(); i++){
318  fill(points.at(i), points.at(i).getWeight());
319  }
320 
321 }
const HyperPoint & at(int i) const
int fill(const HyperPoint &coords, double weight)
unsigned int size() const
double getWeight(int i=0) const
Definition: Weights.cpp:40

◆ getBinning()

const BinningBase& HyperHistogram::getBinning ( ) const
inline

get the HyperVolumeBinning

Definition at line 137 of file HyperHistogram.h.

◆ getBinningType()

TString HyperHistogram::getBinningType ( TString  filename)

Get binning type from file

Definition at line 1187 of file HyperHistogram.cpp.

1187  {
1188 
1189  TFile* file = new TFile(filename, "READ");
1190 
1191  if (file == 0){
1192  ERROR_LOG << "Could not open TFile in HyperBinning::load(" << filename << ")";
1193  return "";
1194  }
1195 
1196  TTree* tree = (TTree*)file->Get("HyperBinning");
1197 
1198  if (tree != 0){
1199  file->Close();
1200  return "HyperBinning";
1201  }
1202 
1203  file->Close();
1204  return "";
1205 
1206 }
#define ERROR_LOG

◆ getBinVolume()

double HyperHistogram::getBinVolume ( int  bin) const
virtual

Get the volume of a HyperVolume bin

Reimplemented from HistogramBase.

Definition at line 1262 of file HyperHistogram.cpp.

1262  {
1263  return _binning->getBinHyperVolume(bin).volume();
1264 }
BinningBase * _binning
double volume() const
virtual HyperVolume getBinHyperVolume(int binNumber) const =0

◆ getDimension()

int HyperHistogram::getDimension ( ) const

Definition at line 705 of file HyperHistogram.cpp.

705  {
706 
707  if (_binning == 0){
708  ERROR_LOG << "HyperHistogram::getDimension - cannot get dimension, binning not set." << std::endl;
709  return 0;
710  }
711  return _binning->getDimension();
712 
713 }
#define ERROR_LOG
BinningBase * _binning
const int & getDimension() const
Definition: BinningBase.cpp:30

◆ getLimits()

HyperCuboid HyperHistogram::getLimits ( ) const

Get the limits of the histogram

Definition at line 327 of file HyperHistogram.cpp.

327  {
328  return _binning->getLimits();
329 }
BinningBase * _binning
virtual HyperCuboid getLimits() const =0

◆ getNames()

HyperName HyperHistogram::getNames ( ) const

Get the HyperName (mainly used for axis labels)

Definition at line 251 of file HyperHistogram.cpp.

◆ getVal() [1/2]

double HyperHistogram::getVal ( const HyperPoint point) const
virtual

Get the bin content where the given HyperPoint lies

Implements HyperFunction.

Definition at line 283 of file HyperHistogram.cpp.

283  {
284 
285  int binNumber = _binning->getBinNum(point);
286  return this->getBinContent(binNumber);
287 
288 }
double getBinContent(int bin) const
virtual int getBinNum(const HyperPoint &coords) const =0
BinningBase * _binning

◆ getVal() [2/2]

std::vector< double > HyperHistogram::getVal ( const HyperPointSet points) const

Get the bin content for each HyperPoint in a HyperPointSet

Definition at line 293 of file HyperHistogram.cpp.

293  {
294 
295  int nPoints = points.size();
296 
297  std::vector<int> binNums = _binning->getBinNum(points);
298 
299  std::vector<double> conts(nPoints);
300 
301  for (int i = 0; i < nPoints; i++){
302  conts.at(i) = this->getBinContent( binNums.at(i) );
303  }
304 
305  return conts;
306 
307 }
double getBinContent(int bin) const
virtual int getBinNum(const HyperPoint &coords) const =0
BinningBase * _binning
unsigned int size() const

◆ load()

void HyperHistogram::load ( TString  filename,
TString  option = "MEMRES READ" 
)

Load the HyperHistogram from a TFile

Definition at line 1211 of file HyperHistogram.cpp.

1211  {
1212 
1213  //If loading from a file, we first need to figure out what
1214  //type of binning is saved in that file.
1215 
1216  TString binningType = getBinningType(filename);
1217 
1218  if (binningType.Contains("HyperBinning")){
1219 
1220  if (_binning != 0) {
1221  delete _binning;
1222  _binning = 0;
1223  }
1224 
1225  if (option.Contains("DISK")) _binning = new HyperBinningDiskRes();
1226  else{
1227  _binning = new HyperBinningMemRes();
1228  }
1229 
1230  }
1231 
1232  if (binningType == ""){
1233  ERROR_LOG << "HyperHistogram::load - I could not find any binning scheme in this file" << std::endl;
1234  }
1235 
1236  _binning->load(filename, "READ");
1237  this->loadBase(filename);
1238 
1239 }
void loadBase(TString filename)
virtual void load(TString filename, TString option="READ")=0
#define ERROR_LOG
BinningBase * _binning
TString getBinningType(TString filename)

◆ loadEmpty()

void HyperHistogram::loadEmpty ( TString  filename,
TString  option = "MEMRES READ",
TString  binningType = "HyperBinning" 
)

Definition at line 1241 of file HyperHistogram.cpp.

1241  {
1242 
1243  if (binningType.Contains("HyperBinning")){
1244 
1245  if (option.Contains("DISK")) _binning = new HyperBinningDiskRes();
1246  else{
1247  _binning = new HyperBinningMemRes();
1248  }
1249 
1250  }
1251  if (binningType == ""){
1252  ERROR_LOG << "HyperHistogram::loadEmpty - I could not find any binning scheme in this file" << std::endl;
1253  }
1254 
1255  _binning->load(filename, "RECREATE");
1256  this->resetBinContents(0);
1257 
1258 }
virtual void load(TString filename, TString option="READ")=0
#define ERROR_LOG
void resetBinContents(int nBins)
BinningBase * _binning

◆ merge() [1/2]

void HyperHistogram::merge ( const HistogramBase other)
virtual

Merge two HyperHistograms

Reimplemented from HistogramBase.

Definition at line 335 of file HyperHistogram.cpp.

335  {
336 
337  //INFO_LOG << "Starting HyperHistogram::merge" <<std::endl;
338 
339  const HyperHistogram* histOther = dynamic_cast<const HyperHistogram*>(&other);
340 
341  if (histOther == 0){
342  ERROR_LOG << "The object passed to HyperHistogram::merge is not of type ";
343  ERROR_LOG << "HyperHistogram, so cannot merge";
344  return;
345  }
346 
347  _binning->mergeBinnings( histOther->getBinning() );
348  HistogramBase::merge( other );
349 
350  //I do this so that getLimits doesn't have to be
351  //called (which will update the cache). Will speed
352  //things up I hope
353 
355  HyperCuboid b = histOther->getFuncLimits();
356 
358  if (a.getDimension() != 0) vol.push_back(a);
359  if (b.getDimension() != 0) vol.push_back(b);
360 
361  setFuncLimits( vol.getLimits() );
362 
363  //INFO_LOG << "Ending HyperHistogram::merge" <<std::endl;
364 
365 
366 }
virtual void mergeBinnings(const BinningBase &other)=0
virtual void merge(const HistogramBase &other)
void push_back(const HyperCuboid &hyperCuboid)
#define ERROR_LOG
const int & getDimension() const
Definition: HyperCuboid.h:45
BinningBase * _binning
const int & getDimension() const
Definition: BinningBase.cpp:30
const BinningBase & getBinning() const
void setFuncLimits(const HyperCuboid &limits)
const HyperCuboid & getFuncLimits() const
Definition: HyperFunction.h:44

◆ merge() [2/2]

void HyperHistogram::merge ( TString  filenameother)

Merge this histogram with another in a file

Definition at line 371 of file HyperHistogram.cpp.

371  {
372 
373  HyperHistogram other(filenameother, "DISK");
374  merge( other );
375 
376 }
virtual void merge(const HistogramBase &other)

◆ mergeBinsWithSameContent()

void HyperHistogram::mergeBinsWithSameContent ( )

Definition at line 403 of file HyperHistogram.cpp.

403  {
404 
405 
406  if ( _binning->getBinningType() != "HyperBinning" ){
407  ERROR_LOG << "It is only possible to merge bins when using HyperBinning. Doing nothing." << std::endl;
408  return;
409  }
410 
411 
412 
413  const HyperBinning& hyperBinning = dynamic_cast<const HyperBinning&>( getBinning() );
414 
415  std::map<int, bool> volumeKept;
416  for (int i = 0; i < hyperBinning.getNumHyperVolumes(); i++){
417  volumeKept[i] = true;
418  }
419 
420  //Loop over all HyperVolumes and see if there are any linked bins.
421  //If there are, see if these linked bins are actually bins, and not
422  //just part of the binning hierarchy. If they are actually bins,
423  //See if they al have the same bin content. If they do, mark them
424  //to be removed.
425 
426  for (int i = 0; i < hyperBinning.getNumHyperVolumes(); i++){
427 
428  std::vector<int> linkedVols = hyperBinning.getLinkedHyperVolumes(i);
429  if (linkedVols.size() == 0) continue;
430 
431  bool linksLeadToBins = true;
432 
433  for (unsigned j = 0; j < linkedVols.size(); j++){
434  int volNum = linkedVols.at(j);
435  if (hyperBinning.getLinkedHyperVolumes(volNum).size() != 0){
436  linksLeadToBins = false;
437  break;
438  }
439  }
440 
441  if (linksLeadToBins == false) continue;
442 
443  double binCont0 = -99999.9;
444 
445  bool binsHaveSameContent = true;
446 
447  for (unsigned j = 0; j < linkedVols.size(); j++){
448  int volNum = linkedVols.at(j);
449  int binNum = hyperBinning.getBinNum(volNum);
450  double binContent = getBinContent(binNum);
451  if (j == 0) binCont0 = binContent;
452 
453  if (binContent != binCont0){
454  binsHaveSameContent = false;
455  break;
456  }
457 
458  }
459 
460  if (binsHaveSameContent == false) continue;
461 
462  for (unsigned j = 0; j < linkedVols.size(); j++){
463  int volNum = linkedVols.at(j);
464  volumeKept[volNum] = false;
465  }
466 
467 
468  }
469 
470  //Make a map of the old volume numbers to the new ones (once the removed bins have
471  //actually been removed).
472 
473  std::map<int, int> oldToNewVolumeNum;
474  int newVolNum = 0;
475 
476  for (int i = 0; i < hyperBinning.getNumHyperVolumes(); i++){
477  bool exists = volumeKept[i];
478 
479  if (exists == true){
480  oldToNewVolumeNum[i] = newVolNum;
481  newVolNum++;
482  }
483  else{
484  oldToNewVolumeNum[i] = -1;
485  }
486 
487  }
488 
489  volumeKept.clear();
490 
491  //Create a new HyperVolumeBinning with the bins removed
492  HyperBinning* binningNew;
493 
494 
495  if (hyperBinning.isDiskResident() == true){
496  binningNew = new HyperBinningDiskRes();
497  binningNew->load( hyperBinning.filename().ReplaceAll(".root", "_temp.root"), "RECREATE" );
498  }
499  else{
500  binningNew = new HyperBinningMemRes();
501  }
502 
503  INFO_LOG << "Created a new HyperBinning for the reduced binning" << std::endl;
504 
505  int count = 0;
506 
507  for (int i = 0; i < hyperBinning.getNumHyperVolumes(); i++){
508  HyperVolume vol = hyperBinning.getHyperVolume(i);
509 
510  int newVolNum = oldToNewVolumeNum[i];
511 
512  if (newVolNum == -1) continue;
513 
514  if (newVolNum != count){
515  ERROR_LOG << "Something has gone wrong in mergeBinsWithSameContent()" << std::endl;
516  }
517 
518 
519  std::vector<int> linkedVols = hyperBinning.getLinkedHyperVolumes( i );
520  std::vector<int> newLinkedVols;
521 
522  int nLinked = 0;
523  for (unsigned j = 0; j < linkedVols.size(); j++){
524  int linkVolNum = linkedVols.at(j);
525  int newLinkVolNum = oldToNewVolumeNum[linkVolNum];
526 
527  if (newLinkVolNum != -1){
528  newLinkedVols.push_back(newLinkVolNum);
529  nLinked++;
530  }
531 
532 
533  }
534 
535  binningNew->addHyperVolume(vol, newLinkedVols);
536 
537 
538  if (nLinked == 1){
539  INFO_LOG << "This should never be one" << std::endl;
540  }
541 
542  count++;
543  }
544 
545 
546  std::vector<int> primVolNums = hyperBinning.getPrimaryVolumeNumbers();
547 
548  for (unsigned i = 0; i < primVolNums.size(); i++){
549  int oldVolNum = primVolNums.at(i);
550  int newVolNum = oldToNewVolumeNum[oldVolNum];
551  binningNew->addPrimaryVolumeNumber(newVolNum);
552  }
553 
554  oldToNewVolumeNum.clear();
555 
556  INFO_LOG << "Filled the new binning with reduced bins" << std::endl;
557 
558 
559  HyperHistogram newHist(*binningNew);
560 
561  INFO_LOG << "Made a new hitogram which will clone the binning" << std::endl;
562 
563  newHist.setContentsFromFunc(*this);
564 
565  INFO_LOG << "You have managed to remove " << hyperBinning.getNumBins() - binningNew->getNumBins() << " bins with the same content" << std::endl;
566 
567  bool moreBinsToMerge = false;
568  if ( hyperBinning.getNumBins() - binningNew->getNumBins() > 0 ) moreBinsToMerge = true;
569 
570  if (hyperBinning.isDiskResident() == true){
571  TString filenm = hyperBinning.filename();
572  delete _binning;
573  _binning = 0;
574  newHist.save( filenm );
575 
576  load(filenm, "DISK");
577  }
578  else{
579  *this = newHist;
580  }
581 
582  delete binningNew;
583 
584  if (moreBinsToMerge) this->mergeBinsWithSameContent();
585 
586 }
virtual bool isDiskResident() const
Definition: BinningBase.cpp:54
virtual TString filename() const
Definition: BinningBase.cpp:57
virtual void load(TString filename, TString option="READ")=0
#define INFO_LOG
virtual bool addHyperVolume(const HyperVolume &hyperVolume, std::vector< int > linkedVolumes=std::vector< int >(0, 0))=0
bool exists(const string &fname)
virtual HyperVolume getHyperVolume(int volumeNumber) const =0
#define ERROR_LOG
virtual std::vector< int > getLinkedHyperVolumes(int volumeNumber) const =0
double getBinContent(int bin) const
virtual int getNumHyperVolumes() const =0
TString getBinningType() const
Definition: BinningBase.cpp:50
void mergeBinsWithSameContent()
virtual void addPrimaryVolumeNumber(int volumeNumber)=0
virtual std::vector< int > getPrimaryVolumeNumbers() const
virtual void save(TString filename) const =0
BinningBase * _binning
virtual int getNumBins() const
const BinningBase & getBinning() const
int getBinNum(int volumeNumber) const
void load(TString filename, TString option="MEMRES READ")

◆ operator=()

HyperHistogram & HyperHistogram::operator= ( const HyperHistogram other)

Definition at line 95 of file HyperHistogram.cpp.

95  {
96 
97  HistogramBase::operator=(other);
98  HyperFunction::operator=(other);
99 
100  _binning = other._binning->clone();
101  return *this;
102 
103 }
virtual BinningBase * clone() const =0
BinningBase * _binning

◆ printFull()

void HyperHistogram::printFull ( ) const

Print all info about the HyperHistogram

Definition at line 641 of file HyperHistogram.cpp.

641  {
642 
643  for(int i = 0; i < _binning->getNumBins(); i++){
644  INFO_LOG << "Bin Content " << i << ": " << _binContents[i] << " SumW2: " << _sumW2[i];
646  }
647 
648  INFO_LOG << "Overflow: " << _binContents[_nBins] << std::endl;
649 
650 }
#define INFO_LOG
std::vector< double > _binContents
Definition: HistogramBase.h:35
std::vector< double > _sumW2
Definition: HistogramBase.h:36
virtual int getNumBins() const =0
BinningBase * _binning
void print(std::ostream &os=std::cout, int endline=1) const
Definition: HyperCuboid.cpp:99
const HyperCuboid & getHyperCuboid(int i) const
Definition: HyperVolume.h:54
virtual HyperVolume getBinHyperVolume(int binNumber) const =0

◆ project() [1/3]

void HyperHistogram::project ( TH1D *  histogram,
const HyperCuboid cuboid,
double  content,
int  dimension 
) const
Todo:
remember how this works

Definition at line 656 of file HyperHistogram.cpp.

656  {
657 
658  double hyperLowEdge = cuboid.getLowCorner() .at(dimension);
659  double hyperHighEdge = cuboid.getHighCorner().at(dimension);
660  double totWidth = hyperHighEdge - hyperLowEdge;
661  int lowBin = histogram->GetXaxis()->FindFixBin(hyperLowEdge);
662  int highBin = histogram->GetXaxis()->FindFixBin(hyperHighEdge);
663 
664  if (lowBin==highBin) histogram->Fill(hyperLowEdge, content);
665  else{
666 
667  //first deal with the highest and lowest bin as there will be a fractional overlap with the HyperCuboid
668 
669  double widthInLowBin = histogram->GetXaxis()->GetBinUpEdge(lowBin) - hyperLowEdge;
670  double widthInHighBin = hyperHighEdge - histogram->GetXaxis()->GetBinLowEdge(highBin);
671  double eventsInLowBin = (widthInLowBin /totWidth)*content;
672  double eventsInHighBin = (widthInHighBin/totWidth)*content;
673  histogram->Fill(hyperLowEdge , eventsInLowBin);
674  histogram->Fill(hyperHighEdge, eventsInHighBin);
675 
676  //now do the bins in the middle
677 
678  for(int bin = (lowBin + 1); bin <= (highBin - 1); bin++){
679  double lowEdge = histogram->GetXaxis()->GetBinLowEdge(bin);
680  double highEdge = histogram->GetXaxis()->GetBinUpEdge (bin);
681  double events = ((highEdge - lowEdge)/totWidth)*content;
682  histogram->Fill( histogram->GetXaxis()->GetBinCenter(bin) , events);
683  }
684 
685  }
686 
687 }
const HyperPoint & getHighCorner() const
Definition: HyperCuboid.h:89
const double & at(int i) const
Definition: HyperPoint.cpp:433
const HyperPoint & getLowCorner() const
Definition: HyperCuboid.h:87

◆ project() [2/3]

void HyperHistogram::project ( TH1D *  histogram,
const HyperVolume hyperVolume,
double  content,
int  dimension 
) const
Todo:
remember how this works

Definition at line 692 of file HyperHistogram.cpp.

692  {
693 
694  double volume = hyperVolume.volume();
695  for(int i = 0; i < hyperVolume.size(); i++){
696  const HyperCuboid& cuboid = hyperVolume.getHyperCuboid(i);
697  double cuboidVolume = cuboid.volume();
698  double cuboidContent = (content*cuboidVolume)/volume;
699  project(histogram, cuboid, cuboidContent, dimension);
700  }
701 
702 }
int size() const
Definition: HyperVolume.h:60
double volume() const
double volume() const
void project(TH1D *histogram, const HyperCuboid &cuboid, double content, int dimension) const
const HyperCuboid & getHyperCuboid(int i) const
Definition: HyperVolume.h:54

◆ project() [3/3]

TH1D HyperHistogram::project ( int  dim = 0,
int  bins = 100,
TString  name = "projection" 
) const
Todo:
remember how this works

Definition at line 989 of file HyperHistogram.cpp.

989  {
990 
991  double lowEdge = _binning->getMin(dim);
992  double highEdge = _binning->getMax(dim);
993 
994  TH1D projection(name, name, bins, lowEdge, highEdge);
995  projection.GetXaxis()->SetTitle(_binning->getNames().at(dim));
996 
997  for(int i = 0; i < _binning->getNumBins(); i++){
998  project(&projection, _binning->getBinHyperVolume(i), this->getBinContent(i), dim);
999  }
1000 
1001  for (int i = 1; i <= projection.GetNbinsX(); i++){
1002  projection.SetBinError(i, 0.0);
1003  }
1004 
1005  return projection;
1006 
1007 }
double getMax(int dimension) const
Definition: BinningBase.cpp:46
HyperName getNames() const
Definition: BinningBase.cpp:25
double getMin(int dimension) const
Definition: BinningBase.cpp:42
double getBinContent(int bin) const
virtual int getNumBins() const =0
BinningBase * _binning
const TString & at(int dim) const
Definition: HyperName.h:50
void project(TH1D *histogram, const HyperCuboid &cuboid, double content, int dimension) const
virtual HyperVolume getBinHyperVolume(int binNumber) const =0

◆ save()

void HyperHistogram::save ( TString  filename)

Save the HyperHistogram to a TFile

Definition at line 1058 of file HyperHistogram.cpp.

1058  {
1059 
1060  TFile* file = new TFile(filename, "RECREATE");
1061 
1062  if (file == 0){
1063  ERROR_LOG << "Could not open TFile in HyperHistogram::save(" << filename << ")";
1064  return;
1065  }
1066 
1067  //save the bin contents
1068  this->saveBase();
1069  //save the binning
1070  _binning->save();
1071 
1072  file->Write();
1073  file->Close();
1074 
1075 }
#define ERROR_LOG
virtual void save(TString filename) const =0
BinningBase * _binning

◆ saveToTxtFile()

void HyperHistogram::saveToTxtFile ( TString  filename,
bool  incError = true 
) const

Save the HyperHistogram to a .txt file

Definition at line 1080 of file HyperHistogram.cpp.

1080  {
1081 
1082  if ( _binning->getBinningType() != "HyperBinning" ){
1083  ERROR_LOG << "It is only possible to saveToTxtFile when using HyperBinning. Doing nothing." << std::endl;
1084  return;
1085  }
1086 
1087  const HyperBinning& hyperBinning = dynamic_cast<const HyperBinning&>( getBinning() );
1088 
1089  int nVolumes = hyperBinning.getNumHyperVolumes();
1090 
1091  std::ofstream myfile;
1092  myfile.open (filename);
1093 
1094  int widthForVolNum = ceil(log10(nVolumes)) + 2;
1095 
1096  for (int i = 0; i < nVolumes; i++ ){
1097  HyperVolume vol = hyperBinning.getHyperVolume(i);
1098  HyperCuboid cube = vol.getHyperCuboid(0);
1099  int binNumber = hyperBinning.getBinNum(i);
1100  double content = -1.0;
1101  double error = -1.0;
1102 
1103  bool isPrimary = hyperBinning.isPrimaryVolume(i);
1104  bool isBin = false;
1105 
1106  std::vector<int> linkedBins = hyperBinning.getLinkedHyperVolumes(i);
1107 
1108  if (binNumber != -1){
1109  content = getBinContent(binNumber);
1110  error = getBinError (binNumber);
1111  isBin = true;
1112  }
1113 
1114  myfile << std::setw(widthForVolNum) << std::left << i;
1115 
1116  TString binType = "";
1117  if (isPrimary ) myfile << "P";
1118  //if (!isPrimary) myfile << "L";
1119  if (isBin ) myfile << "B";
1120  if (!isBin ) myfile << "V";
1121 
1122  myfile << std::setw(4) << std::left << binType;
1123 
1124 
1125  int width = vol.getDimension()*10 + 10;
1126 
1127  myfile << std::setw(width) << std::left << cube.getLowCorner() << std::setw(width) << std::left << cube.getHighCorner();
1128 
1129  if (isBin){
1130  myfile << std::setw(10) << std::left << content;
1131  if (incError) myfile << std::setw(10) << std::left << error;
1132  }
1133 
1134  if (!isBin){
1135  for (unsigned j = 0; j < linkedBins.size(); j++){
1136  myfile << std::setw(10) << std::left << linkedBins.at(j);
1137  }
1138  }
1139 
1140  myfile << std::endl;
1141  }
1142 
1143  myfile.close();
1144 
1145 
1146 }
bool isPrimaryVolume(int volumeNumber) const
double getBinError(int bin) const
virtual HyperVolume getHyperVolume(int volumeNumber) const =0
#define ERROR_LOG
virtual std::vector< int > getLinkedHyperVolumes(int volumeNumber) const =0
double getBinContent(int bin) const
virtual int getNumHyperVolumes() const =0
TString getBinningType() const
Definition: BinningBase.cpp:50
const int & getDimension() const
Definition: HyperVolume.h:44
const HyperPoint & getHighCorner() const
Definition: HyperCuboid.h:89
BinningBase * _binning
const double & at(int i) const
Definition: HyperPoint.cpp:433
const HyperPoint & getLowCorner() const
Definition: HyperCuboid.h:87
const BinningBase & getBinning() const
int getBinNum(int volumeNumber) const
const HyperCuboid & getHyperCuboid(int i) const
Definition: HyperVolume.h:54

◆ saveToTxtFileNoLinks()

void HyperHistogram::saveToTxtFileNoLinks ( TString  filename,
bool  incError = true 
) const

Definition at line 1148 of file HyperHistogram.cpp.

1148  {
1149 
1150  if ( _binning->getBinningType() != "HyperBinning" ){
1151  ERROR_LOG << "It is only possible to saveToTxtFile when using HyperBinning. Doing nothing." << std::endl;
1152  return;
1153  }
1154 
1155  const HyperBinning& hyperBinning = dynamic_cast<const HyperBinning&>( getBinning() );
1156 
1157  int nBins = hyperBinning.getNumBins();
1158 
1159  std::ofstream myfile;
1160  myfile.open (filename);
1161 
1162  for (int i = 0; i < nBins; i++ ){
1163  HyperVolume vol = hyperBinning.getBinHyperVolume(i);
1164  HyperCuboid cube = vol.getHyperCuboid(0);
1165  double content = getBinContent(i);
1166  double error = getBinError (i);
1167 
1168  int width = vol.getDimension()*10 + 10;
1169 
1170  myfile << std::setw(width) << std::left << cube.getLowCorner() << " " << std::setw(width) << std::left << cube.getHighCorner();
1171  myfile << " ";
1172  myfile << std::setw(10) << std::left << content;
1173  if (incError) myfile << std::setw(10) << std::left << error;
1174 
1175  myfile << std::endl;
1176  }
1177 
1178  myfile.close();
1179 
1180 
1181 }
double getBinError(int bin) const
#define ERROR_LOG
double getBinContent(int bin) const
TString getBinningType() const
Definition: BinningBase.cpp:50
virtual HyperVolume getBinHyperVolume(int binNumber) const
const int & getDimension() const
Definition: HyperVolume.h:44
const HyperPoint & getHighCorner() const
Definition: HyperCuboid.h:89
BinningBase * _binning
const HyperPoint & getLowCorner() const
Definition: HyperCuboid.h:87
virtual int getNumBins() const
const BinningBase & getBinning() const
const HyperCuboid & getHyperCuboid(int i) const
Definition: HyperVolume.h:54

◆ setContentsFromFunc()

void HyperHistogram::setContentsFromFunc ( const HyperFunction func)

Set the bin contents of the histogram using parsed function. Will set bin errors to zero and use bin centers for evaluating function

Definition at line 384 of file HyperHistogram.cpp.

384  {
385 
386  int nbins = getNBins();
387 
388  for (int i = 0; i < nbins; i++){
390  double funcVal = func.getVal(binCenter);
391  setBinContent(i, funcVal);
392  setBinError (i, 0 );
393  }
394 
395 
396 }
int getNBins() const
Definition: HistogramBase.h:63
void setBinError(int bin, double val)
void setBinContent(int bin, double val)
BinningBase * _binning
virtual double getVal(const HyperPoint &point) const =0
HyperPoint getAverageCenter() const
virtual HyperVolume getBinHyperVolume(int binNumber) const =0

◆ setNames()

void HyperHistogram::setNames ( HyperName  names)

Set the HyperName (mainly used for axis labels)

Definition at line 246 of file HyperHistogram.cpp.

◆ slice() [1/2]

HyperHistogram HyperHistogram::slice ( std::vector< int >  sliceDims,
const HyperPoint slicePoint 
) const

Take a slice of the HyperHistogram and return it as HyperHistogram. The slice is taken in the given slice dimesions i.e. if we had a 5D space with dims [0 1 2 3 4] we could slice through dimensions 2 3 and 4 to return a 2D histogram in 0 vs. 1. The slice in dimensions 2 3 and 4 is taken from the given slicePoint

Definition at line 722 of file HyperHistogram.cpp.

722  {
723 
724  HyperBinningMemRes temp;
725 
726  std::vector<double> binContents;
727  std::vector<double> binErrors ;
728 
729  //Loop over all the bins and slice them. If the slice doesn't
730  //pass through the bin volume, the returned volume with have dimesnion 0.
731  //In these cases, skip the bin as it will not show up in the slice.
732 
733  for (int i = 0; i < getNBins(); i++){
734 
736 
737  HyperVolume slicedVol = vol.slice(slicePoint, sliceDims);
738 
739  if (slicedVol.size() == 0) continue;
740 
741  temp.addHyperVolume(slicedVol);
742 
743  binContents.push_back( getBinContent(i) );
744  binErrors .push_back( getBinError (i) );
745  }
746 
747  //Set the bin contents and errors
748 
749  HyperHistogram slicedHist(temp);
750 
751  for (unsigned i = 0; i < binContents.size(); i++){
752  slicedHist.setBinContent(i, binContents.at(i) );
753  slicedHist.setBinError (i, binErrors .at(i) );
754  }
755 
756  //Set the names and the min/max from this histogram
757 
758  slicedHist.setNames( getNames().slice(sliceDims) );
759 
760  slicedHist.setMin( getMin() );
761  slicedHist.setMax( getMax() );
762 
763  return slicedHist;
764 
765 }
double getBinError(int bin) const
int getNBins() const
Definition: HistogramBase.h:63
int size() const
Definition: HyperVolume.h:60
double getBinContent(int bin) const
double getMin() const
double getMax() const
BinningBase * _binning
HyperVolume slice(const HyperPoint &coords, std::vector< int > dims) const
HyperName getNames() const
virtual bool addHyperVolume(const HyperVolume &hyperVolume, std::vector< int > linkedVolumes=std::vector< int >(0, 0))
virtual HyperVolume getBinHyperVolume(int binNumber) const =0
HyperHistogram slice(std::vector< int > sliceDims, const HyperPoint &slicePoint) const

◆ slice() [2/2]

std::vector< HyperHistogram > HyperHistogram::slice ( std::vector< int >  sliceDims,
const HyperPointSet slicePoints 
) const

Definition at line 768 of file HyperHistogram.cpp.

768  {
769 
770  int nSlices = slicePoints.size();
771 
772  std::vector< HyperBinningMemRes > slicedBinnings(nSlices, HyperBinningMemRes() );
773  std::vector< std::vector<double> > binContents (nSlices, std::vector<double>(0, 0.0) );
774  std::vector< std::vector<double> > binErrors (nSlices, std::vector<double>(0, 0.0) );
775 
776  //Loop over all the bins and slice them. If the slice doesn't
777  //pass through the bin volume, the returned volume with have dimesnion 0.
778  //In these cases, skip the bin as it will not show up in the slice.
779 
780  for (int i = 0; i < getNBins(); i++){
781 
783  double content = getBinContent(i);
784  double error = getBinContent(i);
785 
786  for (int sli = 0; sli < nSlices; sli++){
787  HyperVolume slicedVol = vol.slice( slicePoints.at(sli), sliceDims);
788 
789  if (slicedVol.size() == 0) continue;
790 
791  slicedBinnings.at(sli).addHyperVolume( slicedVol );
792  binContents .at(sli).push_back ( content );
793  binErrors .at(sli).push_back ( error );
794  }
795 
796  }
797 
798  //Set the bin contents and errors
799  //Set the names and the min/max from this histogram
800 
801  std::vector< HyperHistogram > slicedHist;
802  slicedHist.reserve(nSlices);
803 
804  for (int sli = 0; sli < nSlices; sli++){
805 
806  slicedHist.push_back( HyperHistogram( slicedBinnings.at(sli) ) );
807 
808  int nBins = binContents.at(sli).size();
809  for (int bin = 0; bin < nBins; bin++){
810  slicedHist.at(sli).setBinContent(bin, binContents.at(sli).at(bin) );
811  slicedHist.at(sli).setBinError (bin, binErrors .at(sli).at(bin) );
812  }
813 
814  slicedHist.at(sli).setNames( getNames().slice(sliceDims) );
815  slicedHist.at(sli).setMin ( getMin() );
816  slicedHist.at(sli).setMax ( getMax() );
817  }
818 
819  return slicedHist;
820 
821 }
int getNBins() const
Definition: HistogramBase.h:63
const HyperCuboid & at(int i) const
Definition: HyperVolume.h:55
int size() const
Definition: HyperVolume.h:60
double getBinContent(int bin) const
double getMin() const
const HyperPoint & at(int i) const
double getMax() const
BinningBase * _binning
unsigned int size() const
HyperVolume slice(const HyperPoint &coords, std::vector< int > dims) const
HyperName getNames() const
virtual HyperVolume getBinHyperVolume(int binNumber) const =0
HyperHistogram slice(std::vector< int > sliceDims, const HyperPoint &slicePoint) const

Member Data Documentation

◆ _binning

BinningBase* HyperHistogram::_binning
protected

The HyperVolumeBinning used for the HyperHistogram

Definition at line 57 of file HyperHistogram.h.


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