class THnSparse: public TNamed



    Efficient multidimensional histogram.

 Use a THnSparse instead of TH1 / TH2 / TH3 / array for histogramming when
 only a small fraction of bins is filled. A 10-dimensional histogram with 10
 bins per dimension has 10^10 bins; in a naive implementation this will not
 fit in memory. THnSparse only allocates memory for the bins that have
 non-zero bin content instead, drastically reducing both the memory usage
 and the access time.

 To construct a THnSparse object you must use one of its templated, derived
 classes:
 THnSparseD (typedef for THnSparse<ArrayD>): bin content held by a Double_t,
 THnSparseF (typedef for THnSparse<ArrayF>): bin content held by a Float_t,
 THnSparseL (typedef for THnSparse<ArrayL>): bin content held by a Long_t,
 THnSparseI (typedef for THnSparse<ArrayI>): bin content held by an Int_t,
 THnSparseS (typedef for THnSparse<ArrayS>): bin content held by a Short_t,
 THnSparseC (typedef for THnSparse<ArrayC>): bin content held by a Char_t,

 They take name and title, the number of dimensions, and for each dimension
 the number of bins, the minimal, and the maximal value on the dimension's
 axis. A TH2 h("h","h",10, 0., 10., 20, -5., 5.) would correspond to
   Int_t bins[2] = {10, 20};
   Double_t xmin[2] = {0., -5.};
   Double_t xmax[2] = {10., 5.};
   THnSparse hs("hs", "hs", 2, bins, min, max);

 * Filling
 A THnSparse is filled just like a regular histogram, using
 THnSparse::Fill(x, weight), where x is a n-dimensional Double_t value.
 To take errors into account, Sumw2() must be called before filling the
 histogram.
 Bins are allocated as needed; the status of the allocation can be observed
 by GetSparseFractionBins(), GetSparseFractionMem().

 * Fast Bin Content Access
 When iterating over a THnSparse one should only look at filled bins to save
 processing time. The number of filled bins is returned by
 THnSparse::GetNbins(); the bin content for each (linear) bin number can
 be retrieved by THnSparse::GetBinContent(linidx, (Int_t*)coord).
 After the call, coord will contain the bin coordinate of each axis for the bin
 with linear index linidx. A possible call would be
   cout << hs.GetBinContent(0, coord);
   cout <<" is the content of bin [x = " << coord[0] "
        << " | y = " << coord[1] << "]" << endl;

 * Efficiency
 TH1 and TH2 are generally faster than THnSparse for one and two dimensional
 distributions. THnSparse becomes competitive for a sparsely filled TH3
 with large numbers of bins per dimension. The tutorial hist/sparsehist.C
 shows the turning point. On a AMD64 with 8GB memory, THnSparse "wins"
 starting with a TH3 with 30 bins per dimension. Using a THnSparse for a
 one-dimensional histogram is only reasonable if it has a huge number of bins.

 * Projections
 The dimensionality of a THnSparse can be reduced by projecting it to
 1, 2, 3, or n dimensions, which can be represented by a TH1, TH2, TH3, or
 a THnSparse. See the Projection() members. To only project parts of the
 histogram, call
   THnSparse::GetAxis(12)->SetRange(from_bin, to_bin);
 See the important remark in THnSparse::IsInRange() when excluding under-
 and overflow bins!

 * Internal Representation
 An entry for a filled bin consists of its n-dimensional coordinates and
 its bin content. The coordinates are compacted to use as few bits as
 possible; e.g. a histogram with 10 bins in x and 20 bins in y will only
 use 4 bits for the x representation and 5 bits for the y representation.
 This is handled by the internal class THnSparseCompactBinCoord.
 Bin data (content and coordinates) are allocated in chunks of size
 fChunkSize; this parameter can be set when constructing a THnSparse. Each
 chunk is represented by an object of class THnSparseArrayChunk.

 Translation from an n-dimensional bin coordinate to the linear index within
 the chunks is done by GetBin(). It creates a hash from the compacted bin
 coordinates (the hash of a bin coordinate is the compacted coordinate itself
 if it takes less than 4 bytes, the minimal supported size of a Long_t).
 This hash is used to lookup the linear index in the TExMap member fBins;
 the coordinates of the entry fBins points to is compared to the coordinates
 passed to GetBin(). If they do not match, these two coordinates have the same
 hash - which is extremely unlikely but (for the case where the compact bin
 coordinates are larger than 4 bytes) possible. In this case, fBinsContinued
 contains a chain of linear indexes with the same hash. Iterating through this
 chain and comparing each bin coordinates with the one passed to GetBin() will
 retrieve the matching bin.

