THnSparse.cxx

Go to the documentation of this file.
00001 // @(#)root/hist:$Id: THnSparse.cxx 36903 2010-11-24 14:13:17Z moneta $
00002 // Author: Axel Naumann (2007-09-11)
00003 
00004 /*************************************************************************
00005  * Copyright (C) 1995-2007, Rene Brun and Fons Rademakers.               *
00006  * All rights reserved.                                                  *
00007  *                                                                       *
00008  * For the licensing terms see $ROOTSYS/LICENSE.                         *
00009  * For the list of contributors see $ROOTSYS/README/CREDITS.             *
00010  *************************************************************************/
00011 
00012 #include "THnSparse.h"
00013 
00014 #include "TArrayI.h"
00015 #include "TAxis.h"
00016 #include "TClass.h"
00017 #include "TCollection.h"
00018 #include "TDataMember.h"
00019 #include "TDataType.h"
00020 #include "TH1D.h"
00021 #include "TH2D.h"
00022 #include "TH3D.h"
00023 #include "TF1.h"
00024 #include "TInterpreter.h"
00025 #include "TMath.h"
00026 #include "TRandom.h"
00027 
00028 #include "TError.h"
00029 
00030 #include "HFitInterface.h"
00031 #include "Fit/SparseData.h"
00032 #include "Math/MinimizerOptions.h"
00033 #include "Math/WrappedMultiTF1.h"
00034 
00035 //______________________________________________________________________________
00036 //
00037 // THnSparseCoordCompression is a class used by THnSparse internally. It
00038 // represents a compacted n-dimensional array of bin coordinates (indices).
00039 // As the total number of bins in each dimension is known by THnSparse, bin
00040 // indices can be compacted to only use the amount of bins needed by the total
00041 // number of bins in each dimension. E.g. for a THnSparse with
00042 // {15, 100, 2, 20, 10, 100} bins per dimension, a bin index will only occupy
00043 // 28 bits (4+7+1+5+4+7), i.e. less than a 32bit integer. The tricky part is
00044 // the fast compression and decompression, the platform-independent storage
00045 // (think of endianness: the bits of the number 0x123456 depend on the
00046 // platform), and the hashing needed by THnSparseArrayChunk.
00047 //______________________________________________________________________________
00048 
00049 class THnSparseCoordCompression {
00050 public:
00051    THnSparseCoordCompression(Int_t dim, const Int_t* nbins);
00052    THnSparseCoordCompression(const THnSparseCoordCompression& other);
00053    ~THnSparseCoordCompression();
00054 
00055    ULong64_t GetHashFromBuffer(const Char_t* buf) const;
00056    Int_t     GetBufferSize() const { return fCoordBufferSize; }
00057    Int_t     GetNdimensions() const { return fNdimensions; }
00058    void      SetCoordFromBuffer(const Char_t* buf_in, Int_t* coord_out) const;
00059    ULong64_t SetBufferFromCoord(const Int_t* coord_in, Char_t* buf_out) const;
00060 
00061 protected:
00062    Int_t     GetNumBits(Int_t n) const {
00063       // return the number of bits allocated by the number "n"
00064       Int_t r = (n > 0);
00065       while (n/=2) ++r;
00066       return r;
00067    }
00068 private:
00069    Int_t  fNdimensions;     // number of dimensions
00070    Int_t  fCoordBufferSize; // size of coordbuf
00071    Int_t *fBitOffsets;      //[fNdimensions + 1] bit offset of each axis index
00072 };
00073 
00074 
00075 //______________________________________________________________________________
00076 //______________________________________________________________________________
00077 
00078 
00079 //______________________________________________________________________________
00080 THnSparseCoordCompression::THnSparseCoordCompression(Int_t dim, const Int_t* nbins):
00081    fNdimensions(dim), fCoordBufferSize(0), fBitOffsets(0)
00082 {
00083    // Initialize a THnSparseCoordCompression object with "dim" dimensions
00084    // and "bins" holding the number of bins for each dimension; it
00085    // stores the 
00086 
00087    fBitOffsets = new Int_t[dim + 1];
00088 
00089    int shift = 0;
00090    for (Int_t i = 0; i < dim; ++i) {
00091       fBitOffsets[i] = shift;
00092       shift += GetNumBits(nbins[i] + 2);
00093    }
00094    fBitOffsets[dim] = shift;
00095    fCoordBufferSize = (shift + 7) / 8;
00096 }
00097 
00098 
00099 //______________________________________________________________________________
00100 THnSparseCoordCompression::THnSparseCoordCompression(const THnSparseCoordCompression& other)
00101 {
00102    // Construct a THnSparseCoordCompression from another one
00103    fNdimensions = other.fNdimensions;
00104    fCoordBufferSize = other.fCoordBufferSize;
00105    fBitOffsets = new Int_t[fNdimensions + 1];
00106    memcpy(fBitOffsets, other.fBitOffsets, sizeof(Int_t) * fNdimensions);
00107 }
00108 
00109 
00110 //______________________________________________________________________________
00111 THnSparseCoordCompression::~THnSparseCoordCompression()
00112 {
00113    // destruct a THnSparseCoordCompression
00114    delete [] fBitOffsets;
00115 }
00116 
00117 
00118 //______________________________________________________________________________
00119 void THnSparseCoordCompression::SetCoordFromBuffer(const Char_t* buf_in,
00120                                                   Int_t* coord_out) const
00121 {
00122    // Given the compressed coordinate buffer buf_in, calculate ("decompact")
00123    // the bin coordinates and return them in coord_out.
00124 
00125    for (Int_t i = 0; i < fNdimensions; ++i) {
00126       const Int_t offset = fBitOffsets[i] / 8;
00127       Int_t shift = fBitOffsets[i] % 8;
00128       Int_t nbits = fBitOffsets[i + 1] - fBitOffsets[i];
00129       const UChar_t* pbuf = (const UChar_t*) buf_in + offset;
00130       coord_out[i] = *pbuf >> shift;
00131       Int_t subst = (Int_t) -1;
00132       subst = subst << nbits;
00133       nbits -= (8 - shift);
00134       shift = 8 - shift;
00135       for (Int_t n = 0; n * 8 < nbits; ++n) {
00136          ++pbuf;
00137          coord_out[i] += *pbuf << shift;
00138          shift += 8;
00139       }
00140       coord_out[i] &= ~subst;
00141    }
00142 }
00143 
00144 
00145 //______________________________________________________________________________
00146 ULong64_t THnSparseCoordCompression::SetBufferFromCoord(const Int_t* coord_in,
00147                                                        Char_t* buf_out) const
00148 {
00149    // Given the cbin coordinates coord_in, calculate ("compact")
00150    // the bin coordinates and return them in buf_in.
00151    // Return the hash value.
00152 
00153    if (fCoordBufferSize <= 8) {
00154       ULong64_t l64buf = 0;
00155       for (Int_t i = 0; i < fNdimensions; ++i) {
00156          l64buf += ((ULong64_t)((UInt_t)coord_in[i])) << fBitOffsets[i];
00157       }
00158       memcpy(buf_out, &l64buf, sizeof(Long64_t));
00159       return l64buf;
00160    }
00161 
00162    // else: doesn't fit into a Long64_t:
00163    memset(buf_out, 0, fCoordBufferSize);
00164    for (Int_t i = 0; i < fNdimensions; ++i) {
00165       const Int_t offset = fBitOffsets[i] / 8;
00166       const Int_t shift = fBitOffsets[i] % 8;
00167       ULong64_t val = coord_in[i];
00168 
00169       Char_t* pbuf = buf_out + offset;
00170       *pbuf += 0xff & (val << shift);
00171       val = val >> (8 - shift);
00172       while (val) {
00173          ++pbuf;
00174          *pbuf += 0xff & val;
00175          val = val >> 8;
00176       }
00177    }
00178 
00179    return GetHashFromBuffer(buf_out);
00180 }
00181 
00182 /*
00183 //______________________________________________________________________________
00184 ULong64_t THnSparseCoordCompression::GetHashFromCoords(const Int_t* coord) const
00185 {
00186    // Calculate hash from bin indexes.
00187 
00188    // Bins are addressed in two different modes, depending
00189    // on whether the compact bin index fits into a Long64_t or not.
00190    // If it does, we can use it as a "perfect hash" for the TExMap.
00191    // If not we build a hash from the compact bin index, and use that
00192    // as the TExMap's hash.
00193 
00194    if (fCoordBufferSize <= 8) {
00195       // fits into a Long64_t
00196       ULong64_t hash1 = 0;
00197       for (Int_t i = 0; i < fNdimensions; ++i) {
00198          hash1 += coord[i] << fBitOffsets[i];
00199       }
00200       return hash1;
00201    }
00202 
00203    // else: doesn't fit into a Long64_t:
00204    memset(coord, 0, fCoordBufferSize);
00205    for (Int_t i = 0; i < fNdimensions; ++i) {
00206       const Int_t offset = fBitOffsets[i] / 8;
00207       const Int_t shift = fBitOffsets[i] % 8;
00208       ULong64_t val = coord[i];
00209 
00210       Char_t* pbuf = fCoordBuffer + offset;
00211       *pbuf += 0xff & (val << shift);
00212       val = val >> (8 - shift);
00213       while (val) {
00214          ++pbuf;
00215          *pbuf += 0xff & val;
00216          val = val >> 8;
00217       }
00218    }
00219 
00220    ULong64_t hash = 5381;
00221    Char_t* str = fCoordBuffer;
00222    while (str - fCoordBuffer < fCoordBufferSize) {
00223       hash *= 5;
00224       hash += *(str++);
00225    }
00226    return hash;
00227 }
00228 */
00229 
00230 
00231 //______________________________________________________________________________
00232 ULong64_t THnSparseCoordCompression::GetHashFromBuffer(const Char_t* buf) const
00233 {
00234    // Calculate hash from compact bin index.
00235 
00236    // Bins are addressed in two different modes, depending
00237    // on whether the compact bin index fits into a Long64_t or not.
00238    // If it does, we can use it as a "perfect hash" for the TExMap.
00239    // If not we build a hash from the compact bin index, and use that
00240    // as the TExMap's hash.
00241 
00242    if (fCoordBufferSize <= 8) {
00243       // fits into a Long64_t
00244       ULong64_t hash1 = 0;
00245       memcpy(&hash1, buf, fCoordBufferSize);
00246       return hash1;
00247    }
00248 
00249    // else: doesn't fit into a Long64_t:
00250    ULong64_t hash = 5381;
00251    const Char_t* str = buf;
00252    while (str - buf < fCoordBufferSize) {
00253       hash *= 5;
00254       hash += *(str++);
00255    }
00256    return hash;
00257 }
00258 
00259 
00260 
00261 
00262 //______________________________________________________________________________
00263 //
00264 // THnSparseCompactBinCoord is a class used by THnSparse internally. It
00265 // maps between an n-dimensional array of bin coordinates (indices) and
00266 // its compact version, the THnSparseCoordCompression.
00267 //______________________________________________________________________________
00268 
00269 class THnSparseCompactBinCoord: public THnSparseCoordCompression {
00270 public:
00271    THnSparseCompactBinCoord(Int_t dim, const Int_t* nbins);
00272    ~THnSparseCompactBinCoord();
00273    Int_t*    GetCoord() { return fCurrentBin; }
00274    const Char_t*   GetBuffer() const { return fCoordBuffer; }
00275    ULong64_t GetHash() const { return fHash; }
00276    void UpdateCoord() {
00277       fHash = SetBufferFromCoord(fCurrentBin, fCoordBuffer);
00278    }
00279    void      SetCoord(const Int_t* coord) {
00280       memcpy(fCurrentBin, coord, sizeof(Int_t) * GetNdimensions());
00281       fHash = SetBufferFromCoord(coord, fCoordBuffer);
00282    }
00283    void      SetBuffer(const Char_t* buf) {
00284       memcpy(fCoordBuffer, buf, GetBufferSize());
00285       fHash = GetHashFromBuffer(fCoordBuffer);
00286    }
00287 
00288 private:
00289    ULong64_t fHash;      // hash for current coordinates; 0 if not calculated
00290    Char_t *fCoordBuffer; // compact buffer of coordinates
00291    Int_t *fCurrentBin;   // current coordinates
00292 };
00293 
00294 
00295 //______________________________________________________________________________
00296 //______________________________________________________________________________
00297 
00298 
00299 //______________________________________________________________________________
00300 THnSparseCompactBinCoord::THnSparseCompactBinCoord(Int_t dim, const Int_t* nbins):
00301    THnSparseCoordCompression(dim, nbins),
00302    fHash(0), fCoordBuffer(0), fCurrentBin(0)
00303 {
00304    // Initialize a THnSparseCompactBinCoord object with "dim" dimensions
00305    // and "bins" holding the number of bins for each dimension.
00306 
00307    fCurrentBin = new Int_t[dim];
00308    size_t bufAllocSize = GetBufferSize();
00309    if (bufAllocSize < sizeof(Long64_t))
00310       bufAllocSize = sizeof(Long64_t);
00311    fCoordBuffer = new Char_t[bufAllocSize];
00312 }
00313 
00314 
00315 //______________________________________________________________________________
00316 THnSparseCompactBinCoord::~THnSparseCompactBinCoord()
00317 {
00318    // destruct a THnSparseCompactBinCoord
00319 
00320    delete [] fCoordBuffer;
00321    delete [] fCurrentBin;
00322 }
00323 
00324 //______________________________________________________________________________
00325 //
00326 // THnSparseArrayChunk is used internally by THnSparse.
00327 //
00328 // THnSparse stores its (dynamic size) array of bin coordinates and their
00329 // contents (and possibly errors) in a TObjArray of THnSparseArrayChunk. Each
00330 // of the chunks holds an array of THnSparseCompactBinCoord and the content
00331 // (a TArray*), which is created outside (by the templated derived classes of
00332 // THnSparse) and passed in at construction time.
00333 //______________________________________________________________________________
00334 
00335 
00336 ClassImp(THnSparseArrayChunk);
00337 
00338 //______________________________________________________________________________
00339 THnSparseArrayChunk::THnSparseArrayChunk(Int_t coordsize, bool errors, TArray* cont):
00340       fCoordinateAllocationSize(-1), fSingleCoordinateSize(coordsize), fCoordinatesSize(0),
00341       fCoordinates(0), fContent(cont),
00342       fSumw2(0)
00343 {
00344    // (Default) initialize a chunk. Takes ownership of cont (~THnSparseArrayChunk deletes it),
00345    // and create an ArrayF for errors if "errors" is true.
00346 
00347    fCoordinateAllocationSize = fSingleCoordinateSize * cont->GetSize();
00348    fCoordinates = new Char_t[fCoordinateAllocationSize];
00349    if (errors) Sumw2();
00350 }
00351 
00352 //______________________________________________________________________________
00353 THnSparseArrayChunk::~THnSparseArrayChunk()
00354 {
00355    // Destructor
00356    delete fContent;
00357    delete [] fCoordinates;
00358    delete fSumw2;
00359 }
00360 
00361 //______________________________________________________________________________
00362 void THnSparseArrayChunk::AddBin(Int_t idx, const Char_t* coordbuf)
00363 {
00364    // Create a new bin in this chunk
00365 
00366    // When streaming out only the filled chunk is saved.
00367    // When reading back only the memory needed for that filled part gets
00368    // allocated. We need to check whether the allowed chunk size is
00369    // bigger than the allocated size. If fCoordinateAllocationSize is
00370    // set to -1 this chunk has been allocated by the  streamer and the
00371    // buffer allocation size is defined by [fCoordinatesSize]. In that
00372    // case we need to compare fCoordinatesSize to
00373    // fSingleCoordinateSize * fContent->GetSize()
00374    // to determine whether we need to expand the buffer.
00375    if (fCoordinateAllocationSize == -1 && fContent) {
00376       Int_t chunksize = fSingleCoordinateSize * fContent->GetSize();
00377       if (fCoordinatesSize < chunksize) {
00378          // need to re-allocate:
00379          Char_t *newcoord = new Char_t[chunksize];
00380          memcpy(newcoord, fCoordinates, fCoordinatesSize);
00381          delete [] fCoordinates;
00382          fCoordinates = newcoord;
00383       }
00384       fCoordinateAllocationSize = chunksize;
00385    }
00386 
00387    memcpy(fCoordinates + idx * fSingleCoordinateSize, coordbuf, fSingleCoordinateSize);
00388    fCoordinatesSize += fSingleCoordinateSize;
00389 }
00390 
00391 //______________________________________________________________________________
00392 void THnSparseArrayChunk::Sumw2()
00393 {
00394    // Turn on support of errors
00395    if (!fSumw2)
00396       fSumw2 = new TArrayD(fContent->GetSize());
00397 }
00398 
00399 
00400 
00401 //______________________________________________________________________________
00402 //
00403 //
00404 //    Efficient multidimensional histogram.
00405 //
00406 // Use a THnSparse instead of TH1 / TH2 / TH3 / array for histogramming when
00407 // only a small fraction of bins is filled. A 10-dimensional histogram with 10
00408 // bins per dimension has 10^10 bins; in a naive implementation this will not
00409 // fit in memory. THnSparse only allocates memory for the bins that have
00410 // non-zero bin content instead, drastically reducing both the memory usage
00411 // and the access time.
00412 //
00413 // To construct a THnSparse object you must use one of its templated, derived
00414 // classes:
00415 // THnSparseD (typedef for THnSparse<ArrayD>): bin content held by a Double_t,
00416 // THnSparseF (typedef for THnSparse<ArrayF>): bin content held by a Float_t,
00417 // THnSparseL (typedef for THnSparse<ArrayL>): bin content held by a Long_t,
00418 // THnSparseI (typedef for THnSparse<ArrayI>): bin content held by an Int_t,
00419 // THnSparseS (typedef for THnSparse<ArrayS>): bin content held by a Short_t,
00420 // THnSparseC (typedef for THnSparse<ArrayC>): bin content held by a Char_t,
00421 //
00422 // They take name and title, the number of dimensions, and for each dimension
00423 // the number of bins, the minimal, and the maximal value on the dimension's
00424 // axis. A TH2 h("h","h",10, 0., 10., 20, -5., 5.) would correspond to
00425 //   Int_t bins[2] = {10, 20};
00426 //   Double_t xmin[2] = {0., -5.};
00427 //   Double_t xmax[2] = {10., 5.};
00428 //   THnSparse hs("hs", "hs", 2, bins, min, max);
00429 //
00430 // * Filling
00431 // A THnSparse is filled just like a regular histogram, using
00432 // THnSparse::Fill(x, weight), where x is a n-dimensional Double_t value.
00433 // To take errors into account, Sumw2() must be called before filling the
00434 // histogram.
00435 // Bins are allocated as needed; the status of the allocation can be observed
00436 // by GetSparseFractionBins(), GetSparseFractionMem().
00437 //
00438 // * Fast Bin Content Access
00439 // When iterating over a THnSparse one should only look at filled bins to save
00440 // processing time. The number of filled bins is returned by
00441 // THnSparse::GetNbins(); the bin content for each (linear) bin number can
00442 // be retrieved by THnSparse::GetBinContent(linidx, (Int_t*)coord).
00443 // After the call, coord will contain the bin coordinate of each axis for the bin
00444 // with linear index linidx. A possible call would be
00445 //   cout << hs.GetBinContent(0, coord);
00446 //   cout <<" is the content of bin [x = " << coord[0] "
00447 //        << " | y = " << coord[1] << "]" << endl;
00448 //
00449 // * Efficiency
00450 // TH1 and TH2 are generally faster than THnSparse for one and two dimensional
00451 // distributions. THnSparse becomes competitive for a sparsely filled TH3
00452 // with large numbers of bins per dimension. The tutorial hist/sparsehist.C
00453 // shows the turning point. On a AMD64 with 8GB memory, THnSparse "wins"
00454 // starting with a TH3 with 30 bins per dimension. Using a THnSparse for a
00455 // one-dimensional histogram is only reasonable if it has a huge number of bins.
00456 //
00457 // * Projections
00458 // The dimensionality of a THnSparse can be reduced by projecting it to
00459 // 1, 2, 3, or n dimensions, which can be represented by a TH1, TH2, TH3, or
00460 // a THnSparse. See the Projection() members. To only project parts of the
00461 // histogram, call
00462 //   THnSparse::GetAxis(12)->SetRange(from_bin, to_bin);
00463 // See the important remark in THnSparse::IsInRange() when excluding under-
00464 // and overflow bins!
00465 //
00466 // * Internal Representation
00467 // An entry for a filled bin consists of its n-dimensional coordinates and
00468 // its bin content. The coordinates are compacted to use as few bits as
00469 // possible; e.g. a histogram with 10 bins in x and 20 bins in y will only
00470 // use 4 bits for the x representation and 5 bits for the y representation.
00471 // This is handled by the internal class THnSparseCompactBinCoord.
00472 // Bin data (content and coordinates) are allocated in chunks of size
00473 // fChunkSize; this parameter can be set when constructing a THnSparse. Each
00474 // chunk is represented by an object of class THnSparseArrayChunk.
00475 //
00476 // Translation from an n-dimensional bin coordinate to the linear index within
00477 // the chunks is done by GetBin(). It creates a hash from the compacted bin
00478 // coordinates (the hash of a bin coordinate is the compacted coordinate itself
00479 // if it takes less than 8 bytes, the size of a Long64_t.
00480 // This hash is used to lookup the linear index in the TExMap member fBins;
00481 // the coordinates of the entry fBins points to is compared to the coordinates
00482 // passed to GetBin(). If they do not match, these two coordinates have the same
00483 // hash - which is extremely unlikely but (for the case where the compact bin
00484 // coordinates are larger than 4 bytes) possible. In this case, fBinsContinued
00485 // contains a chain of linear indexes with the same hash. Iterating through this
00486 // chain and comparing each bin coordinates with the one passed to GetBin() will
00487 // retrieve the matching bin.
00488 
00489 
00490 ClassImp(THnSparse);
00491 
00492 //______________________________________________________________________________
00493 THnSparse::THnSparse():
00494    fNdimensions(0), fChunkSize(1024), fFilledBins(0), fEntries(0),
00495    fTsumw(0), fTsumw2(-1.), fCompactCoord(0), fIntegral(0), fIntegralStatus(kNoInt)
00496 {
00497    // Construct an empty THnSparse.
00498    fBinContent.SetOwner();
00499 }
00500 
00501 //______________________________________________________________________________
00502 THnSparse::THnSparse(const char* name, const char* title, Int_t dim,
00503                      const Int_t* nbins, const Double_t* xmin, const Double_t* xmax,
00504                      Int_t chunksize):
00505    TNamed(name, title), fNdimensions(dim), fChunkSize(chunksize), fFilledBins(0),
00506    fAxes(dim), fEntries(0), fTsumw(0), fTsumw2(-1.), fTsumwx(dim), fTsumwx2(dim),
00507    fCompactCoord(0), fIntegral(0), fIntegralStatus(kNoInt)
00508 {
00509    // Construct a THnSparse with "dim" dimensions,
00510    // with chunksize as the size of the chunks.
00511    // "nbins" holds the number of bins for each dimension;
00512    // "xmin" and "xmax" the minimal and maximal value for each dimension.
00513    // The arrays "xmin" and "xmax" can be NULL; in that case SetBinEdges()
00514    // must be called for each dimension.
00515 
00516    for (Int_t i = 0; i < fNdimensions; ++i) {
00517       TAxis* axis = new TAxis(nbins[i], xmin ? xmin[i] : 0., xmax ? xmax[i] : 1.);
00518       fAxes.AddAtAndExpand(axis, i);
00519    }
00520    SetTitle(title);
00521    fAxes.SetOwner();
00522 
00523    fCompactCoord = new THnSparseCompactBinCoord(dim, nbins);
00524    fBinContent.SetOwner();
00525 }
00526 
00527 //______________________________________________________________________________
00528 THnSparse::~THnSparse() {
00529    // Destruct a THnSparse
00530 
00531    delete fCompactCoord;
00532    if (fIntegralStatus != kNoInt) delete [] fIntegral;
00533 }
00534 
00535 //______________________________________________________________________________
00536 void THnSparse::AddBinContent(const Int_t* coord, Double_t v)
00537 {
00538    // Add "v" to the content of bin with coordinates "coord"
00539 
00540    GetCompactCoord()->SetCoord(coord);
00541    Long64_t bin = GetBinIndexForCurrentBin(kTRUE);
00542    return AddBinContent(bin, v);
00543 }
00544 
00545 //______________________________________________________________________________
00546 void THnSparse::AddBinContent(Long64_t bin, Double_t v)
00547 {
00548    // Add "v" to the content of bin with index "bin"
00549 
00550    THnSparseArrayChunk* chunk = GetChunk(bin / fChunkSize);
00551    bin %= fChunkSize;
00552    v += chunk->fContent->GetAt(bin);
00553    return chunk->fContent->SetAt(v, bin);
00554 }
00555 
00556 //______________________________________________________________________________
00557 THnSparseArrayChunk* THnSparse::AddChunk()
00558 {
00559    //Create a new chunk of bin content
00560    THnSparseArrayChunk* chunk =
00561       new THnSparseArrayChunk(GetCompactCoord()->GetBufferSize(),
00562                               GetCalculateErrors(), GenerateArray());
00563    fBinContent.AddLast(chunk);
00564    return chunk;
00565 }
00566 
00567 //______________________________________________________________________________
00568 THnSparse* THnSparse::CloneEmpty(const char* name, const char* title,
00569                                  const TObjArray* axes, Int_t chunksize,
00570                                  Bool_t keepTargetAxis) const
00571 {
00572    // Create a new THnSparse object that is of the same type as *this,
00573    // but with dimensions and bins given by axes.
00574    // If keepTargetAxis is true, the axes will keep their original xmin / xmax,
00575    // else they will be restricted to the range selected (first / last).
00576 
00577    THnSparse* ret = (THnSparse*)IsA()->New();
00578    ret->SetNameTitle(name, title);
00579    ret->fNdimensions = axes->GetEntriesFast();
00580    ret->fChunkSize = chunksize;
00581 
00582    TIter iAxis(axes);
00583    const TAxis* axis = 0;
00584    Int_t pos = 0;
00585    Int_t *nbins = new Int_t[axes->GetEntriesFast()];
00586    while ((axis = (TAxis*)iAxis())) {
00587       TAxis* reqaxis = (TAxis*)axis->Clone();
00588       if (!keepTargetAxis && axis->TestBit(TAxis::kAxisRange)) {
00589          Int_t binFirst = axis->GetFirst();
00590          Int_t binLast = axis->GetLast();
00591          Int_t nBins = binLast - binFirst + 1;
00592          if (axis->GetXbins()->GetSize()) {
00593             // non-uniform bins:
00594             reqaxis->Set(nBins, axis->GetXbins()->GetArray() + binFirst - 1);
00595          } else {
00596             // uniform bins:
00597             reqaxis->Set(nBins, axis->GetBinLowEdge(binFirst), axis->GetBinUpEdge(binLast));
00598          }
00599          reqaxis->ResetBit(TAxis::kAxisRange);
00600       }
00601 
00602       nbins[pos] = reqaxis->GetNbins();
00603       ret->fAxes.AddAtAndExpand(reqaxis->Clone(), pos++);
00604    }
00605    ret->fAxes.SetOwner();
00606 
00607    ret->fCompactCoord = new THnSparseCompactBinCoord(pos, nbins);
00608    delete [] nbins;
00609 
00610    return ret;
00611 }
00612 
00613 //______________________________________________________________________________
00614 TH1* THnSparse::CreateHist(const char* name, const char* title,
00615                            const TObjArray* axes,
00616                            Bool_t keepTargetAxis ) const {
00617    // Create an empty histogram with name and title with a given
00618    // set of axes. Create a TH1D/TH2D/TH3D, depending on the number
00619    // of elements in axes.
00620 
00621    const int ndim = axes->GetSize();
00622 
00623    TH1* hist = 0;
00624    // create hist with dummy axes, we fix them later.
00625    if (ndim == 1)
00626       hist = new TH1D(name, title, 1, 0., 1.);
00627    else if (ndim == 2)
00628       hist = new TH2D(name, title, 1, 0., 1., 1, 0., 1.);
00629    else if (ndim == 3)
00630       hist = new TH3D(name, title, 1, 0., 1., 1, 0., 1., 1, 0., 1.);
00631    else {
00632       Error("CreateHist", "Cannot create histogram %s with %d dimensions!", name, ndim);
00633       return 0;
00634    }
00635 
00636    TAxis* hax[3] = {hist->GetXaxis(), hist->GetYaxis(), hist->GetZaxis()};
00637    for (Int_t d = 0; d < ndim; ++d) {
00638       TAxis* reqaxis = (TAxis*)(*axes)[d];
00639       hax[d]->SetTitle(reqaxis->GetTitle());
00640       if (!keepTargetAxis && reqaxis->TestBit(TAxis::kAxisRange)) {
00641          Int_t binFirst = reqaxis->GetFirst();
00642          if (binFirst == 0) binFirst = 1;
00643          Int_t binLast = reqaxis->GetLast();
00644          Int_t nBins = binLast - binFirst + 1;
00645          if (reqaxis->GetXbins()->GetSize()) {
00646             // non-uniform bins:
00647             hax[d]->Set(nBins, reqaxis->GetXbins()->GetArray() + binFirst - 1);
00648          } else {
00649             // uniform bins:
00650             hax[d]->Set(nBins, reqaxis->GetBinLowEdge(binFirst), reqaxis->GetBinUpEdge(binLast));
00651          }
00652       } else {
00653          if (reqaxis->GetXbins()->GetSize()) {
00654             // non-uniform bins:
00655             hax[d]->Set(reqaxis->GetNbins(), reqaxis->GetXbins()->GetArray());
00656          } else {
00657             // uniform bins:
00658             hax[d]->Set(reqaxis->GetNbins(), reqaxis->GetXmin(), reqaxis->GetXmax());
00659          }
00660       }
00661    }
00662 
00663    hist->Rebuild();
00664 
00665    return hist;
00666 }
00667 
00668 //______________________________________________________________________________
00669 THnSparse* THnSparse::CreateSparse(const char* name, const char* title,
00670                                    const TH1* h, Int_t chunkSize)
00671 {
00672    // Create a THnSparse object from a histogram deriving from TH1.
00673 
00674    // Get the dimension of the TH1
00675    int ndim = h->GetDimension();
00676 
00677    // Axis properties
00678    int nbins[3] = {0,0,0};
00679    double minRange[3] = {0.,0.,0.};
00680    double maxRange[3] = {0.,0.,0.};
00681    TAxis* axis[3] = { h->GetXaxis(), h->GetYaxis(), h->GetZaxis() };
00682    for (int i = 0; i < ndim; ++i) {
00683       nbins[i]    = axis[i]->GetNbins();
00684       minRange[i] = axis[i]->GetXmin();
00685       maxRange[i] = axis[i]->GetXmax();
00686    }
00687 
00688    // Create the corresponding THnSparse, depending on the storage
00689    // type of the TH1. The class name will be "TH??\0" where the first
00690    // ? is 1,2 or 3 and the second ? indicates the stograge as C, S,
00691    // I, F or D.
00692    THnSparse* s = 0;
00693    const char* cname( h->ClassName() );
00694    if (cname[0] == 'T' && cname[1] == 'H' && cname[2] >= '1' && cname[2] <= '3' && cname[4] == 0) {
00695       if (cname[3] == 'F') {
00696          s = new THnSparseF(name, title, ndim, nbins, minRange, maxRange, chunkSize);
00697       } else if (cname[3] == 'D') {
00698          s = new THnSparseD(name, title, ndim, nbins, minRange, maxRange, chunkSize);
00699       } else if (cname[3] == 'I') {
00700          s = new THnSparseI(name, title, ndim, nbins, minRange, maxRange, chunkSize);
00701       } else if (cname[3] == 'S') {
00702          s = new THnSparseS(name, title, ndim, nbins, minRange, maxRange, chunkSize);
00703       } else if (cname[3] == 'C') {
00704          s = new THnSparseC(name, title, ndim, nbins, minRange, maxRange, chunkSize);
00705       }
00706    }
00707    if (!s) {
00708       ::Warning("THnSparse::CreateSparse", "Unknown Type of Histogram");
00709       return 0;
00710    }
00711 
00712    for (int i = 0; i < ndim; ++i) {
00713       s->GetAxis(i)->SetTitle(axis[i]->GetTitle());
00714    }
00715 
00716    // Get the array to know the number of entries of the TH1
00717    const TArray *array = dynamic_cast<const TArray*>(h);
00718    if (!array) {
00719       ::Warning("THnSparse::CreateSparse", "Unknown Type of Histogram");
00720       return 0;
00721    }
00722 
00723    // Fill the THnSparse with the bins that have content.
00724    for (int i = 0; i < array->GetSize(); ++i) {
00725       double value = h->GetBinContent(i);
00726       double error = h->GetBinError(i);
00727       if (!value && !error) continue;
00728       int x[3] = {0,0,0};
00729       h->GetBinXYZ(i, x[0], x[1], x[2]);
00730       s->SetBinContent(x, value);
00731       s->SetBinError(x, error);
00732    }
00733 
00734    return s;
00735 }
00736 
00737 //______________________________________________________________________________
00738 void THnSparse::FillExMap()
00739 {
00740    //We have been streamed; set up fBins
00741    TIter iChunk(&fBinContent);
00742    THnSparseArrayChunk* chunk = 0;
00743    THnSparseCoordCompression compactCoord(*GetCompactCoord());
00744    Long64_t idx = 0;
00745    if (2 * GetNbins() > fBins.Capacity())
00746       fBins.Expand(3 * GetNbins());
00747    while ((chunk = (THnSparseArrayChunk*) iChunk())) {
00748       const Int_t chunkSize = chunk->GetEntries();
00749       Char_t* buf = chunk->fCoordinates;
00750       const Int_t singleCoordSize = chunk->fSingleCoordinateSize;
00751       const Char_t* endbuf = buf + singleCoordSize * chunkSize;
00752       for (; buf < endbuf; buf += singleCoordSize, ++idx) {
00753          Long64_t hash = compactCoord.GetHashFromBuffer(buf);
00754          Long64_t linidx = fBins.GetValue(hash);
00755          if (linidx) {
00756             Long64_t nextidx = linidx;
00757             while (nextidx) {
00758                // must be a collision, so go to fBinsContinued.
00759                linidx = nextidx;
00760                nextidx = fBinsContinued.GetValue(linidx);
00761             }
00762             fBinsContinued.Add(linidx, idx + 1);
00763          } else {
00764             fBins.Add(hash, idx + 1);
00765          }
00766       }
00767    }
00768 }
00769 
00770 //______________________________________________________________________________
00771 TFitResultPtr THnSparse::Fit(TF1 *f ,Option_t *option ,Option_t *goption)
00772 {
00773 //   Fit a THnSparse with function f
00774 // 
00775 //   since the data is sparse by default a likelihood fit is performed 
00776 //   merging all the regions with empty bins for betetr performance efficiency
00777 // 
00778 //  Since the THnSparse is not drawn no graphics options are passed 
00779 //  Here is the list of possible options 
00780 // 
00781 //                = "I"  Use integral of function in bin instead of value at bin center
00782 //                = "X"  Use chi2 method (default is log-likelihood method)
00783 //                = "U"  Use a User specified fitting algorithm (via SetFCN)
00784 //                = "Q"  Quiet mode (minimum printing)
00785 //                = "V"  Verbose mode (default is between Q and V)
00786 //                = "E"  Perform better Errors estimation using Minos technique
00787 //                = "B"  Use this option when you want to fix one or more parameters
00788 //                       and the fitting function is like "gaus", "expo", "poln", "landau".
00789 //                = "M"  More. Improve fit results
00790 //                = "R"  Use the Range specified in the function range
00791 
00792 
00793    Foption_t fitOption;
00794 
00795    if (!TH1::FitOptionsMake(option,fitOption)) return 0;
00796 
00797    // The function used to fit cannot be stored in a THnSparse. It
00798    // cannot be drawn either. Perhaps in the future.
00799    fitOption.Nostore = true;
00800    // Use likelihood fit if not specified
00801    if (!fitOption.Chi2) fitOption.Like = true;
00802    // create range and minimizer options with default values 
00803    ROOT::Fit::DataRange range(GetNdimensions()); 
00804    for ( int i = 0; i < GetNdimensions(); ++i ) {
00805       TAxis *axis = GetAxis(i);
00806       range.AddRange(i, axis->GetXmin(), axis->GetXmax());
00807    }
00808    ROOT::Math::MinimizerOptions minOption; 
00809    
00810    return ROOT::Fit::FitObject(this, f , fitOption , minOption, goption, range); 
00811 }
00812 
00813 //______________________________________________________________________________
00814 Long64_t THnSparse::GetBin(const Double_t* x, Bool_t allocate /* = kTRUE */)
00815 {
00816    // Get the bin index for the n dimensional tuple x,
00817    // allocate one if it doesn't exist yet and "allocate" is true.
00818 
00819    THnSparseCompactBinCoord* cc = GetCompactCoord();
00820    Int_t *coord = cc->GetCoord();
00821    for (Int_t i = 0; i < fNdimensions; ++i)
00822       coord[i] = GetAxis(i)->FindBin(x[i]);
00823    cc->UpdateCoord();
00824 
00825    return GetBinIndexForCurrentBin(allocate);
00826 }
00827 
00828 
00829 //______________________________________________________________________________
00830 Long64_t THnSparse::GetBin(const char* name[], Bool_t allocate /* = kTRUE */)
00831 {
00832    // Get the bin index for the n dimensional tuple addressed by "name",
00833    // allocate one if it doesn't exist yet and "allocate" is true.
00834 
00835    THnSparseCompactBinCoord* cc = GetCompactCoord();
00836    Int_t *coord = cc->GetCoord();
00837    for (Int_t i = 0; i < fNdimensions; ++i)
00838       coord[i] = GetAxis(i)->FindBin(name[i]);
00839    cc->UpdateCoord();
00840 
00841    return GetBinIndexForCurrentBin(allocate);
00842 }
00843 
00844 //______________________________________________________________________________
00845 Long64_t THnSparse::GetBin(const Int_t* coord, Bool_t allocate /*= kTRUE*/)
00846 {
00847    // Get the bin index for the n dimensional coordinates coord,
00848    // allocate one if it doesn't exist yet and "allocate" is true.
00849    GetCompactCoord()->SetCoord(coord);
00850    return GetBinIndexForCurrentBin(allocate);
00851 }
00852 
00853 //______________________________________________________________________________
00854 Double_t THnSparse::GetBinContent(const Int_t *coord) const {
00855    // Get content of bin with coordinates "coord"
00856    GetCompactCoord()->SetCoord(coord);
00857    Long64_t idx = const_cast<THnSparse*>(this)->GetBinIndexForCurrentBin(kFALSE);
00858    if (idx < 0) return 0.;
00859    THnSparseArrayChunk* chunk = GetChunk(idx / fChunkSize);
00860    return chunk->fContent->GetAt(idx % fChunkSize);
00861 }
00862 
00863 //______________________________________________________________________________
00864 Double_t THnSparse::GetBinContent(Long64_t idx, Int_t* coord /* = 0 */) const
00865 {
00866    // Return the content of the filled bin number "idx".
00867    // If coord is non-null, it will contain the bin's coordinates for each axis
00868    // that correspond to the bin.
00869 
00870    if (idx >= 0) {
00871       THnSparseArrayChunk* chunk = GetChunk(idx / fChunkSize);
00872       idx %= fChunkSize;
00873       if (chunk && chunk->fContent->GetSize() > idx) {
00874          if (coord) {
00875             THnSparseCompactBinCoord* cc = GetCompactCoord();
00876             Int_t sizeCompact = cc->GetBufferSize();
00877             cc->SetCoordFromBuffer(chunk->fCoordinates + idx * sizeCompact,
00878                                                   coord);
00879                                                                                
00880          }
00881          return chunk->fContent->GetAt(idx);
00882       }
00883    }
00884    if (coord)
00885       memset(coord, -1, sizeof(Int_t) * fNdimensions);
00886    return 0.;
00887 }
00888 
00889 //______________________________________________________________________________
00890 Double_t THnSparse::GetBinError(const Int_t *coord) const {
00891    // Get error of bin with coordinates "coord" as
00892    // BEGIN_LATEX #sqrt{#sum weight^{2}}
00893    // END_LATEX
00894    // If errors are not enabled (via Sumw2() or CalculateErrors())
00895    // return sqrt(contents).
00896 
00897    if (!GetCalculateErrors())
00898       return TMath::Sqrt(GetBinContent(coord));
00899 
00900    GetCompactCoord()->SetCoord(coord);
00901    Long64_t idx = const_cast<THnSparse*>(this)->GetBinIndexForCurrentBin(kFALSE);
00902    if (idx < 0) return 0.;
00903 
00904    THnSparseArrayChunk* chunk = GetChunk(idx / fChunkSize);
00905    return TMath::Sqrt(chunk->fSumw2->GetAt(idx % fChunkSize));
00906 }
00907 
00908 //______________________________________________________________________________
00909 Double_t THnSparse::GetBinError(Long64_t linidx) const {
00910    // Get error of bin addressed by linidx as
00911    // BEGIN_LATEX #sqrt{#sum weight^{2}}
00912    // END_LATEX
00913    // If errors are not enabled (via Sumw2() or CalculateErrors())
00914    // return sqrt(contents).
00915 
00916    if (!GetCalculateErrors())
00917       return TMath::Sqrt(GetBinContent(linidx));
00918 
00919    if (linidx < 0) return 0.;
00920    THnSparseArrayChunk* chunk = GetChunk(linidx / fChunkSize);
00921    linidx %= fChunkSize;
00922    if (!chunk || chunk->fContent->GetSize() < linidx)
00923       return 0.;
00924 
00925    return TMath::Sqrt(chunk->fSumw2->GetAt(linidx));
00926 }
00927 
00928 
00929 //______________________________________________________________________________
00930 Long64_t THnSparse::GetBinIndexForCurrentBin(Bool_t allocate)
00931 {
00932    // Return the index for fCurrentBinIndex.
00933    // If it doesn't exist then return -1, or allocate a new bin if allocate is set
00934 
00935    THnSparseCompactBinCoord* cc = GetCompactCoord();
00936    ULong64_t hash = cc->GetHash();
00937    if (fBinContent.GetSize() && !fBins.GetSize())
00938       FillExMap();
00939    Long64_t linidx = (Long64_t) fBins.GetValue(hash);
00940    while (linidx) {
00941       // fBins stores index + 1!
00942       THnSparseArrayChunk* chunk = GetChunk((linidx - 1)/ fChunkSize);
00943       if (chunk->Matches((linidx - 1) % fChunkSize, cc->GetBuffer()))
00944          return linidx - 1; // we store idx+1, 0 is "TExMap: not found"
00945 
00946       Long64_t nextlinidx = fBinsContinued.GetValue(linidx);
00947       if (!nextlinidx) break;
00948 
00949       linidx = nextlinidx;
00950    }
00951    if (!allocate) return -1;
00952 
00953    ++fFilledBins;
00954 
00955    // allocate bin in chunk
00956    THnSparseArrayChunk *chunk = (THnSparseArrayChunk*) fBinContent.Last();
00957    Long64_t newidx = chunk ? ((Long64_t) chunk->GetEntries()) : -1;
00958    if (!chunk || newidx == (Long64_t)fChunkSize) {
00959       chunk = AddChunk();
00960       newidx = 0;
00961    }
00962    chunk->AddBin(newidx, cc->GetBuffer());
00963 
00964    // store translation between hash and bin
00965    newidx += (fBinContent.GetEntriesFast() - 1) * fChunkSize;
00966    if (!linidx) {
00967       // fBins didn't find it
00968       if (2 * GetNbins() > fBins.Capacity())
00969          fBins.Expand(3 * GetNbins());
00970       fBins.Add(hash, newidx + 1);
00971    } else {
00972       // fBins contains one, but it's the wrong one;
00973       // add entry to fBinsContinued.
00974       fBinsContinued.Add(linidx, newidx + 1);
00975    }
00976    return newidx;
00977 }
00978 
00979 //______________________________________________________________________________
00980 THnSparseCompactBinCoord* THnSparse::GetCompactCoord() const
00981 {
00982    // Return THnSparseCompactBinCoord object.
00983 
00984    if (!fCompactCoord) {
00985       Int_t *bins = new Int_t[fNdimensions];
00986       for (Int_t d = 0; d < fNdimensions; ++d)
00987          bins[d] = GetAxis(d)->GetNbins();
00988       const_cast<THnSparse*>(this)->fCompactCoord
00989          = new THnSparseCompactBinCoord(fNdimensions, bins);
00990       delete [] bins;
00991    }
00992    return fCompactCoord;
00993 }
00994 
00995 //______________________________________________________________________________
00996 void THnSparse::GetRandom(Double_t *rand, Bool_t subBinRandom /* = kTRUE */)
00997 {
00998    // Generate an n-dimensional random tuple based on the histogrammed
00999    // distribution. If subBinRandom, the returned tuple will be additionally
01000    // randomly distributed within the randomized bin, using a flat
01001    // distribution.
01002 
01003    // check whether the integral array is valid
01004    if (fIntegralStatus != kValidInt)
01005       ComputeIntegral();
01006 
01007    // generate a random bin
01008    Double_t p = gRandom->Rndm();
01009    Long64_t idx = TMath::BinarySearch(GetNbins() + 1, fIntegral, p);
01010    const Int_t nStaticBins = 40;
01011    Int_t bin[nStaticBins];
01012    Int_t* pBin = bin;
01013    if (GetNdimensions() > nStaticBins) {
01014       pBin = new Int_t[GetNdimensions()];
01015    }
01016    GetBinContent(idx, pBin);
01017 
01018    // convert bin coordinates to real values
01019    for (Int_t i = 0; i < fNdimensions; i++) {
01020       rand[i] = GetAxis(i)->GetBinCenter(pBin[i]);
01021 
01022       // randomize the vector withing a bin
01023       if (subBinRandom)
01024          rand[i] += (gRandom->Rndm() - 0.5) * GetAxis(i)->GetBinWidth(pBin[i]);
01025    }
01026    if (pBin != bin) {
01027       delete [] pBin;
01028    }
01029 
01030    return;
01031 }
01032 
01033 //______________________________________________________________________________
01034 Double_t THnSparse::GetSparseFractionBins() const {
01035    // Return the amount of filled bins over all bins
01036 
01037    Double_t nbinsTotal = 1.;
01038    for (Int_t d = 0; d < fNdimensions; ++d)
01039       nbinsTotal *= GetAxis(d)->GetNbins() + 2;
01040    return fFilledBins / nbinsTotal;
01041 }
01042 
01043 //______________________________________________________________________________
01044 Double_t THnSparse::GetSparseFractionMem() const {
01045    // Return the amount of used memory over memory that would be used by a
01046    // non-sparse n-dimensional histogram. The value is approximate.
01047 
01048    Int_t arrayElementSize = 0;
01049    if (fFilledBins) {
01050       TClass* clArray = GetChunk(0)->fContent->IsA();
01051       TDataMember* dm = clArray ? clArray->GetDataMember("fArray") : 0;
01052       arrayElementSize = dm ? dm->GetDataType()->Size() : 0;
01053    }
01054    if (!arrayElementSize) {
01055       Warning("GetSparseFractionMem", "Cannot determine type of elements!");
01056       return -1.;
01057    }
01058 
01059    Double_t sizePerChunkElement = arrayElementSize + GetCompactCoord()->GetBufferSize();
01060    if (fFilledBins && GetChunk(0)->fSumw2)
01061       sizePerChunkElement += sizeof(Double_t); /* fSumw2 */
01062 
01063    Double_t size = 0.;
01064    size += fBinContent.GetEntries() * (GetChunkSize() * sizePerChunkElement + sizeof(THnSparseArrayChunk));
01065    size += + 3 * sizeof(Long64_t) * fBins.GetSize() /* TExMap */;
01066 
01067    Double_t nbinsTotal = 1.;
01068    for (Int_t d = 0; d < fNdimensions; ++d)
01069       nbinsTotal *= GetAxis(d)->GetNbins() + 2;
01070 
01071    return size / nbinsTotal / arrayElementSize;
01072 }
01073 
01074 //______________________________________________________________________________
01075 Bool_t THnSparse::IsInRange(Int_t *coord) const
01076 {
01077    // Check whether bin coord is in range, as defined by TAxis::SetRange().
01078    // Currently, TAxis::SetRange() does not allow to select all but over- and
01079    // underflow bins (it instead resets the axis to "no range selected").
01080    // Instead, simply call
01081    //    TAxis* axis12 = hsparse.GetAxis(12);
01082    //    axis12->SetRange(1, axis12->GetNbins());
01083    //    axis12->SetBit(TAxis::kAxisRange);
01084    // to deselect the under- and overflow bins in the 12th dimension.
01085 
01086    Int_t min = 0;
01087    Int_t max = 0;
01088    for (Int_t i = 0; i < fNdimensions; ++i) {
01089       TAxis *axis = GetAxis(i);
01090       if (!axis->TestBit(TAxis::kAxisRange)) continue;
01091       min = axis->GetFirst();
01092       max = axis->GetLast();
01093       if (min == 0 && max == 0) {
01094          // special case where TAxis::SetBit(kAxisRange) and
01095          // over- and underflow bins are de-selected.
01096          // first and last are == 0 due to axis12->SetRange(1, axis12->GetNbins());
01097          min = 1;
01098          max = axis->GetNbins();
01099       }
01100       if (coord[i] < min || coord[i] > max)
01101          return kFALSE;
01102    }
01103    return kTRUE;
01104 }
01105 
01106 //______________________________________________________________________________
01107 TH1D* THnSparse::Projection(Int_t xDim, Option_t* option /*= ""*/) const
01108 {
01109    // Project all bins into a 1-dimensional histogram,
01110    // keeping only axis "xDim".
01111    // If "option" contains "E" errors will be calculated.
01112    //                      "A" ranges of the taget axes will be ignored.
01113    //                      "O" original axis range of the taget axes will be
01114    //                          kept, but only bins inside the selected range
01115    //                          will be filled.
01116 
01117    return (TH1D*) ProjectionAny(1, &xDim, false, option);
01118 }
01119 
01120 //______________________________________________________________________________
01121 TH2D* THnSparse::Projection(Int_t yDim, Int_t xDim, Option_t* option /*= ""*/) const
01122 {
01123    // Project all bins into a 2-dimensional histogram,
01124    // keeping only axes "xDim" and "yDim".
01125    //
01126    // WARNING: just like TH3::Project3D("yx") and TTree::Draw("y:x"),
01127    // Projection(y,x) uses the first argument to define the y-axis and the
01128    // second for the x-axis!
01129    //
01130    // If "option" contains "E" errors will be calculated.
01131    //                      "A" ranges of the taget axes will be ignored.
01132 
01133    const Int_t dim[2] = {xDim, yDim};
01134    return (TH2D*) ProjectionAny(2, dim, false, option);
01135 }
01136 
01137 //______________________________________________________________________________
01138 TH3D* THnSparse::Projection(Int_t xDim, Int_t yDim, Int_t zDim,
01139                             Option_t* option /*= ""*/) const
01140 {
01141    // Project all bins into a 3-dimensional histogram,
01142    // keeping only axes "xDim", "yDim", and "zDim".
01143    // If "option" contains "E" errors will be calculated.
01144    //                      "A" ranges of the taget axes will be ignored.
01145    //                      "O" original axis range of the taget axes will be
01146    //                          kept, but only bins inside the selected range
01147    //                          will be filled.
01148 
01149    const Int_t dim[3] = {xDim, yDim, zDim};
01150    return (TH3D*) ProjectionAny(3, dim, false, option);
01151 }
01152 
01153 //______________________________________________________________________________
01154 THnSparse* THnSparse::Projection(Int_t ndim, const Int_t* dim,
01155                                  Option_t* option /*= ""*/) const
01156 {
01157    // Project all bins into a ndim-dimensional THnSparse histogram,
01158    // keeping only axes in dim (specifying ndim dimensions)
01159    // If "option" contains "E" errors will be calculated.
01160    //                      "A" ranges of the taget axes will be ignored.
01161    //                      "O" original axis range of the taget axes will be
01162    //                          kept, but only bins inside the selected range
01163    //                          will be filled.
01164 
01165    return (THnSparse*) ProjectionAny(ndim, dim, true, option);
01166 }
01167 
01168 
01169 //______________________________________________________________________________
01170 TObject* THnSparse::ProjectionAny(Int_t ndim, const Int_t* dim,
01171                                   Bool_t wantSparse, Option_t* option /*= ""*/) const
01172 {
01173    // Project all bins into a ndim-dimensional THnSparse histogram,
01174    // keeping only axes in dim (specifying ndim dimensions)
01175    // If "option" contains "E" errors will be calculated.
01176    //                      "A" ranges of the taget axes will be ignored.
01177    //                      "O" original axis range of the taget axes will be
01178    //                          kept, but only bins inside the selected range
01179    //                          will be filled.
01180 
01181    TString name(GetName());
01182    name +="_proj";
01183    
01184    for (Int_t d = 0; d < ndim; ++d) {
01185      name += "_";
01186      name += dim[d];
01187    }
01188 
01189    TString title(GetTitle());
01190    Ssiz_t posInsert = title.First(';');
01191    if (posInsert == kNPOS) {
01192       title += " projection ";
01193       for (Int_t d = 0; d < ndim; ++d)
01194          title += GetAxis(dim[d])->GetTitle();
01195    } else {
01196       for (Int_t d = ndim - 1; d >= 0; --d) {
01197          title.Insert(posInsert, GetAxis(d)->GetTitle());
01198          if (d)
01199             title.Insert(posInsert, ", ");
01200       }
01201       title.Insert(posInsert, " projection ");
01202    }
01203 
01204    TObjArray newaxes(ndim);
01205    for (Int_t d = 0; d < ndim; ++d) {
01206       newaxes.AddAt(GetAxis(dim[d]),d);
01207    }
01208 
01209    THnSparse* sparse = 0;
01210    TH1* hist = 0;
01211    TObject* ret = 0;
01212 
01213    Bool_t* hadRange = 0;
01214    Bool_t ignoreTargetRange = (option && (strchr(option, 'A') || strchr(option, 'a')));
01215    Bool_t keepTargetAxis = ignoreTargetRange || (option && (strchr(option, 'O') || strchr(option, 'o')));
01216    if (ignoreTargetRange) {
01217       hadRange = new Bool_t[ndim];
01218       for (Int_t d = 0; d < ndim; ++d){
01219          TAxis *axis = GetAxis(dim[d]);
01220          hadRange[d] = axis->TestBit(TAxis::kAxisRange);
01221          axis->SetBit(TAxis::kAxisRange, kFALSE);
01222       }
01223    }
01224 
01225    if (wantSparse)
01226       ret = sparse = CloneEmpty(name, title, &newaxes, fChunkSize, keepTargetAxis);
01227    else
01228       ret = hist = CreateHist(name, title, &newaxes, keepTargetAxis); 
01229 
01230    if (keepTargetAxis) {
01231       // make the whole axes visible, i.e. unset the range
01232       if (wantSparse) {
01233          for (Int_t d = 0; d < ndim; ++d) {
01234             sparse->GetAxis(d)->SetRange(0, 0);
01235          }
01236       } else {
01237          hist->GetXaxis()->SetRange(0, 0);
01238          hist->GetYaxis()->SetRange(0, 0);
01239          hist->GetZaxis()->SetRange(0, 0);
01240       }
01241    }
01242 
01243    Bool_t haveErrors = GetCalculateErrors();
01244    Bool_t wantErrors = (option && (strchr(option, 'E') || strchr(option, 'e'))) || haveErrors;
01245 
01246    Int_t* bins  = new Int_t[ndim];
01247    Int_t  linbin = 0;
01248    Int_t* coord = new Int_t[fNdimensions];
01249    memset(coord, 0, sizeof(Int_t) * fNdimensions);
01250 
01251    Double_t err = 0.;
01252    Double_t preverr = 0.;
01253    Double_t v = 0.;
01254 
01255    for (Long64_t i = 0; i < GetNbins(); ++i) {
01256       v = GetBinContent(i, coord);
01257 
01258       if (!IsInRange(coord)) continue;
01259 
01260       for (Int_t d = 0; d < ndim; ++d) {
01261          bins[d] = coord[dim[d]];
01262          if (!keepTargetAxis && GetAxis(dim[d])->TestBit(TAxis::kAxisRange)) {
01263             bins[d] -= GetAxis(dim[d])->GetFirst() - 1;
01264          }
01265       }
01266 
01267       if (!wantSparse) {
01268          if (ndim == 1) linbin = bins[0];
01269          else if (ndim == 2) linbin = hist->GetBin(bins[0], bins[1]);
01270          else if (ndim == 3) linbin = hist->GetBin(bins[0], bins[1], bins[2]);
01271       }
01272 
01273       if (wantErrors) {
01274          if (haveErrors) {
01275             err = GetBinError(i);
01276             err *= err;
01277          } else err = v;
01278          if (wantSparse) {
01279             preverr = sparse->GetBinError(bins);
01280             sparse->SetBinError(bins, TMath::Sqrt(preverr * preverr + err));
01281          } else {
01282             preverr = hist->GetBinError(linbin);
01283             hist->SetBinError(linbin, TMath::Sqrt(preverr * preverr + err));
01284          }
01285       }
01286 
01287       // only _after_ error calculation, or sqrt(v) is taken into account!
01288       if (wantSparse)
01289          sparse->AddBinContent(bins, v);
01290       else
01291          hist->AddBinContent(linbin, v);
01292    }
01293 
01294    delete [] bins;
01295    delete [] coord;
01296 
01297    if (wantSparse)
01298       // need to check also when producing a THnSparse how to reset the number of entries
01299       sparse->SetEntries(fEntries);
01300    else  
01301       // need to reset the statistics which will set also the number of entries
01302       // according to the bins filled
01303       hist->ResetStats(); 
01304 
01305    if (hadRange) {
01306       // reset kAxisRange bit:
01307       for (Int_t d = 0; d < ndim; ++d)
01308          GetAxis(dim[d])->SetBit(TAxis::kAxisRange, hadRange[d]);
01309 
01310       delete [] hadRange;
01311    }
01312 
01313    return ret;
01314 }
01315 
01316 //______________________________________________________________________________
01317 void THnSparse::Scale(Double_t c)
01318 {
01319    // Scale contents and errors of this histogram by c:
01320    // this = this * c
01321    // It does not modify the histogram's number of entries.
01322 
01323 
01324    Double_t nEntries = GetEntries();
01325    // Scale the contents & errors
01326    Bool_t haveErrors = GetCalculateErrors();
01327    for (Long64_t i = 0; i < GetNbins(); ++i) {
01328       // Get the content of the bin from the current histogram
01329       Double_t v = GetBinContent(i);
01330       SetBinContent(i, c * v);
01331       if (haveErrors) {
01332          Double_t err = GetBinError(i);
01333          SetBinError(i, c * err);
01334       }
01335    }
01336    SetEntries(nEntries);
01337 }
01338 
01339 //______________________________________________________________________________
01340 void THnSparse::AddInternal(const THnSparse* h, Double_t c, Bool_t rebinned)
01341 {
01342    // Add() implementation for both rebinned histograms and those with identical
01343    // binning. See THnSparse::Add().
01344 
01345    if (fNdimensions != h->GetNdimensions()) {
01346       Warning("RebinnedAdd", "Different number of dimensions, cannot carry out operation on the histograms");
01347       return;
01348    }
01349 
01350    // Trigger error calculation if h has it
01351    if (!GetCalculateErrors() && h->GetCalculateErrors())
01352       Sumw2();
01353    Bool_t haveErrors = GetCalculateErrors();
01354 
01355    // Now add the contents: in this case we have the union of the sets of bins
01356    Int_t* coord = new Int_t[fNdimensions];
01357    memset(coord, 0, sizeof(Int_t) * fNdimensions);
01358    Double_t* x = 0;
01359    if (rebinned) {
01360       x = new Double_t[fNdimensions];
01361    }
01362 
01363    if (!fBins.GetSize() && fBinContent.GetSize()) {
01364       FillExMap();
01365    }
01366 
01367    // Expand the exmap if needed, to reduce collisions
01368    Long64_t numTargetBins = GetNbins() + h->GetNbins();
01369    if (2 * numTargetBins > fBins.Capacity()) {
01370       fBins.Expand(3 * numTargetBins);
01371    }
01372 
01373    // Add to this whatever is found inside the other histogram
01374    for (Long64_t i = 0; i < h->GetNbins(); ++i) {
01375       // Get the content of the bin from the second histogram
01376       Double_t v = h->GetBinContent(i, coord);
01377 
01378       Long64_t binidx = -1;
01379       if (rebinned) {
01380          // Get the bin center given a coord
01381          for (Int_t j = 0; j < fNdimensions; ++j)
01382             x[j] = h->GetAxis(j)->GetBinCenter(coord[j]);
01383 
01384          binidx = GetBin(x);
01385       } else {
01386          binidx = GetBin(coord);
01387       }
01388 
01389       if (haveErrors) {
01390          Double_t err1 = GetBinError(binidx);
01391          Double_t err2 = h->GetBinError(i) * c;
01392          SetBinError(binidx, TMath::Sqrt(err1 * err1 + err2 * err2));
01393       }
01394       // only _after_ error calculation, or sqrt(v) is taken into account!
01395       AddBinContent(binidx, c * v);
01396    }
01397 
01398 
01399    delete [] coord;
01400    delete [] x;
01401 
01402    Double_t nEntries = GetEntries() + c * h->GetEntries();
01403    SetEntries(nEntries);
01404 }
01405 
01406 //______________________________________________________________________________
01407 void THnSparse::Add(const THnSparse* h, Double_t c)
01408 {
01409    // Add contents of h scaled by c to this histogram:
01410    // this = this + c * h
01411    // Note that if h has Sumw2 set, Sumw2 is automatically called for this
01412    // if not already set.
01413 
01414    // Check consistency of the input
01415    if (!CheckConsistency(h, "Add")) return;
01416 
01417    AddInternal(h, c, kFALSE);
01418 }
01419 
01420 //______________________________________________________________________________
01421 void THnSparse::RebinnedAdd(const THnSparse* h, Double_t c)
01422 {
01423    // Add contents of h scaled by c to this histogram:
01424    // this = this + c * h
01425    // Note that if h has Sumw2 set, Sumw2 is automatically called for this
01426    // if not already set.
01427    // In contrast to Add(), RebinnedAdd() does not require consist binning of
01428    // this and h; instead, each bin's center is used to determine the target bin.
01429 
01430    AddInternal(h, c, kTRUE);
01431 }
01432 
01433 
01434 //______________________________________________________________________________
01435 Long64_t THnSparse::Merge(TCollection* list)
01436 {
01437    // Merge this with a list of THnSparses. All THnSparses provided
01438    // in the list must have the same bin layout!
01439 
01440    if (!list) return 0;
01441    if (list->IsEmpty()) return (Long64_t)GetEntries();
01442 
01443    Long64_t sumNbins = GetNbins();
01444    TIter iter(list);
01445    const TObject* addMeObj = 0;
01446    while ((addMeObj = iter())) {
01447       const THnSparse* addMe = dynamic_cast<const THnSparse*>(addMeObj);
01448       if (addMe) {
01449          sumNbins += addMe->GetNbins();
01450       }
01451    }
01452    if (!fBins.GetSize() && fBinContent.GetSize()) {
01453       FillExMap();
01454    }
01455    if (2 * sumNbins > fBins.Capacity()) {
01456       fBins.Expand(3 * sumNbins);
01457    }
01458 
01459    iter.Reset();
01460    while ((addMeObj = iter())) {
01461       const THnSparse* addMe = dynamic_cast<const THnSparse*>(addMeObj);
01462       if (!addMe)
01463          Error("Merge", "Object named %s is not THnSpase! Skipping it.",
01464                addMeObj->GetName());
01465       else
01466          Add(addMe);
01467    }
01468    return (Long64_t)GetEntries();
01469 }
01470 
01471 
01472 //______________________________________________________________________________
01473 void THnSparse::Multiply(const THnSparse* h)
01474 {
01475    // Multiply this histogram by histogram h
01476    // this = this * h
01477    // Note that if h has Sumw2 set, Sumw2 is automatically called for this
01478    // if not already set.
01479 
01480    // Check consistency of the input
01481    if(!CheckConsistency(h, "Multiply"))return;
01482 
01483    // Trigger error calculation if h has it
01484    Bool_t wantErrors = kFALSE;
01485    if (GetCalculateErrors() || h->GetCalculateErrors())
01486       wantErrors = kTRUE;
01487 
01488    if (wantErrors) Sumw2();
01489 
01490    Double_t nEntries = GetEntries();
01491    // Now multiply the contents: in this case we have the intersection of the sets of bins
01492    Int_t* coord = new Int_t[fNdimensions];
01493    for (Long64_t i = 0; i < GetNbins(); ++i) {
01494       // Get the content of the bin from the current histogram
01495       Double_t v1 = GetBinContent(i, coord);
01496       // Now look at the bin with the same coordinates in h
01497       Double_t v2 = h->GetBinContent(coord);
01498       SetBinContent(coord, v1 * v2);;
01499       if (wantErrors) {
01500          Double_t err1 = GetBinError(coord) * v2;
01501          Double_t err2 = h->GetBinError(coord) * v1;
01502          SetBinError(coord,TMath::Sqrt((err2 * err2 + err1 * err1)));
01503       }
01504    }
01505    SetEntries(nEntries);
01506 
01507    delete [] coord;
01508 }
01509 
01510 //______________________________________________________________________________
01511 void THnSparse::Multiply(TF1* f, Double_t c)
01512 {
01513     // Performs the operation: this = this*c*f1
01514     // if errors are defined, errors are also recalculated.
01515     //
01516     // Only bins inside the function range are recomputed.
01517     // IMPORTANT NOTE: If you intend to use the errors of this histogram later
01518     // you should call Sumw2 before making this operation.
01519     // This is particularly important if you fit the histogram after THnSparse::Multiply
01520 
01521    Int_t* coord = new Int_t[fNdimensions];
01522    memset(coord, 0, sizeof(Int_t) * fNdimensions);
01523    Double_t* points = new Double_t[fNdimensions];
01524 
01525 
01526    Bool_t wantErrors = GetCalculateErrors();
01527    if (wantErrors) Sumw2();
01528 
01529    ULong64_t nEntries = GetNbins();
01530    for (ULong64_t i = 0; i < nEntries; ++i) {
01531       Double_t value = GetBinContent(i, coord);
01532 
01533       // Get the bin co-ordinates given an coord
01534       for (Int_t j = 0; j < fNdimensions; ++j)
01535          points[j] = GetAxis(j)->GetBinCenter(coord[j]);
01536 
01537       if (!f->IsInside(points))
01538          continue;
01539       TF1::RejectPoint(kFALSE);
01540 
01541       // Evaulate function at points
01542       Double_t fvalue = f->EvalPar(points, NULL);
01543 
01544       SetBinContent(coord, c * fvalue * value);
01545       if (wantErrors) {
01546          Double_t error(GetBinError(i));
01547          SetBinError(coord, c * fvalue * error);
01548       }
01549    }
01550    
01551    delete [] points;
01552    delete [] coord;
01553 }
01554 
01555 //______________________________________________________________________________
01556 void THnSparse::Divide(const THnSparse *h)
01557 {
01558    // Divide this histogram by h
01559    // this = this/(h)
01560    // Note that if h has Sumw2 set, Sumw2 is automatically called for
01561    // this if not already set.
01562    // The resulting errors are calculated assuming uncorrelated content.
01563 
01564    // Check consistency of the input
01565    if (!CheckConsistency(h, "Divide"))return;
01566 
01567    // Trigger error calculation if h has it
01568    Bool_t wantErrors=GetCalculateErrors();
01569    if (!GetCalculateErrors() && h->GetCalculateErrors())
01570       wantErrors=kTRUE;
01571 
01572    // Remember original histogram statistics
01573    Double_t nEntries = fEntries;
01574 
01575    if (wantErrors) Sumw2();
01576    Bool_t didWarn = kFALSE;
01577 
01578    // Now divide the contents: also in this case we have the intersection of the sets of bins
01579    Int_t* coord = new Int_t[fNdimensions];
01580    memset(coord, 0, sizeof(Int_t) * fNdimensions);
01581    Double_t err = 0.;
01582    Double_t b22 = 0.;
01583    for (Long64_t i = 0; i < GetNbins(); ++i) {
01584       // Get the content of the bin from the first histogram
01585       Double_t v1 = GetBinContent(i, coord);
01586       // Now look at the bin with the same coordinates in h
01587       Double_t v2 = h->GetBinContent(coord);
01588       if (!v2) {
01589          v1 = 0.;
01590          v2 = 1.;
01591          if (!didWarn) {
01592             Warning("Divide(h)", "Histogram h has empty bins - division by zero! Setting bin to 0.");
01593             didWarn = kTRUE;
01594          }
01595       }
01596       SetBinContent(coord, v1 / v2);
01597       if (wantErrors) {
01598          Double_t err1 = GetBinError(coord) * v2;
01599          Double_t err2 = h->GetBinError(coord) * v1;
01600          b22 = v2 * v2;
01601          err = (err1 * err1 + err2 * err2) / (b22 * b22);
01602          SetBinError(coord, TMath::Sqrt(err));
01603       }
01604    }
01605    delete [] coord;
01606    SetEntries(nEntries);
01607 }
01608 
01609 //______________________________________________________________________________
01610 void THnSparse::Divide(const THnSparse *h1, const THnSparse *h2, Double_t c1, Double_t c2, Option_t *option)
01611 {
01612    // Replace contents of this histogram by multiplication of h1 by h2
01613    // this = (c1*h1)/(c2*h2)
01614    // Note that if h1 or h2 have Sumw2 set, Sumw2 is automatically called for
01615    // this if not already set.
01616    // The resulting errors are calculated assuming uncorrelated content.
01617    // However, if option ="B" is specified, Binomial errors are computed.
01618    // In this case c1 and c2 do not make real sense and they are ignored.
01619 
01620 
01621    TString opt = option;
01622    opt.ToLower();
01623    Bool_t binomial = kFALSE;
01624    if (opt.Contains("b")) binomial = kTRUE;
01625 
01626    // Check consistency of the input
01627    if (!CheckConsistency(h1, "Divide") || !CheckConsistency(h2, "Divide"))return;
01628    if (!c2) {
01629       Error("Divide","Coefficient of dividing histogram cannot be zero");
01630       return;
01631    }
01632 
01633    Reset();
01634 
01635    // Trigger error calculation if h1 or h2 have it
01636    if (!GetCalculateErrors() && (h1->GetCalculateErrors()|| h2->GetCalculateErrors() != 0))
01637       Sumw2();
01638 
01639    // Count filled bins
01640    Long64_t nFilledBins=0;
01641 
01642    // Now divide the contents: we have the intersection of the sets of bins
01643 
01644    Int_t* coord = new Int_t[fNdimensions];
01645    memset(coord, 0, sizeof(Int_t) * fNdimensions);
01646    Float_t w = 0.;
01647    Float_t err = 0.;
01648    Float_t b22 = 0.;
01649    Bool_t didWarn = kFALSE;
01650 
01651    for (Long64_t i = 0; i < h1->GetNbins(); ++i) {
01652       // Get the content of the bin from the first histogram
01653       Double_t v1 = h1->GetBinContent(i, coord);
01654       // Now look at the bin with the same coordinates in h2
01655       Double_t v2 = h2->GetBinContent(coord);
01656       if (!v2) {
01657          v1 = 0.;
01658          v2 = 1.;
01659          if (!didWarn) {
01660             Warning("Divide(h1, h2)", "Histogram h2 has empty bins - division by zero! Setting bin to 0.");
01661             didWarn = kTRUE;
01662          }
01663       }
01664       nFilledBins++;
01665       SetBinContent(coord, c1 * v1 / c2 / v2);
01666       if(GetCalculateErrors()){
01667          Double_t err1=h1->GetBinError(coord);
01668          Double_t err2=h2->GetBinError(coord);
01669          if (binomial) {
01670             if (v1 != v2) {
01671                w = v1 / v2;
01672                err2 *= w;
01673                err = TMath::Abs( ( (1. - 2.*w) * err1 * err1 + err2 * err2 ) / (v2 * v2) );
01674             } else {
01675                err = 0;
01676             }
01677          } else {
01678             c1 *= c1;
01679             c2 *= c2;
01680             b22 = v2 * v2 * c2;
01681             err1 *= v2;
01682             err2 *= v1;
01683             err = c1 * c2 * (err1 * err1 + err2 * err2) / (b22 * b22);
01684          }
01685          SetBinError(coord,TMath::Sqrt(err));
01686       }
01687    }
01688 
01689    delete [] coord;
01690    fFilledBins = nFilledBins;
01691 
01692    // Set as entries in the result histogram the entries in the numerator
01693    SetEntries(h1->GetEntries());
01694 }
01695 
01696 //______________________________________________________________________________
01697 Bool_t THnSparse::CheckConsistency(const THnSparse *h, const char *tag) const
01698 {
01699    // Consistency check on (some of) the parameters of two histograms (for operations).
01700 
01701    if (fNdimensions!=h->GetNdimensions()) {
01702       Warning(tag,"Different number of dimensions, cannot carry out operation on the histograms");
01703       return kFALSE;
01704    }
01705    for (Int_t dim = 0; dim < fNdimensions; dim++){
01706       if (GetAxis(dim)->GetNbins()!=h->GetAxis(dim)->GetNbins()) {
01707          Warning(tag,"Different number of bins on axis %i, cannot carry out operation on the histograms", dim);
01708          return kFALSE;
01709       }
01710    }
01711    return kTRUE;
01712 }
01713 
01714 //______________________________________________________________________________
01715 void THnSparse::SetBinEdges(Int_t idim, const Double_t* bins)
01716 {
01717    // Set the axis # of bins and bin limits on dimension idim
01718 
01719    TAxis* axis = (TAxis*) fAxes[idim];
01720    axis->Set(axis->GetNbins(), bins);
01721 }
01722 
01723 //______________________________________________________________________________
01724 void THnSparse::SetBinContent(const Int_t* coord, Double_t v)
01725 {
01726    // Set content of bin with coordinates "coord" to "v"
01727 
01728    GetCompactCoord()->SetCoord(coord);
01729    Long_t bin = GetBinIndexForCurrentBin(kTRUE);
01730    SetBinContent(bin, v);
01731 }
01732 
01733 //______________________________________________________________________________
01734 void THnSparse::SetBinContent(Long64_t bin, Double_t v)
01735 {
01736    // Set content of bin with index "bin" to "v"
01737 
01738    THnSparseArrayChunk* chunk = GetChunk(bin / fChunkSize);
01739    chunk->fContent->SetAt(v, bin % fChunkSize);
01740    ++fEntries;
01741 }
01742 
01743 //______________________________________________________________________________
01744 void THnSparse::SetBinError(const Int_t* coord, Double_t e)
01745 {
01746    // Set error of bin with coordinates "coord" to "e", enable errors if needed
01747 
01748    GetCompactCoord()->SetCoord(coord);
01749    Long_t bin = GetBinIndexForCurrentBin(kTRUE);
01750    SetBinError(bin, e);
01751 }
01752 
01753 //______________________________________________________________________________
01754 void THnSparse::SetBinError(Long64_t bin, Double_t e)
01755 {
01756    // Set error of bin with index "bin" to "e", enable errors if needed
01757 
01758    THnSparseArrayChunk* chunk = GetChunk(bin / fChunkSize);
01759    if (!chunk->fSumw2 ) {
01760       // if fSumw2 is zero GetCalculateErrors should return false
01761       if (GetCalculateErrors()) {
01762          Error("SetBinError", "GetCalculateErrors() logic error!");
01763       }
01764       Sumw2(); // enable error calculation
01765    }
01766 
01767    chunk->fSumw2->SetAt(e * e, bin % fChunkSize);
01768 }
01769 
01770 //______________________________________________________________________________
01771 void THnSparse::Sumw2()
01772 {
01773    // Enable calculation of errors
01774 
01775    if (GetCalculateErrors()) return;
01776 
01777    fTsumw2 = 0.;
01778    TIter iChunk(&fBinContent);
01779    THnSparseArrayChunk* chunk = 0;
01780    while ((chunk = (THnSparseArrayChunk*) iChunk()))
01781       chunk->Sumw2();
01782 }
01783 
01784 //______________________________________________________________________________
01785 THnSparse* THnSparse::Rebin(Int_t group) const
01786 {
01787    // Combine the content of "group" neighboring bins into
01788    // a new bin and return the resulting THnSparse.
01789    // For group=2 and a 3 dimensional histogram, all "blocks"
01790    // of 2*2*2 bins will be put into a bin.
01791 
01792    Int_t* ngroup = new Int_t[GetNdimensions()];
01793    for (Int_t d = 0; d < GetNdimensions(); ++d)
01794       ngroup[d] = group;
01795    THnSparse* ret = Rebin(ngroup);
01796    delete [] ngroup;
01797    return ret;
01798 }
01799 
01800 //______________________________________________________________________________
01801 void THnSparse::SetTitle(const char *title)
01802 {
01803    // Change (i.e. set) the title.
01804    //
01805    // If title is in the form "stringt;string0;string1;string2 ..."
01806    // the histogram title is set to stringt, the title of axis0 to string0,
01807    // of axis1 to string1, of axis2 to string2, etc, just like it is done
01808    // for TH1/TH2/TH3.
01809    // To insert the character ";" in one of the titles, one should use "#;"
01810    // or "#semicolon".
01811 
01812   fTitle = title;
01813   fTitle.ReplaceAll("#;",2,"#semicolon",10);
01814   
01815   Int_t endHistTitle = fTitle.First(';');
01816   if (endHistTitle >= 0) {
01817      // title contains a ';' so parse the axis titles
01818      Int_t posTitle = endHistTitle + 1;
01819      Int_t lenTitle = fTitle.Length();
01820      Int_t dim = 0;
01821      while (posTitle > 0 && posTitle < lenTitle && dim < fNdimensions){
01822         Int_t endTitle = fTitle.Index(";", posTitle);
01823         TString axisTitle = fTitle(posTitle, endTitle - posTitle);
01824         axisTitle.ReplaceAll("#semicolon", 10, ";", 1);
01825         GetAxis(dim)->SetTitle(axisTitle);
01826         dim++;
01827         if (endTitle > 0)
01828            posTitle = endTitle + 1;
01829         else
01830            posTitle = -1;
01831      }
01832      // Remove axis titles from histogram title
01833      fTitle.Remove(endHistTitle, lenTitle - endHistTitle);
01834   }
01835   
01836   fTitle.ReplaceAll("#semicolon", 10, ";", 1);
01837 
01838 }
01839 //______________________________________________________________________________
01840 THnSparse* THnSparse::Rebin(const Int_t* group) const
01841 {
01842    // Combine the content of "group" neighboring bins for each dimension
01843    // into a new bin and return the resulting THnSparse.
01844    // For group={2,1,1} and a 3 dimensional histogram, pairs of x-bins
01845    // will be grouped.
01846 
01847    Int_t ndim = GetNdimensions();
01848    TString name(GetName());
01849    for (Int_t d = 0; d < ndim; ++d)
01850       name += Form("_%d", group[d]);
01851 
01852 
01853    TString title(GetTitle());
01854    Ssiz_t posInsert = title.First(';');
01855    if (posInsert == kNPOS) {
01856       title += " rebin ";
01857       for (Int_t d = 0; d < ndim; ++d)
01858          title += Form("{%d}", group[d]);
01859    } else {
01860       for (Int_t d = ndim - 1; d >= 0; --d)
01861          title.Insert(posInsert, Form("{%d}", group[d]));
01862       title.Insert(posInsert, " rebin ");
01863    }
01864 
01865    TObjArray newaxes(ndim);
01866    newaxes.SetOwner();
01867    for (Int_t d = 0; d < ndim; ++d) {
01868       newaxes.AddAt(GetAxis(d)->Clone(),d);
01869       if (group[d] > 1) {
01870          TAxis* newaxis = (TAxis*) newaxes.At(d);
01871          Int_t newbins = (newaxis->GetNbins() + group[d] - 1) / group[d];
01872          if (newaxis->GetXbins() && newaxis->GetXbins()->GetSize()) {
01873             // variable bins
01874             Double_t *edges = new Double_t[newbins + 1];
01875             for (Int_t i = 0; i < newbins + 1; ++i)
01876                if (group[d] * i <= newaxis->GetNbins())
01877                   edges[i] = newaxis->GetXbins()->At(group[d] * i);
01878                else edges[i] = newaxis->GetXmax();
01879             newaxis->Set(newbins, edges);
01880             delete [] edges;
01881          } else {
01882             newaxis->Set(newbins, newaxis->GetXmin(), newaxis->GetXmax());
01883          }
01884       }
01885    }
01886 
01887    THnSparse* h = CloneEmpty(name.Data(), title.Data(), &newaxes, fChunkSize, kTRUE);
01888    Bool_t haveErrors = GetCalculateErrors();
01889    Bool_t wantErrors = haveErrors;
01890 
01891    Int_t* bins  = new Int_t[ndim];
01892    Int_t* coord = new Int_t[fNdimensions];
01893    memset(coord, 0, sizeof(Int_t) * fNdimensions);
01894    Double_t err = 0.;
01895    Double_t preverr = 0.;
01896    Double_t v = 0.;
01897 
01898    for (Long64_t i = 0; i < GetNbins(); ++i) {
01899       v = GetBinContent(i, coord);
01900       for (Int_t d = 0; d < ndim; ++d)
01901          bins[d] = TMath::CeilNint( (double) coord[d]/group[d] );
01902 
01903       if (wantErrors) {
01904          if (haveErrors) {
01905             err = GetBinError(i);
01906             err *= err;
01907          } else err = v;
01908          preverr = h->GetBinError(bins);
01909          h->SetBinError(bins, TMath::Sqrt(preverr * preverr + err));
01910       }
01911 
01912       // only _after_ error calculation, or sqrt(v) is taken into account!
01913       h->AddBinContent(bins, v);
01914    }
01915 
01916    delete [] bins;
01917    delete [] coord;
01918 
01919    h->SetEntries(fEntries);
01920 
01921    return h;
01922 
01923 }
01924 
01925 //______________________________________________________________________________
01926 void THnSparse::Reset(Option_t * /*option = ""*/)
01927 {
01928    // Clear the histogram
01929    fFilledBins = 0;
01930    fEntries = 0.;
01931    fTsumw = 0.;
01932    fTsumw2 = -1.;
01933    fBins.Delete();
01934    fBinsContinued.Clear();
01935    fBinContent.Delete();
01936    if (fIntegralStatus != kNoInt) {
01937       delete [] fIntegral;
01938       fIntegralStatus = kNoInt;
01939    }
01940 }
01941 
01942 //______________________________________________________________________________
01943 Double_t THnSparse::ComputeIntegral()
01944 {
01945    // Calculate the integral of the histogram
01946 
01947    // delete old integral
01948    if (fIntegralStatus != kNoInt) {
01949       delete [] fIntegral;
01950       fIntegralStatus = kNoInt;
01951    }
01952 
01953    // check number of bins
01954    if (GetNbins() == 0) {
01955       Error("ComputeIntegral", "The histogram must have at least one bin.");
01956       return 0.;
01957    }
01958 
01959    // allocate integral array
01960    fIntegral = new Double_t [GetNbins() + 1];
01961    fIntegral[0] = 0.;
01962 
01963    // fill integral array with contents of regular bins (non over/underflow)
01964    Int_t* coord = new Int_t[fNdimensions];
01965    for (Long64_t i = 0; i < GetNbins(); ++i) {
01966       Double_t v = GetBinContent(i, coord);
01967 
01968       // check whether the bin is regular
01969       bool regularBin = true;
01970       for (Int_t dim = 0; dim < fNdimensions; dim++)
01971          if (coord[dim] < 1 || coord[dim] > GetAxis(dim)->GetNbins()) {
01972             regularBin = false;
01973             break;
01974          }
01975 
01976       // if outlayer, count it with zero weight
01977       if (!regularBin) v = 0.;
01978 
01979       fIntegral[i + 1] = fIntegral[i] + v;
01980    }
01981    delete [] coord;
01982 
01983    // check sum of weights
01984    if (fIntegral[GetNbins()] == 0.) {
01985       Error("ComputeIntegral", "No hits in regular bins (non over/underflow).");
01986       delete [] fIntegral;
01987       return 0.;
01988    }
01989 
01990    // normalize the integral array
01991    for (Long64_t i = 0; i <= GetNbins(); ++i)
01992       fIntegral[i] = fIntegral[i] / fIntegral[GetNbins()];
01993 
01994    // set status to valid
01995    fIntegralStatus = kValidInt;
01996    return fIntegral[GetNbins()];
01997 }
01998 
01999 //______________________________________________________________________________
02000 void THnSparse::PrintBin(Long64_t idx, Option_t* options) const
02001 {
02002    // Print bin with linex index "idx".
02003    // For valid options see PrintBin(Long64_t idx, Int_t* bin, Option_t* options).
02004    Int_t* coord = new Int_t[fNdimensions];
02005    PrintBin(idx, coord, options);
02006    delete [] coord;
02007 }
02008 
02009 //______________________________________________________________________________
02010 Bool_t THnSparse::PrintBin(Long64_t idx, Int_t* bin, Option_t* options) const
02011 {
02012    // Print one bin. If "idx" is != -1 use that to determine the bin,
02013    // otherwise (if "idx" == -1) use the coordinate in "bin".
02014    // If "options" contains:
02015    //   '0': only print bins with an error or content != 0
02016    // Return whether the bin was printed (depends on options)
02017    
02018    Double_t v = -42;
02019    if (idx == -1) {
02020       idx = const_cast<THnSparse*>(this)->GetBin(bin, kFALSE);
02021       v = GetBinContent(idx);
02022    } else {
02023       v = GetBinContent(idx, bin);
02024    }
02025 
02026    if (options && strchr(options, '0') && v == 0.) {
02027       if (GetCalculateErrors()) {
02028          if (GetBinError(idx) == 0.) {
02029             // suppress zeros, and we have one.
02030             return kFALSE;
02031          }
02032       } else {
02033          // suppress zeros, and we have one.
02034          return kFALSE;
02035       }
02036    }
02037 
02038    TString coord;
02039    for (Int_t dim = 0; dim < fNdimensions; ++dim) {
02040       coord += bin[dim];
02041       coord += ',';
02042    }
02043    coord.Remove(coord.Length() - 1);
02044 
02045    if (GetCalculateErrors()) {
02046       Double_t err = 0.;
02047       if (idx != -1) {
02048          err = GetBinError(idx);
02049       }
02050       Printf("Bin at (%s) = %g (+/- %g)", coord.Data(), v, err);
02051    } else {
02052       Printf("Bin at (%s) = %g", coord.Data(), v);
02053    }
02054 
02055    return kTRUE;
02056 }
02057 
02058 //______________________________________________________________________________
02059 void THnSparse::PrintEntries(Long64_t from /*=0*/, Long64_t howmany /*=-1*/,
02060                             Option_t* options /*=0*/) const
02061 {
02062    // Print "howmany" entries starting at "from". If "howmany" is -1, print all.
02063    // If "options" contains:
02064    //   'x': print in the order of axis bins, i.e. (0,0,...,0), (0,0,...,1),...
02065    //   '0': only print bins with content != 0
02066 
02067    if (from < 0) from = 0;
02068    if (howmany == -1) howmany = GetNbins();
02069 
02070    Int_t* bin = new Int_t[fNdimensions];
02071 
02072    if (options && (strchr(options, 'x') || strchr(options, 'X'))) {
02073       Int_t* nbins = new Int_t[fNdimensions];
02074       for (Int_t dim = fNdimensions - 1; dim >= 0; --dim) {
02075          nbins[dim] = GetAxis(dim)->GetNbins();
02076          bin[dim] = from % nbins[dim];
02077          from /= nbins[dim];
02078       }
02079 
02080       for (Long64_t i = 0; i < howmany; ++i) {
02081          if (!PrintBin(-1, bin, options))
02082             ++howmany;
02083          // Advance to next bin:
02084          ++bin[fNdimensions - 1];
02085          for (Int_t dim = fNdimensions - 1; dim >= 0; --dim) {
02086             if (bin[dim] >= nbins[dim]) {
02087                bin[dim] = 0;
02088                if (dim > 0) {
02089                   ++bin[dim - 1];
02090                } else {
02091                   howmany = -1; // aka "global break"
02092                }
02093             }
02094          }
02095       }
02096       delete [] nbins;
02097    } else {
02098       for (Long64_t i = from; i < from + howmany; ++i) {
02099          if (!PrintBin(i, bin, options))
02100             ++howmany;
02101       }
02102    }
02103    delete [] bin;
02104 }
02105 
02106 //______________________________________________________________________________
02107 void THnSparse::Print(Option_t* options) const
02108 {
02109    // Print a THnSparse. If "option" contains:
02110    //   'a': print axis details
02111    //   'm': print memory usage
02112    //   's': print statistics
02113    //   'c': print its content, too (this can generate a LOT of output!)
02114    // Other options are forwarded to PrintEntries().
02115 
02116    Bool_t optAxis    = options && (strchr(options, 'A') || (strchr(options, 'a')));
02117    Bool_t optMem     = options && (strchr(options, 'M') || (strchr(options, 'm')));
02118    Bool_t optStat    = options && (strchr(options, 'S') || (strchr(options, 's')));
02119    Bool_t optContent = options && (strchr(options, 'C') || (strchr(options, 'c')));
02120 
02121    Printf("%s (*0x%lx): \"%s\" \"%s\"", IsA()->GetName(), (unsigned long)this, GetName(), GetTitle());
02122    Printf("  %d dimensions, %g entries in %lld filled bins", GetNdimensions(), GetEntries(), GetNbins());
02123 
02124    if (optAxis) {
02125       for (Int_t dim = 0; dim < fNdimensions; ++dim) {
02126          TAxis* axis = GetAxis(dim);
02127          Printf("    axis %d \"%s\": %d bins (%g..%g), %s bin sizes", dim,
02128                 axis->GetTitle(), axis->GetNbins(), axis->GetXmin(), axis->GetXmax(),
02129                 (axis->GetXbins() ? "variable" : "fixed"));
02130       }
02131    }
02132 
02133    if (optStat) {
02134       Printf("  %s error calculation", (GetCalculateErrors() ? "with" : "without"));
02135       if (GetCalculateErrors()) {
02136          Printf("    Sum(w)=%g, Sum(w^2)=%g", GetSumw(), GetSumw2());
02137          for (Int_t dim = 0; dim < fNdimensions; ++dim) {
02138             Printf("    axis %d: Sum(w*x)=%g, Sum(w*x^2)=%g", dim, GetSumwx(dim), GetSumwx2(dim));
02139          }
02140       }
02141    }
02142 
02143    if (optMem) {
02144       Printf("  coordinates stored in %d chunks of %d entries\n    %g of bins filled using %g of memory compared to an array",
02145              GetNChunks(), GetChunkSize(), GetSparseFractionBins(), GetSparseFractionMem());
02146    }
02147 
02148    if (optContent) {
02149       Printf("  BIN CONTENT:");
02150       PrintEntries(0, -1, options);
02151    }
02152 }
02153 

Generated on Tue Jul 5 14:24:20 2011 for ROOT_528-00b_version by  doxygen 1.5.1