]> git.uio.no Git - u/mrichter/AliRoot.git/blob - CORRFW/AliCFGridSparse.cxx
consistent naming
[u/mrichter/AliRoot.git] / CORRFW / AliCFGridSparse.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
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 //--------------------------------------------------------------------//
16 //                                                                    //
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 //--------------------------------------------------------------------//
25 //
26 //
27 #include "AliCFGridSparse.h"
28 #include "THnSparse.h"
29 #include "AliLog.h"
30 #include "TMath.h"
31 #include "TROOT.h"
32 #include "TH1D.h"
33 #include "TH2D.h"
34 #include "TH3D.h"
35 #include "TAxis.h"
36
37 //____________________________________________________________________
38 ClassImp(AliCFGridSparse)
39
40 //____________________________________________________________________
41 AliCFGridSparse::AliCFGridSparse() : 
42   AliCFVGrid(),
43   fData(0x0)
44 {
45   // default constructor
46 }
47 //____________________________________________________________________
48 AliCFGridSparse::AliCFGridSparse(const Char_t* name, const Char_t* title) : 
49   AliCFVGrid(name,title),
50   fData(0x0)
51 {
52   // default constructor
53 }
54 //____________________________________________________________________
55 AliCFGridSparse::AliCFGridSparse(const Char_t* name, const Char_t* title, const Int_t nVarIn, const Int_t * nBinIn, const Double_t *binLimitsIn) :  
56   AliCFVGrid(name,title,nVarIn,nBinIn,binLimitsIn),
57   fData(0x0)
58 {
59   //
60   // main constructor
61   //
62
63   fData=new THnSparseF(name,title,fNVar,fNVarBins);
64   
65   if(binLimitsIn){
66     for(Int_t ivar=0;ivar<fNVar;ivar++){    
67       Int_t nbins=fNVarBins[ivar]+1;
68       Double_t *array= new Double_t[nbins];
69       for(Int_t i=0;i<nbins;i++){
70         array[i]=fVarBinLimits[fOffset[ivar]+i];
71       } 
72       fData->SetBinEdges(ivar, array);
73       delete array;
74     }
75   }
76 }
77 //____________________________________________________________________
78 AliCFGridSparse::AliCFGridSparse(const AliCFGridSparse& c) : 
79   AliCFVGrid(c),
80   fData(c.fData)
81 {
82   //
83   // copy constructor
84   //
85   ((AliCFGridSparse &)c).Copy(*this);
86 }
87
88 //____________________________________________________________________
89 AliCFGridSparse::~AliCFGridSparse()
90 {
91   //
92   // destructor
93   //
94   if(fData) delete fData;
95 }
96
97 //____________________________________________________________________
98 AliCFGridSparse &AliCFGridSparse::operator=(const AliCFGridSparse &c)
99 {
100   //
101   // assigment operator
102   //
103   if (this != &c)
104     ((AliCFGridSparse &) c).Copy(*this);
105   return *this;
106
107
108 //____________________________________________________________________
109 void AliCFGridSparse::SetBinLimits(Int_t ivar, Double_t *array)
110 {
111   //
112   // setting the arrays containing the bin limits 
113   //
114   fData->SetBinEdges(ivar, array);
115   //then fill the appropriate array in ALICFFrame, to be able to use 
116   //the getter, in case....
117   Int_t nbins=fNVarBins[ivar]+1;
118   for(Int_t i=0;i<nbins;i++){
119     fVarBinLimits[fOffset[ivar]+i] =array[i];
120   } 
121
122
123 //____________________________________________________________________
124 void AliCFGridSparse::Fill(Double_t *var, Double_t weight)
125 {
126   //
127   // Fill the grid,
128   // given a set of values of the input variable, 
129   // with weight (by default w=1)
130   //
131   fData->Fill(var,weight);
132 }
133
134 //___________________________________________________________________
135 TH1D *AliCFGridSparse::Project(Int_t ivar) const
136 {
137   //
138   // Make a 1D projection along variable ivar 
139   //
140
141   TH1D *hist=fData->Projection(ivar);
142   return hist;
143 }
144 //___________________________________________________________________
145 TH2D *AliCFGridSparse::Project(Int_t ivar1, Int_t ivar2) const
146 {
147   //
148   // Make a 2D projection along variables ivar1 & ivar2 
149   //
150
151   TH2D *hist=fData->Projection(ivar2,ivar1); //notice inverted axis (THnSparse uses TH3 2d-projection convention...)
152   return hist;
153
154 }
155 //___________________________________________________________________
156 TH3D *AliCFGridSparse::Project(Int_t ivar1, Int_t ivar2, Int_t ivar3) const
157 {
158   //
159   // Make a 3D projection along variables ivar1 & ivar2 & ivar3 
160   //
161
162   TH3D *hist=fData->Projection(ivar1,ivar2,ivar3); 
163   return hist;
164
165 }
166
167 //___________________________________________________________________
168 AliCFGridSparse* AliCFGridSparse::Project(Int_t nVars, Int_t* vars, Double_t* varMin, Double_t* varMax) const
169 {
170
171   // binning for new grid
172   Int_t* bins = new Int_t[nVars];
173   for (Int_t iVar=0; iVar<nVars; iVar++) {
174     bins[iVar] = fNVarBins[vars[iVar]];
175   }
176   
177   // create new grid sparse
178   AliCFGridSparse* out = new AliCFGridSparse(fName,fTitle,nVars,bins);
179
180   //set the range in the THnSparse to project
181   THnSparse* clone = ((THnSparse*)fData->Clone());
182   for (Int_t iAxis=0; iAxis<fNVar; iAxis++) {
183     clone->GetAxis(iAxis)->SetRangeUser(varMin[iAxis],varMax[iAxis]);
184   }
185   out->SetGrid(clone->Projection(nVars,vars));
186   return out;
187 }
188
189 //____________________________________________________________________
190 Float_t AliCFGridSparse::GetOverFlows(Int_t ivar) const
191 {
192   //
193   // Returns exclusive overflows in variable ivar
194   //
195   Int_t* bin = new Int_t[fNVar];
196   memset(bin, 0, sizeof(Int_t) * fNVar);
197   Float_t ovfl=0.;
198   for (Long64_t i = 0; i < fData->GetNbins(); ++i) {
199     Double_t v = fData->GetBinContent(i, bin);
200     Bool_t add=kTRUE;
201     for(Int_t j=0;j<fNVar;j++){
202       if(ivar==j)continue;
203       if((bin[j]==0) || (bin[j]==fNVarBins[j]+1))add=kFALSE;
204     }
205     if(bin[ivar]==fNVarBins[ivar]+1 && add) ovfl+=v;
206   }
207
208   delete[] bin;
209   return ovfl;
210 }
211
212 //____________________________________________________________________
213 Float_t AliCFGridSparse::GetUnderFlows(Int_t ivar) const
214 {
215   //
216   // Returns exclusive overflows in variable ivar
217   //
218   Int_t* bin = new Int_t[fNVar];
219   memset(bin, 0, sizeof(Int_t) * fNVar);
220   Float_t unfl=0.;
221   for (Long64_t i = 0; i < fData->GetNbins(); ++i) {
222     Double_t v = fData->GetBinContent(i, bin);
223     Bool_t add=kTRUE;
224     for(Int_t j=0;j<fNVar;j++){
225       if(ivar==j)continue;
226       if((bin[j]==0) || (bin[j]==fNVarBins[j]+1))add=kFALSE;
227     }
228     if(bin[ivar]==0 && add) unfl+=v;
229   }
230
231   delete[] bin;
232   return unfl;
233 }
234
235
236 //____________________________________________________________________
237 Float_t AliCFGridSparse::GetEntries() const
238 {
239   //
240   // total entries (including overflows and underflows)
241   //
242
243   return fData->GetEntries();
244 }
245
246 //____________________________________________________________________
247 Float_t AliCFGridSparse::GetElement(Int_t index) const
248 {
249   //
250   // Returns content of grid element index according to the
251   // linear indexing in AliCFFrame
252   //
253   Int_t *bin = new Int_t[fNVar];
254   GetBinIndex(index, bin);
255   for(Int_t i=0;i<fNVar;i++)fIndex[i]=bin[i]+1; //consistency with AliCFGrid
256   Float_t val=GetElement(fIndex);
257   delete [] bin;
258   return val; 
259   
260 }
261 //____________________________________________________________________
262 Float_t AliCFGridSparse::GetElement(Int_t *bin) const
263 {
264   //
265   // Get the content in a bin corresponding to a set of bin indexes
266   //
267   return fData->GetBinContent(bin);
268
269 }  
270 //____________________________________________________________________
271 Float_t AliCFGridSparse::GetElement(Double_t *var) const
272 {
273   //
274   // Get the content in a bin corresponding to a set of input variables
275   //
276
277   Long_t index=fData->GetBin(var,kFALSE); //this is the THnSparse index (do not allocate new cells if content is empty)
278   if(index<0){
279     return 0.;
280   }else{
281     return fData->GetBinContent(index);
282   }
283
284
285 //____________________________________________________________________
286 Float_t AliCFGridSparse::GetElementError(Int_t index) const
287 {
288   //
289   // Returns the error on the content of a bin according to a linear
290   // indexing in AliCFFrame
291   //
292
293   Int_t *bin = new Int_t[fNVar];
294   GetBinIndex(index, bin);
295   for(Int_t i=0;i<fNVar;i++)fIndex[i]=bin[i]+1; //consistency with AliCFGrid
296   Float_t val=GetElementError(fIndex);
297   delete [] bin;
298   return val;
299
300 }
301 //____________________________________________________________________
302 Float_t AliCFGridSparse::GetElementError(Int_t *bin) const
303 {
304  //
305   // Get the error in a bin corresponding to a set of bin indexes
306   //
307   return fData->GetBinError(bin);
308
309 }  
310 //____________________________________________________________________
311 Float_t AliCFGridSparse::GetElementError(Double_t *var) const
312 {
313   //
314   // Get the error in a bin corresponding to a set of input variables
315   //
316
317   Long_t index=fData->GetBin(var,kFALSE); //this is the THnSparse index (do not allocate new cells if content is empy)
318   if(index<0){
319     return 0.;
320   }else{
321     return fData->GetBinError(index);
322   }
323
324
325
326 //____________________________________________________________________
327 void AliCFGridSparse::SetElement(Int_t index, Float_t val)
328 {
329   //
330   // Sets grid element iel to val (linear indexing) in AliCFFrame
331   //
332   Int_t *bin = new Int_t[fNVar];
333   GetBinIndex(index, bin);
334   for(Int_t i=0;i<fNVar;i++)fIndex[i]=bin[i]+1;
335   SetElement(fIndex,val);
336   delete [] bin;
337 }
338 //____________________________________________________________________
339 void AliCFGridSparse::SetElement(Int_t *bin, Float_t val)
340 {
341   //
342   // Sets grid element of bin indeces bin to val
343   //
344   fData->SetBinContent(bin,val);
345 }
346 //____________________________________________________________________
347 void AliCFGridSparse::SetElement(Double_t *var, Float_t val) 
348 {
349   //
350   // Set the content in a bin to value val corresponding to a set of input variables
351   //
352   Long_t index=fData->GetBin(var); //THnSparse index: allocate the cell
353   Int_t *bin = new Int_t[fNVar];
354   fData->GetBinContent(index,bin); //trick to access the array of bins
355   fData->SetBinContent(bin,val);
356   delete [] bin;
357
358 }
359
360 //____________________________________________________________________
361 void AliCFGridSparse::SetElementError(Int_t index, Float_t val)
362 {
363   //
364   // Sets grid element iel error to val (linear indexing) in AliCFFrame
365   //
366   Int_t *bin = new Int_t[fNVar];
367   GetBinIndex(index, bin);
368   for(Int_t i=0;i<fNVar;i++)fIndex[i]=bin[i]+1;
369   SetElementError(fIndex,val);
370   delete [] bin;
371 }
372 //____________________________________________________________________
373 void AliCFGridSparse::SetElementError(Int_t *bin, Float_t val)
374 {
375   //
376   // Sets grid element error of bin indeces bin to val
377   //
378   fData->SetBinError(bin,val);
379 }
380 //____________________________________________________________________
381 void AliCFGridSparse::SetElementError(Double_t *var, Float_t val) 
382 {
383   //
384   // Set the error in a bin to value val corresponding to a set of input variables
385   //
386   Long_t index=fData->GetBin(var); //THnSparse index
387   Int_t *bin = new Int_t[fNVar];
388   fData->GetBinContent(index,bin); //trick to access the array of bins
389   fData->SetBinError(bin,val);
390   delete [] bin;
391 }
392
393 //____________________________________________________________________
394 void AliCFGridSparse::SumW2()
395 {
396   //
397   //set calculation of the squared sum of the weighted entries
398   //
399   if(!fSumW2){
400     fData->CalculateErrors(kTRUE); 
401   }
402
403   fSumW2=kTRUE;
404 }
405
406 //____________________________________________________________________
407 void AliCFGridSparse::Add(AliCFVGrid* aGrid, Double_t c)
408 {
409   //
410   //add aGrid to the current one
411   //
412
413
414   if(aGrid->GetNVar()!=fNVar){
415     AliInfo("Different number of variables, cannot add the grids");
416     return;
417   } 
418   if(aGrid->GetNDim()!=fNDim){
419     AliInfo("Different number of dimensions, cannot add the grids!");
420     return;
421   } 
422   
423   if(!fSumW2  && aGrid->GetSumW2())SumW2();
424
425
426   fData->Add(((AliCFGridSparse*)aGrid)->GetGrid(),c);
427
428 }
429
430 //____________________________________________________________________
431 void AliCFGridSparse::Add(AliCFVGrid* aGrid1, AliCFVGrid* aGrid2, Double_t c1,Double_t c2)
432 {
433   //
434   //Add aGrid1 and aGrid2 and deposit the result into the current one
435   //
436
437   if(fNVar!=aGrid1->GetNVar()|| fNVar!=aGrid2->GetNVar()){
438     AliInfo("Different number of variables, cannot add the grids");
439     return;
440   } 
441   if(fNDim!=aGrid1->GetNDim()|| fNDim!=aGrid2->GetNDim()){
442     AliInfo("Different number of dimensions, cannot add the grids!");
443     return;
444   } 
445   
446   if(!fSumW2  && (aGrid1->GetSumW2() || aGrid2->GetSumW2()))SumW2();
447
448
449   fData->Reset();
450   fData->Add(((AliCFGridSparse*)aGrid1)->GetGrid(),c1);
451   fData->Add(((AliCFGridSparse*)aGrid2)->GetGrid(),c2);
452
453 }
454
455 //____________________________________________________________________
456 void AliCFGridSparse::Multiply(AliCFVGrid* aGrid, Double_t c)
457 {
458   //
459   // Multiply aGrid to the current one
460   //
461
462
463   if(aGrid->GetNVar()!=fNVar){
464     AliInfo("Different number of variables, cannot multiply the grids");
465     return;
466   } 
467   if(aGrid->GetNDim()!=fNDim){
468     AliInfo("Different number of dimensions, cannot multiply the grids!");
469     return;
470   } 
471   
472   if(!fSumW2  && aGrid->GetSumW2())SumW2();
473
474   THnSparse *h= ((AliCFGridSparse*)aGrid)->GetGrid();
475
476   fData->Multiply(h);
477   fData->Scale(c);
478
479 }
480
481 //____________________________________________________________________
482 void AliCFGridSparse::Multiply(AliCFVGrid* aGrid1, AliCFVGrid* aGrid2, Double_t c1,Double_t c2)
483 {
484   //
485   //Multiply aGrid1 and aGrid2 and deposit the result into the current one
486   //
487
488   if(fNVar!=aGrid1->GetNVar()|| fNVar!=aGrid2->GetNVar()){
489     AliInfo("Different number of variables, cannot multiply the grids");
490     return;
491   } 
492   if(fNDim!=aGrid1->GetNDim()|| fNDim!=aGrid2->GetNDim()){
493     AliInfo("Different number of dimensions, cannot multiply the grids!");
494     return;
495   } 
496   
497   if(!fSumW2  && (aGrid1->GetSumW2() || aGrid2->GetSumW2()))SumW2();
498
499
500   fData->Reset();
501   THnSparse *h1= ((AliCFGridSparse*)aGrid1)->GetGrid();
502   THnSparse *h2= ((AliCFGridSparse*)aGrid2)->GetGrid();
503   h2->Multiply(h1);
504   h2->Scale(c1*c2);
505   fData->Add(h2);
506 }
507
508
509
510 //____________________________________________________________________
511 void AliCFGridSparse::Divide(AliCFVGrid* aGrid, Double_t c)
512 {
513   //
514   // Divide aGrid to the current one
515   //
516
517
518   if(aGrid->GetNVar()!=fNVar){
519     AliInfo("Different number of variables, cannot divide the grids");
520     return;
521   } 
522   if(aGrid->GetNDim()!=fNDim){
523     AliInfo("Different number of dimensions, cannot divide the grids!");
524     return;
525   } 
526   
527   if(!fSumW2  && aGrid->GetSumW2())SumW2();
528
529   THnSparse *h= ((AliCFGridSparse*)aGrid)->GetGrid();
530
531   fData->Divide(h);
532   fData->Scale(c);
533
534 }
535
536 //____________________________________________________________________
537 void AliCFGridSparse::Divide(AliCFVGrid* aGrid1, AliCFVGrid* aGrid2, Double_t c1,Double_t c2, Option_t *option)
538 {
539   //
540   //Divide aGrid1 and aGrid2 and deposit the result into the current one
541   //bynomial errors are supported
542   //
543
544   if(fNVar!=aGrid1->GetNVar()|| fNVar!=aGrid2->GetNVar()){
545     AliInfo("Different number of variables, cannot divide the grids");
546     return;
547   } 
548   if(fNDim!=aGrid1->GetNDim()|| fNDim!=aGrid2->GetNDim()){
549     AliInfo("Different number of dimensions, cannot divide the grids!");
550     return;
551   } 
552   
553   if(!fSumW2  && (aGrid1->GetSumW2() || aGrid2->GetSumW2()))SumW2();
554
555
556   THnSparse *h1= ((AliCFGridSparse*)aGrid1)->GetGrid();
557   THnSparse *h2= ((AliCFGridSparse*)aGrid2)->GetGrid();
558   fData->Divide(h1,h2,c1,c2,option);
559 }
560
561
562 //____________________________________________________________________
563 void AliCFGridSparse::Rebin(const Int_t* group)
564 {
565   //
566   // rebin the grid according to Rebin() as in THnSparse
567   // Please notice that the original number of bins on
568   // a given axis has to be divisible by the rebin group.
569   //
570
571   for(Int_t i=0;i<fNVar;i++){
572     if(group[i]!=1)AliInfo(Form(" merging bins along dimension %i in groups of %i bins", i,group[i]));
573   }
574
575   THnSparse *rebinned =fData->Rebin(group);
576   fData->Reset();
577   fData = rebinned;
578
579   //redefine the needed stuff
580
581   Int_t ndimTot=1;
582   Int_t nbinTot=0;
583
584   //number of bins in each dimension, auxiliary variables
585
586   for(Int_t ivar=0;ivar<fNVar;ivar++){
587     Int_t nbins = fData->GetAxis(ivar)->GetNbins();
588     fNVarBins[ivar]=nbins;
589     ndimTot*=fNVarBins[ivar];
590     nbinTot+=(fNVarBins[ivar]+1);
591     Int_t offset=0;
592     for(Int_t i =0;i<ivar;i++)offset+=(fNVarBins[i]+1);      
593     fOffset[ivar]=offset;
594     Int_t prod=1;
595     for(Int_t i=0;i<ivar;i++)prod*=fNVarBins[i];
596     fProduct[ivar]=prod;
597   }
598
599   fNDim=ndimTot;
600
601   //now the array of bin limits
602
603   delete fVarBinLimits;
604   fNVarBinLimits=nbinTot;
605   fVarBinLimits=new Double_t[fNVarBinLimits];
606
607   for(Int_t ivar=0;ivar<fNVar;ivar++){
608     Double_t low = fData->GetAxis(ivar)->GetXmin();
609     Double_t high = fData->GetAxis(ivar)->GetXmax();    
610     const TArrayD *xbins = fData->GetAxis(ivar)->GetXbins();
611     if (xbins->fN == 0){
612       for(Int_t ibin=0;ibin<=fNVarBins[ivar];ibin++){
613         fVarBinLimits[ibin+fOffset[ivar]] = low + ibin*(high-low)/((Double_t) fNVarBins[ivar]);
614       }
615     }
616     else{
617       
618       for(Int_t ibin=0;ibin<=fNVarBins[ivar];ibin++) {
619         fVarBinLimits[ibin+fOffset[ivar]] = xbins->At(ibin);
620       }
621     }
622   }   
623   
624 }
625 //____________________________________________________________________
626 void AliCFGridSparse::Copy(TObject& c) const
627 {
628   //
629   // copy function
630   //
631   AliCFGridSparse& target = (AliCFGridSparse &) c;
632
633   if(fData)target.fData = fData;
634 }
635
636 //____________________________________________________________________
637 TH1D* AliCFGridSparse::Slice(Int_t iVar, Double_t *varMin, Double_t *varMax) const
638 {
639   //
640   // return a slice (1D-projection) on variable iVar while axis ranges are defined with varMin,varMax
641   // arrays varMin and varMax contain the min and max values of each variable.
642   // therefore varMin and varMax must have their dimensions equal to fNVar
643   //
644   
645   THnSparse* clone = (THnSparse*)fData->Clone();
646   for (Int_t iAxis=0; iAxis<fNVar; iAxis++) {
647     clone->GetAxis(iAxis)->SetRangeUser(varMin[iAxis],varMax[iAxis]);
648   }
649   return clone->Projection(iVar);
650 }
651
652 //____________________________________________________________________
653 TH2D* AliCFGridSparse::Slice(Int_t iVar1, Int_t iVar2, Double_t *varMin, Double_t *varMax) const
654 {
655   //
656   // return a slice (2D-projection) on variables iVar1 and iVar2 while axis ranges are defined with varMin,varMax
657   // arrays varMin and varMax contain the min and max values of each variable.
658   // therefore varMin and varMax must have their dimensions equal to fNVar
659   //
660   
661   THnSparse* clone = (THnSparse*)fData->Clone();
662   for (Int_t iAxis=0; iAxis<fNVar; iAxis++) {
663     clone->GetAxis(iAxis)->SetRangeUser(varMin[iAxis],varMax[iAxis]);
664   }
665   return clone->Projection(iVar1,iVar2);
666 }
667
668 //____________________________________________________________________
669 TH3D* AliCFGridSparse::Slice(Int_t iVar1, Int_t iVar2, Int_t iVar3, Double_t *varMin, Double_t *varMax) const
670 {
671   //
672   // return a slice (3D-projection) on variables iVar1, iVar2 and iVar3 while axis ranges are defined with varMin,varMax
673   // arrays varMin and varMax contain the min and max values of each variable.
674   // therefore varMin and varMax must have their dimensions equal to fNVar
675   //
676
677   THnSparse* clone = (THnSparse*)fData->Clone();
678   for (Int_t iAxis=0; iAxis<fNVar; iAxis++) {
679     clone->GetAxis(iAxis)->SetRangeUser(varMin[iAxis],varMax[iAxis]);
680   }
681   return clone->Projection(iVar1,iVar2,iVar3);
682 }
683
684 //____________________________________________________________________
685 void AliCFGridSparse::SetRangeUser(Int_t iVar, Double_t varMin, Double_t varMax) {
686   //
687   // set range of axis iVar. 
688   //
689   fData->GetAxis(iVar)->SetRangeUser(varMin,varMax);
690   AliWarning(Form("THnSparse axis %d range has been modified",iVar));
691 }
692
693 //____________________________________________________________________
694 void AliCFGridSparse::SetRangeUser(Double_t *varMin, Double_t *varMax) {
695   //
696   // set range of every axis. varMin and varMax must be of dimension fNVar
697   //
698   for (Int_t iAxis=0; iAxis<fNVar ; iAxis++) { // set new range for every axis
699     SetRangeUser(iAxis,varMin[iAxis],varMax[iAxis]);
700   }
701   AliWarning("THnSparse axes ranges have been modified");
702 }
703
704 //____________________________________________________________________
705 void AliCFGridSparse::UseAxisRange(Bool_t b) const {
706   for (Int_t iAxis=0; iAxis<fNVar; iAxis++) fData->GetAxis(iAxis)->SetBit(TAxis::kAxisRange,b);
707 }