Function Members (Methods)

 
    This is an abstract class, constructors will not be documented.
    Look at the header to check for available constructors.

public:
virtual~THnSparse()
voidTObject::AbstractMethod(const char* method) const
voidAdd(const THnSparse* h, Double_t c = 1.)
voidAddBinContent(const Int_t* x, Double_t v = 1.)
virtual voidTObject::AppendPad(Option_t* option = "")
virtual voidTObject::Browse(TBrowser* b)
voidCalculateErrors(Bool_t calc = kTRUE)
static TClass*Class()
virtual const char*TObject::ClassName() const
virtual voidTNamed::Clear(Option_t* option = "")
virtual TObject*TNamed::Clone(const char* newname = "") const
virtual Int_tTNamed::Compare(const TObject* obj) const
Double_tComputeIntegral()
virtual voidTNamed::Copy(TObject& named) const
virtual voidTObject::Delete(Option_t* option = "")
virtual Int_tTObject::DistancetoPrimitive(Int_t px, Int_t py)
voidDivide(const THnSparse* h)
voidDivide(const THnSparse* h1, const THnSparse* h2, Double_t c1 = 1., Double_t c2 = 1., Option_t* option = "")
virtual voidTObject::Draw(Option_t* option = "")
virtual voidTObject::DrawClass() const
virtual TObject*TObject::DrawClone(Option_t* option = "") const
virtual voidTObject::Dump() const
virtual voidTObject::Error(const char* method, const char* msgfmt) const
virtual voidTObject::Execute(const char* method, const char* params, Int_t* error = 0)
virtual voidTObject::Execute(TMethod* method, TObjArray* params, Int_t* error = 0)
virtual voidTObject::ExecuteEvent(Int_t event, Int_t px, Int_t py)
virtual voidTObject::Fatal(const char* method, const char* msgfmt) const
Long_tFill(const Double_t* x, Double_t w = 1.)
Long_tFill(const char** name, Double_t w = 1.)
virtual voidTNamed::FillBuffer(char*& buffer)
virtual TObject*TObject::FindObject(const char* name) const
virtual TObject*TObject::FindObject(const TObject* obj) const
TAxis*GetAxis(Int_t dim) const
Long_tGetBin(const Int_t* idx, Bool_t allocate = kTRUE)
Long_tGetBin(const Double_t* x, Bool_t allocate = kTRUE)
Long_tGetBin(const char** name, Bool_t allocate = kTRUE)
Double_tGetBinContent(const Int_t* idx) const
Double_tGetBinContent(Long64_t bin, Int_t* idx = 0) const
Double_tGetBinError(const Int_t* idx) const
Double_tGetBinError(Long64_t linidx) const
Bool_tGetCalculateErrors() const
virtual Option_t*TObject::GetDrawOption() const
static Long_tTObject::GetDtorOnly()
Double_tGetEntries() const
virtual const char*TObject::GetIconName() const
TObjArray*GetListOfAxes()
virtual const char*TNamed::GetName() const
Long64_tGetNbins() const
Int_tGetNChunks() const
Int_tGetNdimensions() const
virtual char*TObject::GetObjectInfo(Int_t px, Int_t py) const
static Bool_tTObject::GetObjectStat()
virtual Option_t*TObject::GetOption() const
voidGetRandom(Double_t* rand, Bool_t subBinRandom = kTRUE)
Double_tGetSparseFractionBins() const
Double_tGetSparseFractionMem() const
Double_tGetSumw() const
Double_tGetSumw2() const
Double_tGetSumwx(Int_t dim) const
Double_tGetSumwx2(Int_t dim) const
virtual const char*TNamed::GetTitle() const
virtual UInt_tTObject::GetUniqueID() const
Double_tGetWeightSum() const
virtual Bool_tTObject::HandleTimer(TTimer* timer)
virtual ULong_tTNamed::Hash() const
virtual voidTObject::Info(const char* method, const char* msgfmt) const
virtual Bool_tTObject::InheritsFrom(const char* classname) const
virtual Bool_tTObject::InheritsFrom(const TClass* cl) const
virtual voidTObject::Inspect() const
voidTObject::InvertBit(UInt_t f)
virtual TClass*IsA() const
virtual Bool_tTObject::IsEqual(const TObject* obj) const
virtual Bool_tTObject::IsFolder() const
Bool_tTObject::IsOnHeap() const
virtual Bool_tTNamed::IsSortable() const
Bool_tTObject::IsZombie() const
virtual voidTNamed::ls(Option_t* option = "") const
voidTObject::MayNotUse(const char* method) const
voidMultiply(const THnSparse* h)
virtual Bool_tTObject::Notify()
static voidTObject::operator delete(void* ptr)
static voidTObject::operator delete(void* ptr, void* vp)
static voidTObject::operator delete[](void* ptr)
static voidTObject::operator delete[](void* ptr, void* vp)
void*TObject::operator new(size_t sz)
void*TObject::operator new(size_t sz, void* vp)
void*TObject::operator new[](size_t sz)
void*TObject::operator new[](size_t sz, void* vp)
TNamed&TNamed::operator=(const TNamed& rhs)
virtual voidTObject::Paint(Option_t* option = "")
virtual voidTObject::Pop()
virtual voidTNamed::Print(Option_t* option = "") const
TH1D*Projection(Int_t xDim, Option_t* option = "") const
TH2D*Projection(Int_t xDim, Int_t yDim, Option_t* option = "") const
THnSparse*Projection(Int_t ndim, const Int_t* dim, Option_t* option = "") const
TH3D*Projection(Int_t xDim, Int_t yDim, Int_t zDim, Option_t* option = "") const
virtual Int_tTObject::Read(const char* name)
THnSparse*Rebin(Int_t group) const
THnSparse*Rebin(const Int_t* group) const
virtual voidTObject::RecursiveRemove(TObject* obj)
voidReset(Option_t* option = "")
voidTObject::ResetBit(UInt_t f)
virtual voidTObject::SaveAs(const char* filename = "", Option_t* option = "") const
virtual voidTObject::SavePrimitive(basic_ostream<char,char_traits<char> >& out, Option_t* option = "")
voidScale(Double_t c)
voidSetBinContent(const Int_t* x, Double_t v)
voidSetBinEdges(Int_t idim, const Double_t* bins)
voidSetBinError(const Int_t* x, Double_t e)
voidTObject::SetBit(UInt_t f)
voidTObject::SetBit(UInt_t f, Bool_t set)
virtual voidTObject::SetDrawOption(Option_t* option = "")
static voidTObject::SetDtorOnly(void* obj)
voidSetEntries(Double_t entries)
virtual voidTNamed::SetName(const char* name)
virtual voidTNamed::SetNameTitle(const char* name, const char* title)
static voidTObject::SetObjectStat(Bool_t stat)
virtual voidTNamed::SetTitle(const char* title = "")
virtual voidTObject::SetUniqueID(UInt_t uid)
virtual voidShowMembers(TMemberInspector& insp, char* parent)
virtual Int_tTNamed::Sizeof() const
virtual voidStreamer(TBuffer& b)
voidStreamerNVirtual(TBuffer& b)
voidSumw2()
virtual voidTObject::SysError(const char* method, const char* msgfmt) const
Bool_tTObject::TestBit(UInt_t f) const
Int_tTObject::TestBits(UInt_t f) const
virtual voidTObject::UseCurrentStyle()
virtual voidTObject::Warning(const char* method, const char* msgfmt) const
virtual Int_tTObject::Write(const char* name = 0, Int_t option = 0, Int_t bufsize = 0)
virtual Int_tTObject::Write(const char* name = 0, Int_t option = 0, Int_t bufsize = 0) const
protected:
THnSparseArrayChunk*AddChunk()
Bool_tCheckConsistency(const THnSparse* h, const char* tag) const
THnSparse*CloneEmpty(const char* name, const char* title, const TObjArray* axes, Int_t chunksize) const
virtual voidTObject::DoError(int level, const char* location, const char* fmt, va_list va) const
Long_tFill(Long_t bin, Double_t w)
virtual TArray*GenerateArray() const
Long_tGetBinIndexForCurrentBin(Bool_t allocate)
THnSparseArrayChunk*GetChunk(Int_t idx) const
Int_tGetChunkSize() const
THnSparseCompactBinCoord*GetCompactCoord() const
Bool_tIsInRange(Int_t* coord) const
voidTObject::MakeZombie()

