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, Double_t min, Double_t max)
104 // set a uniform binning for variable ivar
106 Int_t nBins = GetNBins(ivar);
107 Double_t * array = new Double_t[nBins+1];
108 for (Int_t iEdge=0; iEdge<=nBins; iEdge++) array[iEdge] = min + iEdge * (max-min)/nBins ;
109 fData->SetBinEdges(ivar, array);
113 //____________________________________________________________________
114 void AliCFGridSparse::SetBinLimits(Int_t ivar, const Double_t *array)
117 // setting the arrays containing the bin limits
119 fData->SetBinEdges(ivar, array);
122 //____________________________________________________________________
123 void AliCFGridSparse::Fill(const Double_t *var, Double_t weight)
127 // given a set of values of the input variable,
128 // with weight (by default w=1)
130 fData->Fill(var,weight);
133 //___________________________________________________________________
134 TH1D *AliCFGridSparse::Project(Int_t ivar) const
137 // Make a 1D projection along variable ivar
140 TH1D *hist=fData->Projection(ivar);
141 hist->SetXTitle(GetVarTitle(ivar));
142 hist->SetName(Form("%s_proj-%s",GetName(),GetVarTitle(ivar)));
143 hist->SetTitle(Form("%s: projection on %s",GetTitle(),GetVarTitle(ivar)));
146 //___________________________________________________________________
147 TH2D *AliCFGridSparse::Project(Int_t ivar1, Int_t ivar2) const
150 // Make a 2D projection along variables ivar1 & ivar2
153 TH2D *hist=fData->Projection(ivar2,ivar1); //notice inverted axis (THnSparse uses TH3 2d-projection convention...)
154 hist->SetXTitle(GetVarTitle(ivar1));
155 hist->SetYTitle(GetVarTitle(ivar2));
156 hist->SetName(Form("%s_proj-%s,%s",GetName(),GetVarTitle(ivar1),GetVarTitle(ivar2)));
157 hist->SetTitle(Form("%s: projection on %s-%s",GetTitle(),GetVarTitle(ivar1),GetVarTitle(ivar2)));
160 //___________________________________________________________________
161 TH3D *AliCFGridSparse::Project(Int_t ivar1, Int_t ivar2, Int_t ivar3) const
164 // Make a 3D projection along variables ivar1 & ivar2 & ivar3
167 TH3D *hist=fData->Projection(ivar1,ivar2,ivar3);
168 hist->SetXTitle(GetVarTitle(ivar1));
169 hist->SetYTitle(GetVarTitle(ivar2));
170 hist->SetZTitle(GetVarTitle(ivar3));
171 hist->SetName(Form("%s_proj-%s,%s,%s",GetName(),GetVarTitle(ivar1),GetVarTitle(ivar2),GetVarTitle(ivar3)));
172 hist->SetTitle(Form("%s: projection on %s-%s-%s",GetTitle(),GetVarTitle(ivar1),GetVarTitle(ivar2),GetVarTitle(ivar3)));
176 //___________________________________________________________________
177 AliCFGridSparse* AliCFGridSparse::Project(Int_t nVars, const Int_t* vars, const Double_t* varMin, const Double_t* varMax, Bool_t useBins) const
180 // projects the grid on the nVars dimensions defined in vars.
181 // axis ranges can be defined in arrays varMin, varMax
182 // If useBins=true, varMin and varMax are taken as bin numbers
185 // binning for new grid
186 Int_t* bins = new Int_t[nVars];
187 for (Int_t iVar=0; iVar<nVars; iVar++) {
188 bins[iVar] = GetNBins(vars[iVar]);
191 // create new grid sparse
192 AliCFGridSparse* out = new AliCFGridSparse(fName,fTitle,nVars,bins);
194 //set the range in the THnSparse to project
195 THnSparse* clone = ((THnSparse*)fData->Clone());
196 if (varMin && varMax) {
197 for (Int_t iAxis=0; iAxis<GetNVar(); iAxis++) {
198 if (useBins) clone->GetAxis(iAxis)->SetRange((Int_t)varMin[iAxis],(Int_t)varMax[iAxis]);
199 else clone->GetAxis(iAxis)->SetRangeUser(varMin[iAxis],varMax[iAxis]);
202 else AliInfo("Keeping same axis ranges");
204 out->SetGrid(clone->Projection(nVars,vars));
209 //____________________________________________________________________
210 Float_t AliCFGridSparse::GetBinCenter(Int_t ivar, Int_t ibin) const
213 // Returns the center of specified bin for variable axis ivar
216 return (Float_t) fData->GetAxis(ivar)->GetBinCenter(ibin);
219 //____________________________________________________________________
220 Float_t AliCFGridSparse::GetBinSize(Int_t ivar, Int_t ibin) const
223 // Returns the size of specified bin for variable axis ivar
226 return (Float_t) fData->GetAxis(ivar)->GetBinUpEdge(ibin) - fData->GetAxis(ivar)->GetBinLowEdge(ibin);
229 //____________________________________________________________________
230 Float_t AliCFGridSparse::GetEntries() const
233 // total entries (including overflows and underflows)
236 return fData->GetEntries();
239 //____________________________________________________________________
240 Float_t AliCFGridSparse::GetElement(Long_t index) const
243 // Returns content of grid element index
246 return fData->GetBinContent(index);
248 //____________________________________________________________________
249 Float_t AliCFGridSparse::GetElement(const Int_t *bin) const
252 // Get the content in a bin corresponding to a set of bin indexes
254 return fData->GetBinContent(bin);
257 //____________________________________________________________________
258 Float_t AliCFGridSparse::GetElement(const Double_t *var) const
261 // Get the content in a bin corresponding to a set of input variables
264 Long_t index = fData->GetBin(var,kFALSE);
265 if (index<0) return 0.;
266 return fData->GetBinContent(index);
269 //____________________________________________________________________
270 Float_t AliCFGridSparse::GetElementError(Long_t index) const
273 // Returns the error on the content
276 return fData->GetBinContent(index);
278 //____________________________________________________________________
279 Float_t AliCFGridSparse::GetElementError(const Int_t *bin) const
282 // Get the error in a bin corresponding to a set of bin indexes
284 return fData->GetBinError(bin);
287 //____________________________________________________________________
288 Float_t AliCFGridSparse::GetElementError(const Double_t *var) const
291 // Get the error in a bin corresponding to a set of input variables
294 Long_t index=fData->GetBin(var,kFALSE); //this is the THnSparse index (do not allocate new cells if content is empy)
295 if (index<0) return 0.;
296 return fData->GetBinError(index);
299 //____________________________________________________________________
300 void AliCFGridSparse::SetElement(Long_t index, Float_t val)
303 // Sets grid element value
305 Int_t* bin = new Int_t[GetNVar()];
306 fData->GetBinContent(index,bin); //affects the bin coordinates
311 //____________________________________________________________________
312 void AliCFGridSparse::SetElement(const Int_t *bin, Float_t val)
315 // Sets grid element of bin indeces bin to val
317 fData->SetBinContent(bin,val);
319 //____________________________________________________________________
320 void AliCFGridSparse::SetElement(const Double_t *var, Float_t val)
323 // Set the content in a bin to value val corresponding to a set of input variables
325 Long_t index=fData->GetBin(var,kTRUE); //THnSparse index: allocate the cell
326 Int_t *bin = new Int_t[GetNVar()];
327 fData->GetBinContent(index,bin); //trick to access the array of bins
332 //____________________________________________________________________
333 void AliCFGridSparse::SetElementError(Long_t index, Float_t val)
336 // Sets grid element iel error to val (linear indexing) in AliCFFrame
338 Int_t *bin = new Int_t[GetNVar()];
339 fData->GetBinContent(index,bin);
340 SetElementError(bin,val);
344 //____________________________________________________________________
345 void AliCFGridSparse::SetElementError(const Int_t *bin, Float_t val)
348 // Sets grid element error of bin indeces bin to val
350 fData->SetBinError(bin,val);
352 //____________________________________________________________________
353 void AliCFGridSparse::SetElementError(const Double_t *var, Float_t val)
356 // Set the error in a bin to value val corresponding to a set of input variables
358 Long_t index=fData->GetBin(var); //THnSparse index
359 Int_t *bin = new Int_t[GetNVar()];
360 fData->GetBinContent(index,bin); //trick to access the array of bins
361 SetElementError(bin,val);
365 //____________________________________________________________________
366 void AliCFGridSparse::SumW2()
369 //set calculation of the squared sum of the weighted entries
372 fData->CalculateErrors(kTRUE);
377 //____________________________________________________________________
378 void AliCFGridSparse::Add(const AliCFGridSparse* aGrid, Double_t c)
381 //add aGrid to the current one
384 if (aGrid->GetNVar() != GetNVar()){
385 AliError("Different number of variables, cannot add the grids");
389 if (!fSumW2 && aGrid->GetSumW2()) SumW2();
390 fData->Add(aGrid->GetGrid(),c);
393 //____________________________________________________________________
394 void AliCFGridSparse::Add(const AliCFGridSparse* aGrid1, const AliCFGridSparse* aGrid2, Double_t c1,Double_t c2)
397 //Add aGrid1 and aGrid2 and deposit the result into the current one
400 if (GetNVar() != aGrid1->GetNVar() || GetNVar() != aGrid2->GetNVar()) {
401 AliInfo("Different number of variables, cannot add the grids");
405 if (!fSumW2 && (aGrid1->GetSumW2() || aGrid2->GetSumW2())) SumW2();
408 fData->Add(aGrid1->GetGrid(),c1);
409 fData->Add(aGrid2->GetGrid(),c2);
412 //____________________________________________________________________
413 void AliCFGridSparse::Multiply(const AliCFGridSparse* aGrid, Double_t c)
416 // Multiply aGrid to the current one
419 if (aGrid->GetNVar() != GetNVar()) {
420 AliError("Different number of variables, cannot multiply the grids");
424 if(!fSumW2 && aGrid->GetSumW2()) SumW2();
425 THnSparse *h = aGrid->GetGrid();
430 //____________________________________________________________________
431 void AliCFGridSparse::Multiply(const AliCFGridSparse* aGrid1, const AliCFGridSparse* aGrid2, Double_t c1,Double_t c2)
434 //Multiply aGrid1 and aGrid2 and deposit the result into the current one
437 if (GetNVar() != aGrid1->GetNVar() || GetNVar() != aGrid2->GetNVar()) {
438 AliError("Different number of variables, cannot multiply the grids");
442 if(!fSumW2 && (aGrid1->GetSumW2() || aGrid2->GetSumW2())) SumW2();
445 THnSparse *h1 = aGrid1->GetGrid();
446 THnSparse *h2 = aGrid2->GetGrid();
452 //____________________________________________________________________
453 void AliCFGridSparse::Divide(const AliCFGridSparse* aGrid, Double_t c)
456 // Divide aGrid to the current one
459 if (aGrid->GetNVar() != GetNVar()) {
460 AliError("Different number of variables, cannot divide the grids");
464 if (!fSumW2 && aGrid->GetSumW2()) SumW2();
466 THnSparse *h1 = aGrid->GetGrid();
467 THnSparse *h2 = (THnSparse*)fData->Clone();
468 fData->Divide(h2,h1);
472 //____________________________________________________________________
473 void AliCFGridSparse::Divide(const AliCFGridSparse* aGrid1, const AliCFGridSparse* aGrid2, Double_t c1,Double_t c2, Option_t *option)
476 //Divide aGrid1 and aGrid2 and deposit the result into the current one
477 //binomial errors are supported
480 if (GetNVar() != aGrid1->GetNVar() || GetNVar() != aGrid2->GetNVar()) {
481 AliError("Different number of variables, cannot divide the grids");
485 if (!fSumW2 && (aGrid1->GetSumW2() || aGrid2->GetSumW2())) SumW2();
487 THnSparse *h1= aGrid1->GetGrid();
488 THnSparse *h2= aGrid2->GetGrid();
489 fData->Divide(h1,h2,c1,c2,option);
493 //____________________________________________________________________
494 void AliCFGridSparse::Rebin(const Int_t* group)
497 // rebin the grid according to Rebin() as in THnSparse
498 // Please notice that the original number of bins on
499 // a given axis has to be divisible by the rebin group.
502 for(Int_t i=0;i<GetNVar();i++){
503 if (group[i]!=1) AliInfo(Form(" merging bins along dimension %i in groups of %i bins", i,group[i]));
506 THnSparse *rebinned =fData->Rebin(group);
510 //____________________________________________________________________
511 void AliCFGridSparse::Scale(Long_t index, const Double_t *fact)
514 //scale content of a certain cell by (positive) fact (with error)
517 if (GetElement(index)==0 || fact[0]==0) return;
519 Double_t in[2], out[2];
520 in[0]=GetElement(index);
521 in[1]=GetElementError(index);
522 GetScaledValues(fact,in,out);
523 SetElement(index,out[0]);
524 if (fSumW2) SetElementError(index,out[1]);
526 //____________________________________________________________________
527 void AliCFGridSparse::Scale(const Int_t *bin, const Double_t *fact)
530 //scale content of a certain cell by (positive) fact (with error)
532 if(GetElement(bin)==0 || fact[0]==0)return;
534 Double_t in[2], out[2];
535 in[0]=GetElement(bin);
536 in[1]=GetElementError(bin);
537 GetScaledValues(fact,in,out);
538 SetElement(bin,out[0]);
539 if(fSumW2)SetElementError(bin,out[1]);
542 //____________________________________________________________________
543 void AliCFGridSparse::Scale(const Double_t *var, const Double_t *fact)
546 //scale content of a certain cell by (positive) fact (with error)
548 if(GetElement(var)==0 || fact[0]==0)return;
550 Double_t in[2], out[2];
551 in[0]=GetElement(var);
552 in[1]=GetElementError(var);
553 GetScaledValues(fact,in,out);
554 SetElement(var,out[0]);
555 if(fSumW2)SetElementError(var,out[1]);
558 //____________________________________________________________________
559 void AliCFGridSparse::Scale(const Double_t* fact)
562 //scale contents of the whole grid by fact
565 for (Long_t iel=0; iel<GetNFilledBins(); iel++) {
569 //____________________________________________________________________
570 Long_t AliCFGridSparse::GetEmptyBins() const {
575 return (GetNBinsTotal() - GetNFilledBins()) ;
578 //_____________________________________________________________________
579 // Int_t AliCFGridSparse::GetEmptyBins( Double_t *varMin, Double_t* varMax ) const
582 // // Get empty bins in a range specified by varMin and varMax
585 // Int_t *indexMin = new Int_t[GetNVar()];
586 // Int_t *indexMax = new Int_t[GetNVar()];
588 // //Find out the min and max bins
590 // for (Int_t i=0;i<GetNVar();i++) {
591 // Double_t xmin=varMin[i]; // the min values
592 // Double_t xmax=varMax[i]; // the min values
593 // Int_t nbins = fData->GetAxis(i)->GetNbins()+1;
594 // Double_t *bins=new Double_t[nbins];
595 // for(Int_t ibin =0;ibin<nbins;ibin++){
596 // bins[ibin] = fVarBinLimits[ibin+fOffset[i]];
598 // indexMin[i] = TMath::BinarySearch(nbins,bins,xmin);
599 // indexMax[i] = TMath::BinarySearch(nbins,bins,xmax);
604 // for(Int_t i=0;i<fNDim;i++){
605 // for (Int_t j=0;j<GetNVar();j++)fIndex[j]=GetBinIndex(j,i);
606 // Bool_t isIn=kTRUE;
607 // for (Int_t j=0;j<GetNVar();j++){
608 // if(!(fIndex[j]>=indexMin[j] && fIndex[j]<=indexMax[j]))isIn=kFALSE;
610 // if(isIn && GetElement(i)<=0)val++;
612 // AliInfo(Form(" the empty bins = %i ",val));
614 // delete [] indexMin;
615 // delete [] indexMax;
619 //____________________________________________________________________
620 Int_t AliCFGridSparse::CheckStats(Double_t thr) const
623 // Count the cells below a certain threshold
626 for (Int_t i=0; i<GetNBinsTotal(); i++) {
627 if (GetElement(i)<thr) ncellsLow++;
632 //_____________________________________________________________________
633 Double_t AliCFGridSparse::GetIntegral() const
638 return fData->ComputeIntegral();
641 //_____________________________________________________________________
642 // Double_t AliCFGridSparse::GetIntegral(Int_t *binMin, Int_t* binMax ) const
645 // // Get Integral in a range of bin indeces (extremes included)
650 // for(Int_t i=0;i<GetNVar();i++){
651 // if(binMin[i]<1)binMin[i]=1;
652 // if(binMax[i]>fNVarBins[i])binMax[i]=fNVarBins[i];
653 // if((binMin[i]>binMax[i])){
654 // AliInfo(Form(" Bin indices in variable %i in reverse order, please check!", i));
658 // val=GetSum(0,binMin,binMax);
663 //_____________________________________________________________________
664 // Double_t AliCFGridSparse::GetIntegral(Double_t *varMin, Double_t* varMax ) const
667 // // Get Integral in a range (extremes included)
670 // Int_t *indexMin=new Int_t[GetNVar()];
671 // Int_t *indexMax=new Int_t[GetNVar()];
673 // //Find out the min and max bins
675 // for(Int_t i=0;i<GetNVar();i++){
676 // Double_t xmin=varMin[i]; // the min values
677 // Double_t xmax=varMax[i]; // the min values
678 // Int_t nbins=fNVarBins[i]+1;
679 // Double_t *bins=new Double_t[nbins];
680 // for(Int_t ibin =0;ibin<nbins;ibin++){
681 // bins[ibin] = fVarBinLimits[ibin+fOffset[i]];
683 // indexMin[i] = TMath::BinarySearch(nbins,bins,xmin);
684 // indexMax[i] = TMath::BinarySearch(nbins,bins,xmax);
688 // //move to the TH/THnSparse convention in N-dim bin numbering
689 // for(Int_t i=0;i<GetNVar();i++){
694 // Double_t val=GetIntegral(indexMin,indexMax);
696 // delete [] indexMin;
697 // delete [] indexMax;
702 // //_____________________________________________________________________
703 // Double_t AliCFGridSparse::GetSum(Int_t ivar, Int_t *binMin, Int_t* binMax) const
706 // // recursively add over nested loops....
709 // static Double_t val;
710 // if (ivar==0) val=0.;
711 // for(Int_t ibin=binMin[ivar]-1 ; ibin<=binMax[ivar]-1 ; ibin++){
712 // //-1 is to move from TH/ThnSparse N-dim bin convention to one in AliCFFrame
713 // fIndex[ivar]=ibin;
714 // if(ivar<GetNVar()-1) {
715 // val=GetSum(ivar+1,binMin,binMax);
718 // Int_t iel=GetBinIndex(fIndex);
719 // val+=GetElement(iel);
726 //____________________________________________________________________
727 Long64_t AliCFGridSparse::Merge(TCollection* list)
730 // Merge a list of AliCFGridSparse with this (needed for PROOF).
731 // Returns the number of merged objects (including this).
740 TIterator* iter = list->MakeIterator();
744 while ((obj = iter->Next())) {
745 AliCFGridSparse* entry = dynamic_cast<AliCFGridSparse*> (obj);
755 //____________________________________________________________________
756 void AliCFGridSparse::GetScaledValues(const Double_t *fact, const Double_t *in, Double_t *out) const{
758 // scale input *in and its error by (positive) fact (with error)
759 // and erite it to *out
761 out[0]=in[0]*fact[0];
762 out[1]=TMath::Sqrt(in[1]*in[1]/in[0]/in[0]
763 +fact[1]*fact[1]/fact[0]/fact[0])*out[0];
767 //____________________________________________________________________
768 void AliCFGridSparse::Copy(TObject& c) const
774 AliCFGridSparse& target = (AliCFGridSparse &) c;
775 target.fSumW2 = fSumW2 ;
777 target.fData = (THnSparse*)fData->Clone();
781 //____________________________________________________________________
782 TH1D* AliCFGridSparse::Slice(Int_t iVar, const Double_t *varMin, const Double_t *varMax, Bool_t useBins) const
785 // return a slice (1D-projection) on variable iVar while axis ranges are defined with varMin,varMax
786 // arrays varMin and varMax contain the min and max values of each variable.
787 // therefore varMin and varMax must have their dimensions equal to GetNVar()
788 // If useBins=true, varMin and varMax are taken as bin numbers
790 THnSparse* clone = (THnSparse*)fData->Clone();
791 for (Int_t iAxis=0; iAxis<GetNVar(); iAxis++) {
792 if (useBins) clone->GetAxis(iAxis)->SetRange((Int_t)varMin[iAxis],(Int_t)varMax[iAxis]);
793 else clone->GetAxis(iAxis)->SetRangeUser(varMin[iAxis],varMax[iAxis]);
795 return clone->Projection(iVar);
798 //____________________________________________________________________
799 TH2D* AliCFGridSparse::Slice(Int_t iVar1, Int_t iVar2, const Double_t *varMin, const Double_t *varMax, Bool_t useBins) const
802 // return a slice (2D-projection) on variables iVar1 and iVar2 while axis ranges are defined with varMin,varMax
803 // arrays varMin and varMax contain the min and max values of each variable.
804 // therefore varMin and varMax must have their dimensions equal to GetNVar()
805 // If useBins=true, varMin and varMax are taken as bin numbers
807 THnSparse* clone = (THnSparse*)fData->Clone();
808 for (Int_t iAxis=0; iAxis<GetNVar(); iAxis++) {
809 if (useBins) clone->GetAxis(iAxis)->SetRange((Int_t)varMin[iAxis],(Int_t)varMax[iAxis]);
810 else clone->GetAxis(iAxis)->SetRangeUser(varMin[iAxis],varMax[iAxis]);
812 return clone->Projection(iVar1,iVar2);
815 //____________________________________________________________________
816 TH3D* AliCFGridSparse::Slice(Int_t iVar1, Int_t iVar2, Int_t iVar3, const Double_t *varMin, const Double_t *varMax, Bool_t useBins) const
819 // return a slice (3D-projection) on variables iVar1, iVar2 and iVar3 while axis ranges are defined with varMin,varMax
820 // arrays varMin and varMax contain the min and max values of each variable.
821 // therefore varMin and varMax must have their dimensions equal to GetNVar()
822 // If useBins=true, varMin and varMax are taken as bin numbers
824 THnSparse* clone = (THnSparse*)fData->Clone();
825 for (Int_t iAxis=0; iAxis<GetNVar(); iAxis++) {
826 if (useBins) clone->GetAxis(iAxis)->SetRange((Int_t)varMin[iAxis],(Int_t)varMax[iAxis]);
827 else clone->GetAxis(iAxis)->SetRangeUser(varMin[iAxis],varMax[iAxis]);
829 return clone->Projection(iVar1,iVar2,iVar3);
832 //____________________________________________________________________
833 void AliCFGridSparse::SetRangeUser(Int_t iVar, Double_t varMin, Double_t varMax) {
835 // set range of axis iVar.
837 fData->GetAxis(iVar)->SetRangeUser(varMin,varMax);
838 AliWarning(Form("THnSparse axis %d range has been modified",iVar));
841 //____________________________________________________________________
842 void AliCFGridSparse::SetRangeUser(const Double_t *varMin, const Double_t *varMax) {
844 // set range of every axis. varMin and varMax must be of dimension GetNVar()
846 for (Int_t iAxis=0; iAxis<GetNVar() ; iAxis++) { // set new range for every axis
847 SetRangeUser(iAxis,varMin[iAxis],varMax[iAxis]);
849 AliWarning("THnSparse axes ranges have been modified");
852 //____________________________________________________________________
853 void AliCFGridSparse::UseAxisRange(Bool_t b) const {
854 for (Int_t iAxis=0; iAxis<GetNVar(); iAxis++) fData->GetAxis(iAxis)->SetBit(TAxis::kAxisRange,b);
857 //____________________________________________________________________
858 Float_t AliCFGridSparse::GetOverFlows(Int_t ivar, Bool_t exclusive) const
861 // Returns overflows in variable ivar
862 // Set 'exclusive' to true for an exclusive check on variable ivar
864 Int_t* bin = new Int_t[GetNVar()];
865 memset(bin, 0, sizeof(Int_t) * GetNVar());
867 for (Long64_t i = 0; i < fData->GetNbins(); i++) {
868 Double_t v = fData->GetBinContent(i, bin);
871 for(Int_t j=0;j<GetNVar();j++){
873 if((bin[j]==0) || (bin[j]==GetNBins(j)+1))add=kFALSE;
876 if(bin[ivar]==GetNBins(ivar)+1 && add) ovfl+=v;
883 //____________________________________________________________________
884 Float_t AliCFGridSparse::GetUnderFlows(Int_t ivar, Bool_t exclusive) const
887 // Returns exclusive overflows in variable ivar
888 // Set 'exclusive' to true for an exclusive check on variable ivar
890 Int_t* bin = new Int_t[GetNVar()];
891 memset(bin, 0, sizeof(Int_t) * GetNVar());
893 for (Long64_t i = 0; i < fData->GetNbins(); i++) {
894 Double_t v = fData->GetBinContent(i, bin);
897 for(Int_t j=0;j<GetNVar();j++){
899 if((bin[j]==0) || (bin[j]==GetNBins(j)+1))add=kFALSE;
902 if(bin[ivar]==0 && add) unfl+=v;