1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
15 //--------------------------------------------------------------------//
17 // AliCFGridSparse Class //
18 // Class to accumulate data on an N-dimensional grid, to be used //
19 // as input to get corrections for Reconstruction & Trigger efficiency//
20 // Based on root THnSparse //
21 // -- Author : S.Arcelli //
22 // Still to be done: //
23 // --Interpolate among bins in a range //
24 //--------------------------------------------------------------------//
27 #include "AliCFGridSparse.h"
28 #include "THnSparse.h"
37 //____________________________________________________________________
38 ClassImp(AliCFGridSparse)
40 //____________________________________________________________________
41 AliCFGridSparse::AliCFGridSparse() :
46 // default constructor
48 //____________________________________________________________________
49 AliCFGridSparse::AliCFGridSparse(const Char_t* name, const Char_t* title) :
50 AliCFFrame(name,title),
54 // default constructor
56 //____________________________________________________________________
57 AliCFGridSparse::AliCFGridSparse(const Char_t* name, const Char_t* title, Int_t nVarIn, const Int_t * nBinIn) :
58 AliCFFrame(name,title),
66 fData = new THnSparseF(name,title,nVarIn,nBinIn);
69 //____________________________________________________________________
70 AliCFGridSparse::~AliCFGridSparse()
75 if (fData) delete fData;
78 //____________________________________________________________________
79 AliCFGridSparse::AliCFGridSparse(const AliCFGridSparse& c) :
87 ((AliCFGridSparse &)c).Copy(*this);
90 //____________________________________________________________________
91 AliCFGridSparse& AliCFGridSparse::operator=(const AliCFGridSparse &c)
96 if (this != &c) c.Copy(*this);
100 //____________________________________________________________________
101 void AliCFGridSparse::SetBinLimits(Int_t ivar, const Double_t *array)
104 // setting the arrays containing the bin limits
106 fData->SetBinEdges(ivar, array);
109 //____________________________________________________________________
110 void AliCFGridSparse::Fill(const Double_t *var, Double_t weight)
114 // given a set of values of the input variable,
115 // with weight (by default w=1)
117 fData->Fill(var,weight);
120 //___________________________________________________________________
121 TH1D *AliCFGridSparse::Project(Int_t ivar) const
124 // Make a 1D projection along variable ivar
127 TH1D *hist=fData->Projection(ivar);
128 hist->SetXTitle(GetVarTitle(ivar));
129 hist->SetTitle(Form("%s: projection on %s",GetTitle(),GetVarTitle(ivar)));
132 //___________________________________________________________________
133 TH2D *AliCFGridSparse::Project(Int_t ivar1, Int_t ivar2) const
136 // Make a 2D projection along variables ivar1 & ivar2
139 TH2D *hist=fData->Projection(ivar2,ivar1); //notice inverted axis (THnSparse uses TH3 2d-projection convention...)
140 hist->SetXTitle(GetVarTitle(ivar1));
141 hist->SetYTitle(GetVarTitle(ivar2));
142 hist->SetTitle(Form("%s: projection on %s-%s",GetTitle(),GetVarTitle(ivar1),GetVarTitle(ivar2)));
145 //___________________________________________________________________
146 TH3D *AliCFGridSparse::Project(Int_t ivar1, Int_t ivar2, Int_t ivar3) const
149 // Make a 3D projection along variables ivar1 & ivar2 & ivar3
152 TH3D *hist=fData->Projection(ivar1,ivar2,ivar3);
153 hist->SetXTitle(GetVarTitle(ivar1));
154 hist->SetYTitle(GetVarTitle(ivar2));
155 hist->SetZTitle(GetVarTitle(ivar3));
156 hist->SetTitle(Form("%s: projection on %s-%s-%s",GetTitle(),GetVarTitle(ivar1),GetVarTitle(ivar2),GetVarTitle(ivar3)));
160 //___________________________________________________________________
161 AliCFGridSparse* AliCFGridSparse::Project(Int_t nVars, const Int_t* vars, const Double_t* varMin, const Double_t* varMax) const
164 // projects the grid on the nVars dimensions defined in vars.
165 // axis ranges can be defined in arrays varMin, varMax
168 // binning for new grid
169 Int_t* bins = new Int_t[nVars];
170 for (Int_t iVar=0; iVar<nVars; iVar++) {
171 bins[iVar] = GetNBins(vars[iVar]);
174 // create new grid sparse
175 AliCFGridSparse* out = new AliCFGridSparse(fName,fTitle,nVars,bins);
177 //set the range in the THnSparse to project
178 THnSparse* clone = ((THnSparse*)fData->Clone());
179 if (varMin && varMax) {
180 for (Int_t iAxis=0; iAxis<GetNVar(); iAxis++) {
181 clone->GetAxis(iAxis)->SetRangeUser(varMin[iAxis],varMax[iAxis]);
184 else AliInfo("Keeping same axis ranges");
186 out->SetGrid(clone->Projection(nVars,vars));
191 //____________________________________________________________________
192 Float_t AliCFGridSparse::GetBinCenter(Int_t ivar, Int_t ibin) const
195 // Returns the center of specified bin for variable axis ivar
198 return (Float_t) fData->GetAxis(ivar)->GetBinCenter(ibin);
201 //____________________________________________________________________
202 Float_t AliCFGridSparse::GetBinSize(Int_t ivar, Int_t ibin) const
205 // Returns the size of specified bin for variable axis ivar
208 return (Float_t) fData->GetAxis(ivar)->GetBinUpEdge(ibin) - fData->GetAxis(ivar)->GetBinLowEdge(ibin);
211 //____________________________________________________________________
212 Float_t AliCFGridSparse::GetEntries() const
215 // total entries (including overflows and underflows)
218 return fData->GetEntries();
221 //____________________________________________________________________
222 Float_t AliCFGridSparse::GetElement(Long_t index) const
225 // Returns content of grid element index
228 return fData->GetBinContent(index);
230 //____________________________________________________________________
231 Float_t AliCFGridSparse::GetElement(const Int_t *bin) const
234 // Get the content in a bin corresponding to a set of bin indexes
236 return fData->GetBinContent(bin);
239 //____________________________________________________________________
240 Float_t AliCFGridSparse::GetElement(const Double_t *var) const
243 // Get the content in a bin corresponding to a set of input variables
246 Long_t index = fData->GetBin(var,kFALSE);
247 if (index<0) return 0.;
248 return fData->GetBinContent(index);
251 //____________________________________________________________________
252 Float_t AliCFGridSparse::GetElementError(Long_t index) const
255 // Returns the error on the content
258 return fData->GetBinContent(index);
260 //____________________________________________________________________
261 Float_t AliCFGridSparse::GetElementError(const Int_t *bin) const
264 // Get the error in a bin corresponding to a set of bin indexes
266 return fData->GetBinError(bin);
269 //____________________________________________________________________
270 Float_t AliCFGridSparse::GetElementError(const Double_t *var) const
273 // Get the error in a bin corresponding to a set of input variables
276 Long_t index=fData->GetBin(var,kFALSE); //this is the THnSparse index (do not allocate new cells if content is empy)
277 if (index<0) return 0.;
278 return fData->GetBinError(index);
281 //____________________________________________________________________
282 void AliCFGridSparse::SetElement(Long_t index, Float_t val)
285 // Sets grid element value
287 Int_t* bin = new Int_t[GetNVar()];
288 fData->GetBinContent(index,bin); //affects the bin coordinates
293 //____________________________________________________________________
294 void AliCFGridSparse::SetElement(const Int_t *bin, Float_t val)
297 // Sets grid element of bin indeces bin to val
299 fData->SetBinContent(bin,val);
301 //____________________________________________________________________
302 void AliCFGridSparse::SetElement(const Double_t *var, Float_t val)
305 // Set the content in a bin to value val corresponding to a set of input variables
307 Long_t index=fData->GetBin(var,kTRUE); //THnSparse index: allocate the cell
308 Int_t *bin = new Int_t[GetNVar()];
309 fData->GetBinContent(index,bin); //trick to access the array of bins
314 //____________________________________________________________________
315 void AliCFGridSparse::SetElementError(Long_t index, Float_t val)
318 // Sets grid element iel error to val (linear indexing) in AliCFFrame
320 Int_t *bin = new Int_t[GetNVar()];
321 fData->GetBinContent(index,bin);
322 SetElementError(bin,val);
326 //____________________________________________________________________
327 void AliCFGridSparse::SetElementError(const Int_t *bin, Float_t val)
330 // Sets grid element error of bin indeces bin to val
332 fData->SetBinError(bin,val);
334 //____________________________________________________________________
335 void AliCFGridSparse::SetElementError(const Double_t *var, Float_t val)
338 // Set the error in a bin to value val corresponding to a set of input variables
340 Long_t index=fData->GetBin(var); //THnSparse index
341 Int_t *bin = new Int_t[GetNVar()];
342 fData->GetBinContent(index,bin); //trick to access the array of bins
343 SetElementError(bin,val);
347 //____________________________________________________________________
348 void AliCFGridSparse::SumW2()
351 //set calculation of the squared sum of the weighted entries
354 fData->CalculateErrors(kTRUE);
359 //____________________________________________________________________
360 void AliCFGridSparse::Add(const AliCFGridSparse* aGrid, Double_t c)
363 //add aGrid to the current one
366 if (aGrid->GetNVar() != GetNVar()){
367 AliError("Different number of variables, cannot add the grids");
371 if (!fSumW2 && aGrid->GetSumW2()) SumW2();
372 fData->Add(aGrid->GetGrid(),c);
375 //____________________________________________________________________
376 void AliCFGridSparse::Add(const AliCFGridSparse* aGrid1, const AliCFGridSparse* aGrid2, Double_t c1,Double_t c2)
379 //Add aGrid1 and aGrid2 and deposit the result into the current one
382 if (GetNVar() != aGrid1->GetNVar() || GetNVar() != aGrid2->GetNVar()) {
383 AliInfo("Different number of variables, cannot add the grids");
387 if (!fSumW2 && (aGrid1->GetSumW2() || aGrid2->GetSumW2())) SumW2();
390 fData->Add(aGrid1->GetGrid(),c1);
391 fData->Add(aGrid2->GetGrid(),c2);
394 //____________________________________________________________________
395 void AliCFGridSparse::Multiply(const AliCFGridSparse* aGrid, Double_t c)
398 // Multiply aGrid to the current one
401 if (aGrid->GetNVar() != GetNVar()) {
402 AliError("Different number of variables, cannot multiply the grids");
406 if(!fSumW2 && aGrid->GetSumW2()) SumW2();
407 THnSparse *h = aGrid->GetGrid();
412 //____________________________________________________________________
413 void AliCFGridSparse::Multiply(const AliCFGridSparse* aGrid1, const AliCFGridSparse* aGrid2, Double_t c1,Double_t c2)
416 //Multiply aGrid1 and aGrid2 and deposit the result into the current one
419 if (GetNVar() != aGrid1->GetNVar() || GetNVar() != aGrid2->GetNVar()) {
420 AliError("Different number of variables, cannot multiply the grids");
424 if(!fSumW2 && (aGrid1->GetSumW2() || aGrid2->GetSumW2())) SumW2();
427 THnSparse *h1 = aGrid1->GetGrid();
428 THnSparse *h2 = aGrid2->GetGrid();
434 //____________________________________________________________________
435 void AliCFGridSparse::Divide(const AliCFGridSparse* aGrid, Double_t c)
438 // Divide aGrid to the current one
441 if (aGrid->GetNVar() != GetNVar()) {
442 AliError("Different number of variables, cannot divide the grids");
446 if (!fSumW2 && aGrid->GetSumW2()) SumW2();
448 THnSparse *h = aGrid->GetGrid();
453 //____________________________________________________________________
454 void AliCFGridSparse::Divide(const AliCFGridSparse* aGrid1, const AliCFGridSparse* aGrid2, Double_t c1,Double_t c2, Option_t *option)
457 //Divide aGrid1 and aGrid2 and deposit the result into the current one
458 //binomial errors are supported
461 if (GetNVar() != aGrid1->GetNVar() || GetNVar() != aGrid2->GetNVar()) {
462 AliError("Different number of variables, cannot divide the grids");
466 if (!fSumW2 && (aGrid1->GetSumW2() || aGrid2->GetSumW2())) SumW2();
468 THnSparse *h1= aGrid1->GetGrid();
469 THnSparse *h2= aGrid2->GetGrid();
470 fData->Divide(h1,h2,c1,c2,option);
474 //____________________________________________________________________
475 void AliCFGridSparse::Rebin(const Int_t* group)
478 // rebin the grid according to Rebin() as in THnSparse
479 // Please notice that the original number of bins on
480 // a given axis has to be divisible by the rebin group.
483 for(Int_t i=0;i<GetNVar();i++){
484 if (group[i]!=1) AliInfo(Form(" merging bins along dimension %i in groups of %i bins", i,group[i]));
487 THnSparse *rebinned =fData->Rebin(group);
491 //____________________________________________________________________
492 void AliCFGridSparse::Scale(Long_t index, const Double_t *fact)
495 //scale content of a certain cell by (positive) fact (with error)
498 if (GetElement(index)==0 || fact[0]==0) return;
500 Double_t in[2], out[2];
501 in[0]=GetElement(index);
502 in[1]=GetElementError(index);
503 GetScaledValues(fact,in,out);
504 SetElement(index,out[0]);
505 if (fSumW2) SetElementError(index,out[1]);
507 //____________________________________________________________________
508 void AliCFGridSparse::Scale(const Int_t *bin, const Double_t *fact)
511 //scale content of a certain cell by (positive) fact (with error)
513 if(GetElement(bin)==0 || fact[0]==0)return;
515 Double_t in[2], out[2];
516 in[0]=GetElement(bin);
517 in[1]=GetElementError(bin);
518 GetScaledValues(fact,in,out);
519 SetElement(bin,out[0]);
520 if(fSumW2)SetElementError(bin,out[1]);
523 //____________________________________________________________________
524 void AliCFGridSparse::Scale(const Double_t *var, const Double_t *fact)
527 //scale content of a certain cell by (positive) fact (with error)
529 if(GetElement(var)==0 || fact[0]==0)return;
531 Double_t in[2], out[2];
532 in[0]=GetElement(var);
533 in[1]=GetElementError(var);
534 GetScaledValues(fact,in,out);
535 SetElement(var,out[0]);
536 if(fSumW2)SetElementError(var,out[1]);
539 //____________________________________________________________________
540 void AliCFGridSparse::Scale(const Double_t* fact)
543 //scale contents of the whole grid by fact
546 for (Long_t iel=0; iel<GetNFilledBins(); iel++) {
550 //____________________________________________________________________
551 Long_t AliCFGridSparse::GetEmptyBins() const {
556 return (GetNBinsTotal() - GetNFilledBins()) ;
559 //_____________________________________________________________________
560 // Int_t AliCFGridSparse::GetEmptyBins( Double_t *varMin, Double_t* varMax ) const
563 // // Get empty bins in a range specified by varMin and varMax
566 // Int_t *indexMin = new Int_t[GetNVar()];
567 // Int_t *indexMax = new Int_t[GetNVar()];
569 // //Find out the min and max bins
571 // for (Int_t i=0;i<GetNVar();i++) {
572 // Double_t xmin=varMin[i]; // the min values
573 // Double_t xmax=varMax[i]; // the min values
574 // Int_t nbins = fData->GetAxis(i)->GetNbins()+1;
575 // Double_t *bins=new Double_t[nbins];
576 // for(Int_t ibin =0;ibin<nbins;ibin++){
577 // bins[ibin] = fVarBinLimits[ibin+fOffset[i]];
579 // indexMin[i] = TMath::BinarySearch(nbins,bins,xmin);
580 // indexMax[i] = TMath::BinarySearch(nbins,bins,xmax);
585 // for(Int_t i=0;i<fNDim;i++){
586 // for (Int_t j=0;j<GetNVar();j++)fIndex[j]=GetBinIndex(j,i);
587 // Bool_t isIn=kTRUE;
588 // for (Int_t j=0;j<GetNVar();j++){
589 // if(!(fIndex[j]>=indexMin[j] && fIndex[j]<=indexMax[j]))isIn=kFALSE;
591 // if(isIn && GetElement(i)<=0)val++;
593 // AliInfo(Form(" the empty bins = %i ",val));
595 // delete [] indexMin;
596 // delete [] indexMax;
600 //____________________________________________________________________
601 Int_t AliCFGridSparse::CheckStats(Double_t thr) const
604 // Count the cells below a certain threshold
607 for (Int_t i=0; i<GetNBinsTotal(); i++) {
608 if (GetElement(i)<thr) ncellsLow++;
613 //_____________________________________________________________________
614 Double_t AliCFGridSparse::GetIntegral() const
619 return fData->ComputeIntegral();
622 //_____________________________________________________________________
623 // Double_t AliCFGridSparse::GetIntegral(Int_t *binMin, Int_t* binMax ) const
626 // // Get Integral in a range of bin indeces (extremes included)
631 // for(Int_t i=0;i<GetNVar();i++){
632 // if(binMin[i]<1)binMin[i]=1;
633 // if(binMax[i]>fNVarBins[i])binMax[i]=fNVarBins[i];
634 // if((binMin[i]>binMax[i])){
635 // AliInfo(Form(" Bin indices in variable %i in reverse order, please check!", i));
639 // val=GetSum(0,binMin,binMax);
644 //_____________________________________________________________________
645 // Double_t AliCFGridSparse::GetIntegral(Double_t *varMin, Double_t* varMax ) const
648 // // Get Integral in a range (extremes included)
651 // Int_t *indexMin=new Int_t[GetNVar()];
652 // Int_t *indexMax=new Int_t[GetNVar()];
654 // //Find out the min and max bins
656 // for(Int_t i=0;i<GetNVar();i++){
657 // Double_t xmin=varMin[i]; // the min values
658 // Double_t xmax=varMax[i]; // the min values
659 // Int_t nbins=fNVarBins[i]+1;
660 // Double_t *bins=new Double_t[nbins];
661 // for(Int_t ibin =0;ibin<nbins;ibin++){
662 // bins[ibin] = fVarBinLimits[ibin+fOffset[i]];
664 // indexMin[i] = TMath::BinarySearch(nbins,bins,xmin);
665 // indexMax[i] = TMath::BinarySearch(nbins,bins,xmax);
669 // //move to the TH/THnSparse convention in N-dim bin numbering
670 // for(Int_t i=0;i<GetNVar();i++){
675 // Double_t val=GetIntegral(indexMin,indexMax);
677 // delete [] indexMin;
678 // delete [] indexMax;
683 // //_____________________________________________________________________
684 // Double_t AliCFGridSparse::GetSum(Int_t ivar, Int_t *binMin, Int_t* binMax) const
687 // // recursively add over nested loops....
690 // static Double_t val;
691 // if (ivar==0) val=0.;
692 // for(Int_t ibin=binMin[ivar]-1 ; ibin<=binMax[ivar]-1 ; ibin++){
693 // //-1 is to move from TH/ThnSparse N-dim bin convention to one in AliCFFrame
694 // fIndex[ivar]=ibin;
695 // if(ivar<GetNVar()-1) {
696 // val=GetSum(ivar+1,binMin,binMax);
699 // Int_t iel=GetBinIndex(fIndex);
700 // val+=GetElement(iel);
707 //____________________________________________________________________
708 Long64_t AliCFGridSparse::Merge(TCollection* list)
711 // Merge a list of AliCFGridSparse with this (needed for PROOF).
712 // Returns the number of merged objects (including this).
721 TIterator* iter = list->MakeIterator();
725 while ((obj = iter->Next())) {
726 AliCFGridSparse* entry = dynamic_cast<AliCFGridSparse*> (obj);
736 //____________________________________________________________________
737 void AliCFGridSparse::GetScaledValues(const Double_t *fact, const Double_t *in, Double_t *out) const{
739 // scale input *in and its error by (positive) fact (with error)
740 // and erite it to *out
742 out[0]=in[0]*fact[0];
743 out[1]=TMath::Sqrt(in[1]*in[1]/in[0]/in[0]
744 +fact[1]*fact[1]/fact[0]/fact[0])*out[0];
748 //____________________________________________________________________
749 void AliCFGridSparse::Copy(TObject& c) const
755 AliCFGridSparse& target = (AliCFGridSparse &) c;
756 target.fSumW2 = fSumW2 ;
758 target.fData = (THnSparse*)fData->Clone();
762 //____________________________________________________________________
763 TH1D* AliCFGridSparse::Slice(Int_t iVar, const Double_t *varMin, const Double_t *varMax) const
766 // return a slice (1D-projection) on variable iVar while axis ranges are defined with varMin,varMax
767 // arrays varMin and varMax contain the min and max values of each variable.
768 // therefore varMin and varMax must have their dimensions equal to GetNVar()
771 THnSparse* clone = (THnSparse*)fData->Clone();
772 for (Int_t iAxis=0; iAxis<GetNVar(); iAxis++) {
773 clone->GetAxis(iAxis)->SetRangeUser(varMin[iAxis],varMax[iAxis]);
775 return clone->Projection(iVar);
778 //____________________________________________________________________
779 TH2D* AliCFGridSparse::Slice(Int_t iVar1, Int_t iVar2, const Double_t *varMin, const Double_t *varMax) const
782 // return a slice (2D-projection) on variables iVar1 and iVar2 while axis ranges are defined with varMin,varMax
783 // arrays varMin and varMax contain the min and max values of each variable.
784 // therefore varMin and varMax must have their dimensions equal to GetNVar()
787 THnSparse* clone = (THnSparse*)fData->Clone();
788 for (Int_t iAxis=0; iAxis<GetNVar(); iAxis++) {
789 clone->GetAxis(iAxis)->SetRangeUser(varMin[iAxis],varMax[iAxis]);
791 return clone->Projection(iVar1,iVar2);
794 //____________________________________________________________________
795 TH3D* AliCFGridSparse::Slice(Int_t iVar1, Int_t iVar2, Int_t iVar3, const Double_t *varMin, const Double_t *varMax) const
798 // return a slice (3D-projection) on variables iVar1, iVar2 and iVar3 while axis ranges are defined with varMin,varMax
799 // arrays varMin and varMax contain the min and max values of each variable.
800 // therefore varMin and varMax must have their dimensions equal to GetNVar()
803 THnSparse* clone = (THnSparse*)fData->Clone();
804 for (Int_t iAxis=0; iAxis<GetNVar(); iAxis++) {
805 clone->GetAxis(iAxis)->SetRangeUser(varMin[iAxis],varMax[iAxis]);
807 return clone->Projection(iVar1,iVar2,iVar3);
810 //____________________________________________________________________
811 void AliCFGridSparse::SetRangeUser(Int_t iVar, Double_t varMin, Double_t varMax) {
813 // set range of axis iVar.
815 fData->GetAxis(iVar)->SetRangeUser(varMin,varMax);
816 AliWarning(Form("THnSparse axis %d range has been modified",iVar));
819 //____________________________________________________________________
820 void AliCFGridSparse::SetRangeUser(const Double_t *varMin, const Double_t *varMax) {
822 // set range of every axis. varMin and varMax must be of dimension GetNVar()
824 for (Int_t iAxis=0; iAxis<GetNVar() ; iAxis++) { // set new range for every axis
825 SetRangeUser(iAxis,varMin[iAxis],varMax[iAxis]);
827 AliWarning("THnSparse axes ranges have been modified");
830 //____________________________________________________________________
831 void AliCFGridSparse::UseAxisRange(Bool_t b) const {
832 for (Int_t iAxis=0; iAxis<GetNVar(); iAxis++) fData->GetAxis(iAxis)->SetBit(TAxis::kAxisRange,b);
835 //____________________________________________________________________
836 Float_t AliCFGridSparse::GetOverFlows(Int_t ivar, Bool_t exclusive) const
839 // Returns overflows in variable ivar
840 // Set 'exclusive' to true for an exclusive check on variable ivar
842 Int_t* bin = new Int_t[GetNVar()];
843 memset(bin, 0, sizeof(Int_t) * GetNVar());
845 for (Long64_t i = 0; i < fData->GetNbins(); i++) {
846 Double_t v = fData->GetBinContent(i, bin);
849 for(Int_t j=0;j<GetNVar();j++){
851 if((bin[j]==0) || (bin[j]==GetNBins(j)+1))add=kFALSE;
854 if(bin[ivar]==GetNBins(ivar)+1 && add) ovfl+=v;
861 //____________________________________________________________________
862 Float_t AliCFGridSparse::GetUnderFlows(Int_t ivar, Bool_t exclusive) const
865 // Returns exclusive overflows in variable ivar
866 // Set 'exclusive' to true for an exclusive check on variable ivar
868 Int_t* bin = new Int_t[GetNVar()];
869 memset(bin, 0, sizeof(Int_t) * GetNVar());
871 for (Long64_t i = 0; i < fData->GetNbins(); i++) {
872 Double_t v = fData->GetBinContent(i, bin);
875 for(Int_t j=0;j<GetNVar();j++){
877 if((bin[j]==0) || (bin[j]==GetNBins(j)+1))add=kFALSE;
880 if(bin[ivar]==0 && add) unfl+=v;