Data Members

public:
enum { kNoInt
kValidInt
kInvalidInt
};
enum TObject::EStatusBits { kCanDelete
kMustCleanup
kObjInCanvas
kIsReferenced
kHasUUID
kCannotPick
kNoContextMenu
kInvalidObject
};
enum TObject::[unnamed] { kIsOnHeap
kNotDeleted
kZombie
kBitMask
kSingleKey
kOverwrite
kWriteDelete
};
protected:
TStringTNamed::fNameobject identifier
TStringTNamed::fTitleobject title
private:
TObjArrayfAxesaxes of the histogram
TObjArrayfBinContentarray of THnSparseArrayChunk
TExMapfBinsfilled bins
TExMapfBinsContinuedfilled bins for non-unique hashes, containing pairs of (bin index 0, bin index 1)
Int_tfChunkSizenumber of entries for each chunk
THnSparseCompactBinCoord*fCompactCoord! compact coordinate
Double_tfEntriesnumber of entries, spread over chunks
Long64_tfFilledBinsnumber of filled bins
Double_t*fIntegral! array with bin weight sums
enum THnSparse::fIntegralStatus! status of integral
Int_tfNdimensionsnumber of dimensions
Double_tfTsumwtotal sum of weights
Double_tfTsumw2total sum of weights squared; -1 if no errors are calculated
TArrayDfTsumwxtotal sum of weight*X for each dimension
TArrayDfTsumwx2total sum of weight*X*X for each dimension

