#include "THnSparse.h"
#include "TArrayI.h"
#include "TAxis.h"
#include "TClass.h"
#include "TCollection.h"
#include "TDataMember.h"
#include "TDataType.h"
#include "TH1D.h"
#include "TH2D.h"
#include "TH3D.h"
#include "TInterpreter.h"
#include "TMath.h"
#include "TRandom.h"
class THnSparseCompactBinCoord {
public:
   THnSparseCompactBinCoord(Int_t dim, const Int_t* nbins);
   ~THnSparseCompactBinCoord();
   void SetCoord(const Int_t* coord) { memcpy(fCurrentBin, coord, sizeof(Int_t) * fNdimensions); }
   ULong64_t GetHash();
   Int_t GetSize() const { return fCoordBufferSize; }
   Int_t* GetCoord() const { return fCurrentBin; }
   Char_t* GetBuffer() const { return fCoordBuffer; }
   void GetCoordFromBuffer(Int_t* coord) const;
protected:
   Int_t GetNumBits(Int_t n) const {
      
      Int_t r = (n > 0);
      while (n/=2) ++r;
      return r;
   }
private:
   Int_t  fNdimensions;     
   Int_t *fBitOffsets;      
   Char_t *fCoordBuffer;     
   Int_t  fCoordBufferSize; 
   Int_t *fCurrentBin;      
};
THnSparseCompactBinCoord::THnSparseCompactBinCoord(Int_t dim, const Int_t* nbins):
   fNdimensions(dim), fBitOffsets(0), fCoordBuffer(0), fCoordBufferSize(0)
{
   
   
   fCurrentBin = new Int_t[dim];
   fBitOffsets = new Int_t[dim + 1];
   int shift = 0;
   for (Int_t i = 0; i < dim; ++i) {
      fBitOffsets[i] = shift;
      shift += GetNumBits(nbins[i] + 2);
   }
   fBitOffsets[dim] = shift;
   fCoordBufferSize = (shift + 7) / 8;
   fCoordBuffer = new Char_t[fCoordBufferSize];
}
THnSparseCompactBinCoord::~THnSparseCompactBinCoord()
{
   
   delete [] fBitOffsets;
   delete [] fCoordBuffer;
   delete [] fCurrentBin;
}
void THnSparseCompactBinCoord::GetCoordFromBuffer(Int_t* coord) const
{
   
   
   for (Int_t i = 0; i < fNdimensions; ++i) {
      const Int_t offset = fBitOffsets[i] / 8;
      Int_t shift = fBitOffsets[i] % 8;
      Int_t nbits = fBitOffsets[i + 1] - fBitOffsets[i];
      UChar_t* pbuf = (UChar_t*) fCoordBuffer + offset;
      coord[i] = *pbuf >> shift;
      Int_t subst = (Int_t) -1;
      subst = subst << nbits;
      nbits -= (8 - shift);
      shift = 8 - shift;
      for (Int_t n = 0; n * 8 < nbits; ++n) {
         ++pbuf;
         coord[i] += *pbuf << shift;
         shift += 8;
      }
      coord[i] &= ~subst;
   }
}
ULong64_t THnSparseCompactBinCoord::GetHash()
{
   
   memset(fCoordBuffer, 0, fCoordBufferSize);
   for (Int_t i = 0; i < fNdimensions; ++i) {
      const Int_t offset = fBitOffsets[i] / 8;
      const Int_t shift = fBitOffsets[i] % 8;
      ULong64_t val = fCurrentBin[i];
      Char_t* pbuf = fCoordBuffer + offset;
      *pbuf += 0xff & (val << shift);
      val = val >> (8 - shift);
      while (val) {
         ++pbuf;
         *pbuf += 0xff & val;
         val = val >> 8;
      }
   }
   
   
   
   
   
   
   
   
   
   switch (fCoordBufferSize) {
   case 1: return fCoordBuffer[0];
   case 2: return fCoordBuffer[0] + (fCoordBuffer[1] <<  8l);
   case 3: return fCoordBuffer[0] + (fCoordBuffer[1] <<  8l)
              + (fCoordBuffer[2] << 16l);
   case 4: return fCoordBuffer[0] + (fCoordBuffer[1] <<  8l)
              + (fCoordBuffer[2] << 16l) + (fCoordBuffer[3] << 24l);
   }
   return TMath::Hash(fCoordBuffer, fCoordBufferSize);
}
ClassImp(THnSparseArrayChunk);
THnSparseArrayChunk::THnSparseArrayChunk(Int_t coordsize, bool errors, TArray* cont):
      fCoordinateAllocationSize(-1), fSingleCoordinateSize(coordsize), fCoordinatesSize(0),
      fCoordinates(0), fContent(cont),
      fSumw2(0)
{
   
   
   fCoordinateAllocationSize = fSingleCoordinateSize * cont->GetSize();
   fCoordinates = new Char_t[fCoordinateAllocationSize];
   if (errors) Sumw2();
}
THnSparseArrayChunk::~THnSparseArrayChunk()
{
   
   delete fContent;
   delete [] fCoordinates;
   delete fSumw2;
}
void THnSparseArrayChunk::AddBin(ULong_t idx, const Char_t* coordbuf)
{
   
   
   
   
   
   
   
   
   
   
   if (fCoordinateAllocationSize == -1 && fContent) {
      Int_t chunksize = fSingleCoordinateSize * fContent->GetSize();
      if (fCoordinatesSize < chunksize) {
         
         Char_t *newcoord = new Char_t[chunksize];
         memcpy(newcoord, fCoordinates, fCoordinatesSize);
         delete [] fCoordinates;
         fCoordinates = newcoord;
      }
      fCoordinateAllocationSize = chunksize;
   }
   memcpy(fCoordinates + idx * fSingleCoordinateSize, coordbuf, fSingleCoordinateSize);
   fCoordinatesSize += fSingleCoordinateSize;
}
void THnSparseArrayChunk::Sumw2()
{
   
   if (!fSumw2)
      fSumw2 = new TArrayD(fContent->GetSize());
}
ClassImp(THnSparse);
THnSparse::THnSparse():
   fNdimensions(0), fChunkSize(1024), fFilledBins(0), fEntries(0),
   fTsumw(0), fTsumw2(-1.), fCompactCoord(0), fIntegral(0), fIntegralStatus(kNoInt)
{
   
   fBinContent.SetOwner();
}
THnSparse::THnSparse(const char* name, const char* title, Int_t dim,
                     const Int_t* nbins, const Double_t* xmin, const Double_t* xmax,
                     Int_t chunksize):
   TNamed(name, title), fNdimensions(dim), fChunkSize(chunksize), fFilledBins(0),
   fAxes(dim), fEntries(0), fTsumw(0), fTsumw2(-1.), fTsumwx(dim), fTsumwx2(dim),
   fCompactCoord(0), fIntegral(0), fIntegralStatus(kNoInt)
{
   
   
   
   
   
   
   for (Int_t i = 0; i < fNdimensions; ++i) {
      TAxis* axis = new TAxis(nbins[i], xmin ? xmin[i] : 0., xmax ? xmax[i] : 1.);
      TString aname("axis");
      aname += i;
      axis->SetName(aname);
      axis->SetTitle(aname);
      fAxes.AddAtAndExpand(axis, i);
   }
   fAxes.SetOwner();
   fCompactCoord = new THnSparseCompactBinCoord(dim, nbins);
   fBinContent.SetOwner();
}
THnSparse::~THnSparse() {
   
   delete fCompactCoord;
   if (fIntegralStatus != kNoInt) delete [] fIntegral;
}
void THnSparse::AddBinContent(const Int_t* coord, Double_t v)
{
   
   GetCompactCoord()->SetCoord(coord);
   Long_t bin = GetBinIndexForCurrentBin(kTRUE);
   THnSparseArrayChunk* chunk = GetChunk(bin / fChunkSize);
   bin %= fChunkSize;
   v += chunk->fContent->GetAt(bin);
   return chunk->fContent->SetAt(v, bin);
}
THnSparseArrayChunk* THnSparse::AddChunk()
{
   
   THnSparseArrayChunk* first = 0;
   if (fBinContent.GetEntriesFast() > 0)
      first = GetChunk(0);
   THnSparseArrayChunk* chunk =
      new THnSparseArrayChunk(GetCompactCoord()->GetSize(),
                              GetCalculateErrors(), GenerateArray());
   fBinContent.AddLast(chunk);
   return chunk;
}
THnSparse* THnSparse::CloneEmpty(const char* name, const char* title,
                                 const TObjArray* axes, Int_t chunksize) const
{
   
   
   THnSparse* ret = (THnSparse*)IsA()->New();
   ret->SetNameTitle(name, title);
   ret->fNdimensions = axes->GetEntriesFast();
   ret->fChunkSize = chunksize;
   TIter iAxis(axes);
   const TAxis* axis = 0;
   Int_t pos = 0;
   Int_t *nbins = new Int_t[axes->GetEntriesFast()];
   while ((axis = (TAxis*)iAxis())) {
      nbins[pos] = axis->GetNbins();
      ret->fAxes.AddAtAndExpand(axis->Clone(), pos++);
   }
   ret->fAxes.SetOwner();
   ret->fCompactCoord = new THnSparseCompactBinCoord(pos, nbins);
   delete [] nbins;
   return ret;
}
TH1* THnSparse::CreateHist(const char* name, const char* title,
                           const TObjArray* axes) const {
   
   
   
   const int ndim = axes->GetSize();
   TH1* hist = 0;
   
   if (ndim == 1)
      hist = new TH1D(name, title, 1, 0., 1.);
   else if (ndim == 2)
      hist = new TH2D(name, title, 1, 0., 1., 1, 0., 1.);
   else if (ndim == 3)
      hist = new TH3D(name, title, 1, 0., 1., 1, 0., 1., 1, 0., 1.);
   else {
      Error("CreateHist", "Cannot create histogram %s with %d dimensions!", name, ndim);
      return 0;
   }
   TAxis* hax[3] = {hist->GetXaxis(), hist->GetYaxis(), hist->GetZaxis()};
   for (Int_t d = 0; d < ndim; ++d) {
      TAxis* reqaxis = (TAxis*)(*axes)[d];
      if (reqaxis->GetXbins()->GetSize())
         
         hax[d]->Set(reqaxis->GetNbins(), reqaxis->GetXbins()->GetArray());
      else
         
         hax[d]->Set(reqaxis->GetNbins(), reqaxis->GetXmin(), reqaxis->GetXmax());
   }
   hist->Rebuild();
   return hist;
}
Long_t THnSparse::GetBin(const Double_t* x, Bool_t allocate )
{
   
   
   Int_t *coord = GetCompactCoord()->GetCoord();
   for (Int_t i = 0; i < fNdimensions; ++i)
      coord[i] = GetAxis(i)->FindBin(x[i]);
   return GetBinIndexForCurrentBin(allocate);
}
Long_t THnSparse::GetBin(const char* name[], Bool_t allocate )
{
   
   
   Int_t *coord = GetCompactCoord()->GetCoord();
   for (Int_t i = 0; i < fNdimensions; ++i)
      coord[i] = GetAxis(i)->FindBin(name[i]);
   return GetBinIndexForCurrentBin(allocate);
}
Long_t THnSparse::GetBin(const Int_t* coord, Bool_t allocate )
{
   
   
   GetCompactCoord()->SetCoord(coord);
   return GetBinIndexForCurrentBin(allocate);
}
Double_t THnSparse::GetBinContent(const Int_t *coord) const {
   
   GetCompactCoord()->SetCoord(coord);
   Long_t idx = const_cast<THnSparse*>(this)->GetBinIndexForCurrentBin(kFALSE);
   if (idx < 0) return 0.;
   THnSparseArrayChunk* chunk = GetChunk(idx / fChunkSize);
   return chunk->fContent->GetAt(idx % fChunkSize);
}
Double_t THnSparse::GetBinContent(Long64_t idx, Int_t* coord ) const
{
   
   
   
   if (idx >= 0) {
      THnSparseArrayChunk* chunk = GetChunk(idx / fChunkSize);
      idx %= fChunkSize;
      if (chunk && chunk->fContent->GetSize() > idx) {
         if (coord) {
            Int_t sizeCompact = GetCompactCoord()->GetSize();
            memcpy(GetCompactCoord()->GetBuffer(), chunk->fCoordinates + idx * sizeCompact, sizeCompact);
            GetCompactCoord()->GetCoordFromBuffer(coord);
         }
         return chunk->fContent->GetAt(idx);
      }
   }
   if (coord)
      memset(coord, -1, sizeof(Int_t) * fNdimensions);
   return 0.;
}
Double_t THnSparse::GetBinError(const Int_t *coord) const {
   
   
   // END_LATEX
   
   
   if (!GetCalculateErrors())
      return TMath::Sqrt(GetBinContent(coord));
   GetCompactCoord()->SetCoord(coord);
   Long_t idx = const_cast<THnSparse*>(this)->GetBinIndexForCurrentBin(kFALSE);
   if (idx < 0) return 0.;
   THnSparseArrayChunk* chunk = GetChunk(idx / fChunkSize);
   return TMath::Sqrt(chunk->fSumw2->GetAt(idx % fChunkSize));
}
Double_t THnSparse::GetBinError(Long64_t linidx) const {
   
   
   // END_LATEX
   
   
   if (!GetCalculateErrors())
      return TMath::Sqrt(GetBinContent(linidx));
   if (linidx < 0) return 0.;
   THnSparseArrayChunk* chunk = GetChunk(linidx / fChunkSize);
   linidx %= fChunkSize;
   if (!chunk || chunk->fContent->GetSize() < linidx)
      return 0.;
   return TMath::Sqrt(chunk->fSumw2->GetAt(linidx));
}
Long_t THnSparse::GetBinIndexForCurrentBin(Bool_t allocate)
{
   
   
   Long_t hash = GetCompactCoord()->GetHash();
   Long_t linidx = (Long_t) fBins.GetValue(hash);
   while (linidx) {
      
      THnSparseArrayChunk* chunk = GetChunk((linidx - 1)/ fChunkSize);
      if (chunk->Matches((linidx - 1) % fChunkSize, GetCompactCoord()->GetBuffer()))
         return linidx - 1; 
      Long_t nextlinidx = fBinsContinued.GetValue(linidx);
      if (!nextlinidx) break;
      linidx = nextlinidx;
   }
   if (!allocate) return -1;
   ++fFilledBins;
   
   THnSparseArrayChunk *chunk = (THnSparseArrayChunk*) fBinContent.Last();
   Long_t newidx = chunk ? ((Long_t) chunk->GetEntries()) : -1;
   if (!chunk || newidx == (Long_t)fChunkSize) {
      chunk = AddChunk();
      newidx = 0;
   }
   chunk->AddBin(newidx, GetCompactCoord()->GetBuffer());
   
   newidx += (fBinContent.GetEntriesFast() - 1) * fChunkSize;
   if (!linidx)
      
      fBins.Add(hash, newidx + 1);
   else
      
      
      fBinsContinued.Add(linidx, newidx + 1);
   return newidx;
}
THnSparseCompactBinCoord* THnSparse::GetCompactCoord() const
{
   
   if (!fCompactCoord) {
      Int_t *bins = new Int_t[fNdimensions];
      for (Int_t d = 0; d < fNdimensions; ++d)
         bins[d] = GetAxis(d)->GetNbins();
      const_cast<THnSparse*>(this)->fCompactCoord
         = new THnSparseCompactBinCoord(fNdimensions, bins);
      delete [] bins;
   }
   return fCompactCoord;
}
void THnSparse::GetRandom(Double_t *rand, Bool_t subBinRandom )
{
   
   
   
   
   
   if (fIntegralStatus != kValidInt)
      ComputeIntegral();
   
   Double_t p = gRandom->Rndm();
   Long64_t idx = TMath::BinarySearch(GetNbins() + 1, fIntegral, p);
   Int_t bin[20]; 
   GetBinContent(idx, bin);
   
   for (Int_t i = 0; i < fNdimensions; i++) {
      rand[i] = GetAxis(i)->GetBinCenter(bin[i]);
      
      if (subBinRandom)
         rand[i] += (gRandom->Rndm() - 0.5) * GetAxis(i)->GetBinWidth(bin[i]);
   }
   return;
}
Double_t THnSparse::GetSparseFractionBins() const {
   
   Double_t nbinsTotal = 1.;
   for (Int_t d = 0; d < fNdimensions; ++d)
      nbinsTotal *= GetAxis(d)->GetNbins() + 2;
   return fFilledBins / nbinsTotal;
}
Double_t THnSparse::GetSparseFractionMem() const {
   
   
   Int_t arrayElementSize = 0;
   Double_t size = 0.;
   if (fFilledBins) {
      TClass* clArray = GetChunk(0)->fContent->IsA();
      TDataMember* dm = clArray ? clArray->GetDataMember("fArray") : 0;
      arrayElementSize = dm ? dm->GetDataType()->Size() : 0;
   }
   if (!arrayElementSize) {
      Warning("GetSparseFractionMem", "Cannot determine type of elements!");
      return -1.;
   }
   size += fFilledBins * (GetCompactCoord()->GetSize() + arrayElementSize + 2 * sizeof(Long_t) );
   if (fFilledBins && GetChunk(0)->fSumw2)
      size += fFilledBins * sizeof(Float_t); 
   Double_t nbinsTotal = 1.;
   for (Int_t d = 0; d < fNdimensions; ++d)
      nbinsTotal *= GetAxis(d)->GetNbins() + 2;
   return size / nbinsTotal / arrayElementSize;
}
Bool_t THnSparse::IsInRange(Int_t *coord) const
{
   
   
   
   
   
   
   
   
   Int_t min = 0;
   Int_t max = 0;
   for (Int_t i = 0; i < fNdimensions; ++i) {
      TAxis *axis = GetAxis(i);
      if (!axis->TestBit(TAxis::kAxisRange)) continue;
      min = axis->GetFirst();
      max = axis->GetLast();
      if (min == 0 && max == 0) {
         
         
         
         min = 1;
         max = axis->GetNbins();
      }
      if (coord[i] < min || coord[i] > max)
         return kFALSE;
   }
   return kTRUE;
}
TH1D* THnSparse::Projection(Int_t xDim, Option_t* option ) const
{
   
   
   
   return (TH1D*) ProjectionAny(1, &xDim, false, option);
}
TH2D* THnSparse::Projection(Int_t xDim, Int_t yDim, Option_t* option ) const
{
   
   
   
   
   const Int_t dim[2] = {yDim, xDim};
   return (TH2D*) ProjectionAny(2, dim, false, option);
}
TH3D* THnSparse::Projection(Int_t xDim, Int_t yDim, Int_t zDim,
                            Option_t* option ) const
{
   
   
   
   const Int_t dim[3] = {xDim, yDim, zDim};
   return (TH3D*) ProjectionAny(3, dim, false, option);
}
THnSparse* THnSparse::Projection(Int_t ndim, const Int_t* dim,
                                 Option_t* option ) const
{
   
   
   
   return (THnSparse*) ProjectionAny(ndim, dim, true, option);
}
TObject* THnSparse::ProjectionAny(Int_t ndim, const Int_t* dim,
                                  Bool_t wantSparse, Option_t* option ) const
{
   
   
   
   TString name(GetName());
   name += "_";
   for (Int_t d = 0; d < ndim; ++d)
      name += GetAxis(dim[d])->GetName();
   TString title(GetTitle());
   Ssiz_t posInsert = title.First(';');
   if (posInsert == kNPOS) {
      title += " projection ";
      for (Int_t d = 0; d < ndim; ++d)
         title += GetAxis(dim[d])->GetTitle();
   } else {
      for (Int_t d = ndim - 1; d >= 0; --d) {
         title.Insert(posInsert, GetAxis(d)->GetTitle());
         if (d)
            title.Insert(posInsert, ", ");
      }
      title.Insert(posInsert, " projection ");
   }
   TObjArray newaxes(ndim);
   for (Int_t d = 0; d < ndim; ++d) {
      newaxes.AddAt(GetAxis(dim[d]),d);
   }
   THnSparse* sparse = 0;
   TH1* hist = 0;
   TObject* ret = 0;
   if (wantSparse)
      ret = sparse = CloneEmpty(name, title, &newaxes, fChunkSize);
   else
      ret = hist = CreateHist(name, title, &newaxes); 
   Bool_t* hadRange  = new Bool_t[ndim];
   for (Int_t d = 0; d < ndim; ++d){
      TAxis *axis = GetAxis(dim[d]);
      hadRange[d] = axis->TestBit(TAxis::kAxisRange);
      axis->SetBit(TAxis::kAxisRange, kFALSE);
   }
   Bool_t haveErrors = GetCalculateErrors();
   Bool_t wantErrors = (option && (strchr(option, 'E') || strchr(option, 'e'))) || haveErrors;
   Int_t* bins  = new Int_t[ndim];
   Int_t  linbin = 0;
   Int_t* coord = new Int_t[fNdimensions];
   memset(coord, 0, sizeof(Int_t) * fNdimensions);
   Double_t err = 0.;
   Double_t preverr = 0.;
   Double_t v = 0.;
   for (Long64_t i = 0; i < GetNbins(); ++i) {
      v = GetBinContent(i, coord);
      if (!IsInRange(coord)) continue;
      for (Int_t d = 0; d < ndim; ++d)
         bins[d] = coord[dim[d]];
      if (!wantSparse) {
         if (ndim == 1) linbin = bins[0];
         else if (ndim == 2) linbin = hist->GetBin(bins[0], bins[1]);
         else if (ndim == 3) linbin = hist->GetBin(bins[0], bins[1], bins[2]);
      }
      if (wantErrors) {
         if (haveErrors) {
            err = GetBinError(i);
            err *= err;
         } else err = v;
         if (wantSparse) {
            preverr = sparse->GetBinError(bins);
            sparse->SetBinError(bins, TMath::Sqrt(preverr * preverr + err));
         } else {
            preverr = hist->GetBinError(linbin);
            hist->SetBinError(linbin, TMath::Sqrt(preverr * preverr + err));
         }
      }
      
      if (wantSparse)
         sparse->AddBinContent(bins, v);
      else
         hist->AddBinContent(linbin, v);
   }
   delete [] bins;
   delete [] coord;
   if (wantSparse)
      sparse->SetEntries(fEntries);
   else
      hist->SetEntries(fEntries);
   
   for (Int_t d = 0; d < ndim; ++d)
      GetAxis(dim[d])->SetBit(TAxis::kAxisRange, hadRange[d]);
   delete [] hadRange;
   return ret;
}
void THnSparse::Scale(Double_t c)
{
   
   
   
   Int_t* coord = new Int_t[fNdimensions];
   memset(coord, 0, sizeof(Int_t) * fNdimensions);
   
   Bool_t haveErrors = GetCalculateErrors();
   for (Long64_t i = 0; i < GetNbins(); ++i) {
      
      Double_t v = GetBinContent(i, coord);
      SetBinContent(coord, c * v);
      if (haveErrors) {
         Double_t err = GetBinError(coord);
         SetBinError(coord, c * err);
      }
   }
   delete [] coord;
}
void THnSparse::Add(const THnSparse* h, Double_t c)
{
   
   
   
   
   
   if (!CheckConsistency(h, "Add")) return;
   RebinnedAdd(h, c);
}
void THnSparse::RebinnedAdd(const THnSparse* h, Double_t c)
{
   
   
   
   
   
   
   if (fNdimensions!=h->GetNdimensions()) {
      Warning("RebinnedAdd", "Different number of dimensions, cannot carry out operation on the histograms");
      return;
   }
   
   if (!GetCalculateErrors() && h->GetCalculateErrors())
      Sumw2();
   Bool_t haveErrors = GetCalculateErrors();
   
   Int_t* coord = new Int_t[fNdimensions];
   memset(coord, 0, sizeof(Int_t) * fNdimensions);
   
   for (Long64_t i = 0; i < h->GetNbins(); ++i) {
      
      Double_t v = h->GetBinContent(i, coord);
      AddBinContent(coord, c * v);
      if (haveErrors) {
         Double_t err1 = GetBinError(coord);
         Double_t err2 = h->GetBinError(coord) * c;
         SetBinError(coord, TMath::Sqrt(err1 * err1 + err2 * err2));
      }
   }
   delete [] coord;
   Double_t nEntries = GetEntries() + c * h->GetEntries();
   SetEntries(nEntries);
}
Long64_t THnSparse::Merge(TCollection* list)
{
   
   
   if (!list) return 0;
   if (list->IsEmpty()) return (Long64_t)GetEntries();
   TIter iter(list);
   const TObject* addMeObj = 0;
   while ((addMeObj = iter())) {
      const THnSparse* addMe = dynamic_cast<const THnSparse*>(addMeObj);
      if (!addMe) 
         Error("Merge", "Object named %s is not THnSpase! Skipping it.",
               addMeObj->GetName());
      else
         Add(addMe);
   }
   return (Long64_t)GetEntries();
}
void THnSparse::Multiply(const THnSparse* h)
{
   
   
   
   
   
   if(!CheckConsistency(h, "Multiply"))return;
   
   Bool_t wantErrors = kFALSE;
   if (!GetCalculateErrors() && h->GetCalculateErrors())
      wantErrors = kTRUE;
   
   TObjArray newaxes(fNdimensions);
   for (Int_t d = 0; d < fNdimensions; ++d) {
      newaxes.AddAt(GetAxis(d),d);
   }
   if (wantErrors) Sumw2();
   
   Int_t* coord = new Int_t[fNdimensions];
   memset(coord, 0, sizeof(Int_t) * fNdimensions);
   for (Long64_t i = 0; i < GetNbins(); ++i) {
      
      Double_t v1 = GetBinContent(i, coord);
      
      Double_t v2 = h->GetBinContent(coord);
      SetBinContent(coord, v1 * v2);;
      if (wantErrors) {
         Double_t err1 = GetBinError(coord) * v2;
         Double_t err2 = h->GetBinError(coord) * v1;
         SetBinError(coord,TMath::Sqrt((err2 * err2 + err1 * err1)));
      }
   }
   
   delete [] coord;
}
void THnSparse::Divide(const THnSparse *h)
{
   
   
   
   
   
   
   if (!CheckConsistency(h, "Divide"))return;
   
   Bool_t wantErrors=kFALSE;
   if (!GetCalculateErrors() && h->GetCalculateErrors())
      wantErrors=kTRUE;
   
   Double_t nEntries = fEntries;
   
   TObjArray newaxes(fNdimensions);
   for (Int_t d = 0; d < fNdimensions; ++d) {
      newaxes.AddAt(GetAxis(d),d);
   }
   if (wantErrors) Sumw2();
   Bool_t didWarn = kFALSE;
   
   Int_t* coord = new Int_t[fNdimensions];
   memset(coord, 0, sizeof(Int_t) * fNdimensions);
   Double_t err = 0.;
   Double_t b22 = 0.;
   for (Long64_t i = 0; i < GetNbins(); ++i) {
      
      Double_t v1 = GetBinContent(i, coord);
      
      Double_t v2 = h->GetBinContent(coord);
      if (!v2) {
         v1 = 0.;
         v2 = 1.;
         if (!didWarn) {
            Warning("Divide(h)", "Histogram h has empty bins - division by zero! Setting bin to 0.");
            didWarn = kTRUE;
         }
      }
      SetBinContent(coord, v1 / v2);
      if (wantErrors) {
         Double_t err1 = GetBinError(coord) * v2;
         Double_t err2 = h->GetBinError(coord) * v1;
         b22 = v2 * v2;
         err = (err1 * err1 + err2 * err2) / (b22 * b22);
         SetBinError(coord, TMath::Sqrt(err));
      }
   }
   delete [] coord;
   SetEntries(nEntries);
}
void THnSparse::Divide(const THnSparse *h1, const THnSparse *h2, Double_t c1, Double_t c2, Option_t *option)
{
   
   
   
   
   
   
   
   TString opt = option;
   opt.ToLower();
   Bool_t binomial = kFALSE;
   if (opt.Contains("b")) binomial = kTRUE;
   
   if (!CheckConsistency(h1, "Divide") || !CheckConsistency(h2, "Divide"))return;
   if (!c2) {
      Error("Divide","Coefficient of dividing histogram cannot be zero");
      return;
   }
   Reset();
   
   if (!GetCalculateErrors() && (h1->GetCalculateErrors()|| h2->GetCalculateErrors() != 0))
      Sumw2();
   
   Long64_t nFilledBins=0;
   
   Int_t* coord = new Int_t[fNdimensions];
   memset(coord, 0, sizeof(Int_t) * fNdimensions);
   Float_t w = 0.;
   Float_t err = 0.;
   Float_t b22 = 0.;
   Bool_t didWarn = kFALSE;
   for (Long64_t i = 0; i < h1->GetNbins(); ++i) {
      
      Double_t v1 = h1->GetBinContent(i, coord);
      
      Double_t v2 = h2->GetBinContent(coord);
      if (!v2) {
         v1 = 0.;
         v2 = 1.;
         if (!didWarn) {
            Warning("Divide(h1, h2)", "Histogram h2 has empty bins - division by zero! Setting bin to 0.");
            didWarn = kTRUE;
         }
      }
      nFilledBins++;
      SetBinContent(coord, c1 * v1 / c2 / v2);
      if(GetCalculateErrors()){
         Double_t err1=h1->GetBinError(coord);
         Double_t err2=h2->GetBinError(coord);
         if (binomial) {
            if (v1 != v2) {
               w = v1 / v2;
               err2 *= w;
               err = TMath::Abs( ( (1. - 2.*w) * err1 * err1 + err2 * err2 ) / (v2 * v2) );
            } else {
               err = 0;
            }
         } else {
            c1 *= c1;
            c2 *= c2;
            b22 = v2 * v2 * c2;
            err1 *= v2;
            err2 *= v1;
            err = c1 * c2 * (err1 * err1 + err2 * err2) / (b22 * b22);
         }
         SetBinError(coord,TMath::Sqrt(err));
      }
   }
   delete [] coord;
   fFilledBins = nFilledBins;
   
   SetEntries(h1->GetEntries());
}
Bool_t THnSparse::CheckConsistency(const THnSparse *h, const char *tag) const
{
   
   if (fNdimensions!=h->GetNdimensions()) {
      Warning(tag,"Different number of dimensions, cannot carry out operation on the histograms");
      return kFALSE;
   }
   for (Int_t dim = 0; dim < fNdimensions; dim++){
      if (GetAxis(dim)->GetNbins()!=h->GetAxis(dim)->GetNbins()) {
         Warning(tag,"Different number of bins on axis %i, cannot carry out operation on the histograms", dim);
         return kFALSE;
      }
   }
   return kTRUE;
}
void THnSparse::SetBinEdges(Int_t idim, const Double_t* bins)
{
   
   TAxis* axis = (TAxis*) fAxes[idim];
   axis->Set(axis->GetNbins(), bins);
}
void THnSparse::SetBinContent(const Int_t* coord, Double_t v)
{
   
   GetCompactCoord()->SetCoord(coord);
   Long_t bin = GetBinIndexForCurrentBin(kTRUE);
   THnSparseArrayChunk* chunk = GetChunk(bin / fChunkSize);
   chunk->fContent->SetAt(v, bin % fChunkSize);
}
void THnSparse::SetBinError(const Int_t* coord, Double_t e)
{
   
   GetCompactCoord()->SetCoord(coord);
   Long_t bin = GetBinIndexForCurrentBin(kTRUE);
   THnSparseArrayChunk* chunk = GetChunk(bin / fChunkSize);
   chunk->fSumw2->SetAt(e*e, bin % fChunkSize);
}
void THnSparse::Sumw2()
{
   
   if (GetCalculateErrors()) return;
   fTsumw2 = 0.;
   TIter iChunk(&fBinContent);
   THnSparseArrayChunk* chunk = 0;
   while ((chunk = (THnSparseArrayChunk*) iChunk()))
      chunk->Sumw2();
}
THnSparse* THnSparse::Rebin(Int_t group) const
{
   
   
   
   
   Int_t* ngroup = new Int_t[GetNdimensions()];
   for (Int_t d = 0; d < GetNdimensions(); ++d)
      ngroup[d] = group;
   THnSparse* ret = Rebin(ngroup);
   delete [] ngroup;
   return ret;
}
THnSparse* THnSparse::Rebin(const Int_t* group) const
{
   
   
   
   
   Int_t ndim = GetNdimensions();
   TString name(GetName());
   for (Int_t d = 0; d < ndim; ++d)
      name += Form("_%d", group[d]);
   TString title(GetTitle());
   Ssiz_t posInsert = title.First(';');
   if (posInsert == kNPOS) {
      title += " rebin ";
      for (Int_t d = 0; d < ndim; ++d)
         title += Form("{%d}", group[d]);
   } else {
      for (Int_t d = ndim - 1; d >= 0; --d)
         title.Insert(posInsert, Form("{%d}", group[d]));
      title.Insert(posInsert, " rebin ");
   }
   TObjArray newaxes(ndim);
   newaxes.SetOwner();
   for (Int_t d = 0; d < ndim; ++d) {
      newaxes.AddAt(GetAxis(d)->Clone(),d);
      if (group[d] > 1) {
         TAxis* newaxis = (TAxis*) newaxes.At(d);
         Int_t newbins = (newaxis->GetNbins() + group[d] - 1) / group[d];
         if (newaxis->GetXbins() && newaxis->GetXbins()->GetSize()) {
            
            Double_t *edges = new Double_t[newbins + 1];
            for (Int_t i = 0; i < newbins + 1; ++i)
               if (group[d] * i <= newaxis->GetNbins())
                  edges[i] = newaxis->GetXbins()->At(group[d] * i);
               else edges[i] = newaxis->GetXmax();
            newaxis->Set(newbins, edges);
         } else {
            newaxis->Set(newbins, newaxis->GetXmin(), newaxis->GetXmax());
         }
      }
   }
   THnSparse* h = CloneEmpty(name.Data(), title.Data(), &newaxes, fChunkSize);
   Bool_t haveErrors = GetCalculateErrors();
   Bool_t wantErrors = haveErrors;
   Int_t* bins  = new Int_t[ndim];
   Int_t* coord = new Int_t[fNdimensions];
   memset(coord, 0, sizeof(Int_t) * fNdimensions);
   Double_t err = 0.;
   Double_t preverr = 0.;
   Double_t v = 0.;
   for (Long64_t i = 0; i < GetNbins(); ++i) {
      v = GetBinContent(i, coord);
      for (Int_t d = 0; d < ndim; ++d)
         bins[d] = (coord[d] - 1) / group[d] + 1;
      if (wantErrors) {
         if (haveErrors) {
            err = GetBinError(i);
            err *= err;
         } else err = v;
         preverr = h->GetBinError(bins);
         h->SetBinError(bins, TMath::Sqrt(preverr * preverr + err));
      }
      
      h->AddBinContent(bins, v);
   }
   delete [] bins;
   delete [] coord;
   h->SetEntries(fEntries);
   return h;
}
void THnSparse::Reset(Option_t * )
{
   
   fFilledBins = 0;
   fEntries = 0.;
   fTsumw = 0.;
   fTsumw2 = -1.;
   fBins.Delete();
   fBinsContinued.Clear();
   fBinContent.Delete();
   if (fIntegralStatus != kNoInt) {
      delete [] fIntegral;
      fIntegralStatus = kNoInt;
   }
}
Double_t THnSparse::ComputeIntegral()
{
   
   
   if (fIntegralStatus != kNoInt) {
      delete [] fIntegral;
      fIntegralStatus = kNoInt;
   }
   
   if (GetNbins() == 0) {
      Error("ComputeIntegral", "The histogram must have at least one bin.");
      return 0.;
   }
   
   fIntegral = new Double_t [GetNbins() + 1];
   fIntegral[0] = 0.;
   
   Int_t* coord = new Int_t[fNdimensions];
   for (Long64_t i = 0; i < GetNbins(); ++i) {
      Double_t v = GetBinContent(i, coord);
      
      bool regularBin = true;
      for (Int_t dim = 0; dim < fNdimensions; dim++)
         if (coord[dim] < 1 || coord[dim] > GetAxis(dim)->GetNbins()) {
            regularBin = false;
            break;
         }
      
      if (!regularBin) v = 0.;
      fIntegral[i + 1] = fIntegral[i] + v;
   }
   delete [] coord;
   
   if (fIntegral[GetNbins()] == 0.) {
      Error("ComputeIntegral", "No hits in regular bins (non over/underflow).");
      delete [] fIntegral;
      return 0.;
   }
   
   for (Long64_t i = 0; i <= GetNbins(); ++i)
      fIntegral[i] = fIntegral[i] / fIntegral[GetNbins()];
   
   fIntegralStatus = kValidInt;
   return fIntegral[GetNbins()];
}
Last change: Wed Jun 25 08:47:25 2008
Last generated: 2008-06-25 08:47
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.