Class Charts

Inheritance Inherited Members Includes Libraries
Class Charts

Function documentation

~THnSparse()
 Destruct a THnSparse
void AddBinContent(const Int_t* x, Double_t v = 1.)
 Add "v" to the content of bin with coordinates "coord"
THnSparseArrayChunk* AddChunk()
Create a new chunk of bin content
THnSparse* CloneEmpty(const char* name, const char* title, const TObjArray* axes, Int_t chunksize) const
 Create a new THnSparse object that is of the same type as *this,
 but with dimensions and bins given by axes.
Long_t GetBin(const Double_t* x, Bool_t allocate /* = kTRUE */)
 Get the bin index for the n dimensional tuple x,
 allocate one if it doesn't exist yet and "allocate" is true.
Long_t GetBin(const char* name[], Bool_t allocate /* = kTRUE */)
 Get the bin index for the n dimensional tuple addressed by "name",
 allocate one if it doesn't exist yet and "allocate" is true.
Long_t GetBin(const Int_t* coord, Bool_t allocate /*= kTRUE*/)
 Get the bin index for the n dimensional coordinates coord,
 allocate one if it doesn't exist yet and "allocate" is true.
Double_t GetBinContent(const Int_t* idx) const
 Get content of bin with coordinates "coord"
Double_t GetBinContent(Long64_t bin, Int_t* idx = 0) const
 Return the content of the filled bin number "idx".
 If coord is non-null, it will contain the bin's coordinates for each axis
 that correspond to the bin.
Double_t GetBinError(const Int_t *coord)
 Get error of bin with coordinates "coord" as

#sqrt{#sum weight^{2}}
 If errors are not enabled (via Sumw2() or CalculateErrors())
 return sqrt(contents).
Double_t GetBinError(Long64_t linidx)
 Get error of bin addressed by linidx as

#sqrt{#sum weight^{2}}
 If errors are not enabled (via Sumw2() or CalculateErrors())
 return sqrt(contents).
Long_t GetBinIndexForCurrentBin(Bool_t allocate)
 Return the index for fCurrentBinIndex.
 If it doesn't exist then return -1, or allocate a new bin if allocate is set
THnSparseCompactBinCoord* GetCompactCoord()
 Return THnSparseCompactBinCoord object.
void GetRandom(Double_t* rand, Bool_t subBinRandom = kTRUE)
 Generate an n-dimensional random tuple based on the histogrammed
 distribution. If subBinRandom, the returned tuple will be additionally
 randomly distributed within the randomized bin, using a flat
 distribution.
Double_t GetSparseFractionBins()
 Return the amount of filled bins over all bins
Double_t GetSparseFractionMem()
 Return the amount of used memory over memory that would be used by a
 non-sparse n-dimensional histogram. The value is approximate.
Bool_t IsInRange(Int_t* coord) const
 Check whether bin coord is in range, as defined by TAxis::SetRange().
 Currently, TAxis::SetRange() does not allow to select all but over- and
 underflow bins (it instead resets the axis to "no range selected").
 Instead, simply call
    TAxis* axis12 = hsparse.GetAxis(12);
    axis12->SetRange(1, axis12->GetNbins());
    axis12->SetBit(TAxis::kAxisRange);
 to deselect the under- and overflow bins in the 12th dimension.
TH1D* Projection(Int_t xDim, Option_t* option = "") const
 Project all bins into a 1-dimensional histogram,
 keeping only axis "xDim".
 If "option" contains "E" errors will be calculated.
TH2D* Projection(Int_t xDim, Int_t yDim, Option_t* option /*= ""*/)
 Project all bins into a 2-dimensional histogram,
 keeping only axes "xDim" and "yDim".
 If "option" contains "E" errors will be calculated.
TH3D* Projection(Int_t xDim, Int_t yDim, Int_t zDim, Option_t* option = "") const
 Project all bins into a 3-dimensional histogram,
 keeping only axes "xDim", "yDim", and "zDim".
 If "option" contains "E" errors will be calculated.
THnSparse* Projection(Int_t ndim, const Int_t* dim, Option_t* option /*= ""*/)
 Project all bins into a ndim-dimensional histogram,
 keeping only axes "dim".
 If "option" contains "E" errors will be calculated.
void Scale(Double_t c)
 Scale contents and errors of this histogram by c:
 this = this * c
void Add(const THnSparse* h, Double_t c = 1.)
 Add contents of h scaled by c to this histogram:
 this = this + c * h
 Note that if h has Sumw2 set, Sumw2 is automatically called for this
 if not already set.
void Multiply(const THnSparse* h)
 Multiply this histogram by histogram h
 this = this * h
 Note that if h has Sumw2 set, Sumw2 is automatically called for this
 if not already set.
void Divide(const THnSparse* h)
 Divide this histogram by h
 this = this/(h)
 Note that if h has Sumw2 set, Sumw2 is automatically called for
 this if not already set.
 The resulting errors are calculated assuming uncorrelated content.
void Divide(const THnSparse* h1, const THnSparse* h2, Double_t c1 = 1., Double_t c2 = 1., Option_t* option = "")
 Replace contents of this histogram by multiplication of h1 by h2
 this = (c1*h1)/(c2*h2)
 Note that if h1 or h2 have Sumw2 set, Sumw2 is automatically called for
 this if not already set.
 The resulting errors are calculated assuming uncorrelated content.
 However, if option ="B" is specified, Binomial errors are computed.
 In this case c1 and c2 do not make real sense and they are ignored.
Bool_t CheckConsistency(const THnSparse* h, const char* tag) const
 Consistency check on (some of) the parameters of two histograms (for operations).
void SetBinEdges(Int_t idim, const Double_t* bins)
 Set the axis # of bins and bin limits on dimension idim
void SetBinContent(const Int_t* x, Double_t v)
 Set content of bin with coordinates "coord" to "v"
void SetBinError(const Int_t* x, Double_t e)
 Set error of bin with coordinates "coord" to "e", enable errors if needed
void Sumw2()
 Enable calculation of errors
THnSparse* Rebin(Int_t group)
 Combine the content of "group" neighboring bins into
 a new bin and return the resulting THnSparse.
 For group=2 and a 3 dimensional histogram, all "blocks"
 of 2*2*2 bins will be put into a bin.
THnSparse* Rebin(const Int_t* group)
 Combine the content of "group" neighboring bins for each dimension
 into a new bin and return the resulting THnSparse.
 For group={2,1,1} and a 3 dimensional histogram, pairs of x-bins
 will be grouped.
void Reset(Option_t* option = "")
 Clear the histogram
Double_t ComputeIntegral()
 Calculate the integral of the histogram
ULong64_t GetEntries()
{ return fCoordinatesSize / fSingleCoordinateSize; }
Int_t GetChunkSize()
{ return fChunkSize; }
THnSparseArrayChunk* GetChunk(Int_t idx) const
TArray* GenerateArray()
Long_t Fill(Long_t bin, Double_t w)
 Increment the bin content of "bin" by "w",
 return the bin index.
Int_t GetNChunks()
{ return fBinContent.GetEntriesFast(); }
TObjArray* GetListOfAxes()
{ return &fAxes; }
TAxis* GetAxis(Int_t dim) const
{ return (TAxis*)fAxes[dim]; }
Long_t Fill(const Double_t *x, Double_t w = 1.)
return Fill(GetBin(x), w)
Double_t GetWeightSum()
{ return fTsumw; }
Long64_t GetNbins()
{ return fFilledBins; }
Int_t GetNdimensions()
{ return fNdimensions; }
Bool_t GetCalculateErrors()
{ return fTsumw2 >= 0.; }
void CalculateErrors(Bool_t calc = kTRUE)
 Calculate errors (or not if "calc" == kFALSE)
void SetEntries(Double_t entries)
{ fEntries = entries; }
Double_t GetSumw()
{ return fTsumw; }
Double_t GetSumw2()
{ return fTsumw2; }
Double_t GetSumwx(Int_t dim) const
{ return fTsumwx[dim]; }
Double_t GetSumwx2(Int_t dim) const
{ return fTsumwx2[dim]; }

Author: Axel Naumann (2007-09-11)
Last update: root/hist:$Id: THnSparse.h 21401 2007-12-17 10:02:40Z brun $
Copyright (C) 1995-2007, Rene Brun and Fons Rademakers. *

This page has been automatically generated. If you have any comments or suggestions about the page layout send a mail to ROOT support, or contact the developers with any questions or problems regarding ROOT.