]> git.uio.no Git - u/mrichter/AliRoot.git/blob - CORRFW/AliCFGridSparse.cxx
fix of memory potential leak
[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   AliCFFrame(),
43   fSumW2(kFALSE),
44   fData(0x0)
45 {
46   // default constructor
47 }
48 //____________________________________________________________________
49 AliCFGridSparse::AliCFGridSparse(const Char_t* name, const Char_t* title) : 
50   AliCFFrame(name,title),
51   fSumW2(kFALSE),
52   fData(0x0)
53 {
54   // default constructor
55 }
56 //____________________________________________________________________
57 AliCFGridSparse::AliCFGridSparse(const Char_t* name, const Char_t* title, Int_t nVarIn, const Int_t * nBinIn) :  
58   AliCFFrame(name,title),
59   fSumW2(kFALSE),
60   fData(0x0)
61 {
62   //
63   // main constructor
64   //
65
66   fData = new THnSparseF(name,title,nVarIn,nBinIn);
67 }
68
69 //____________________________________________________________________
70 AliCFGridSparse::~AliCFGridSparse()
71 {
72   //
73   // destructor
74   //
75   if (fData) delete fData;
76 }
77
78 //____________________________________________________________________
79 AliCFGridSparse::AliCFGridSparse(const AliCFGridSparse& c) :
80   AliCFFrame(c),
81   fSumW2(kFALSE),
82   fData(0x0)
83 {
84   //
85   // copy constructor
86   //
87   ((AliCFGridSparse &)c).Copy(*this);
88 }
89
90 //____________________________________________________________________
91 AliCFGridSparse& AliCFGridSparse::operator=(const AliCFGridSparse &c)
92 {
93   //
94   // assigment operator
95   //
96   if (this != &c) c.Copy(*this);
97   return *this;
98
99
100 //____________________________________________________________________
101 void AliCFGridSparse::SetBinLimits(Int_t ivar, Double_t min, Double_t max)
102 {
103   //
104   // set a uniform binning for variable ivar
105   //
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);
110   delete [] array ;
111
112
113 //____________________________________________________________________
114 void AliCFGridSparse::SetBinLimits(Int_t ivar, const Double_t *array)
115 {
116   //
117   // setting the arrays containing the bin limits 
118   //
119   fData->SetBinEdges(ivar, array);
120
121
122 //____________________________________________________________________
123 void AliCFGridSparse::Fill(const Double_t *var, Double_t weight)
124 {
125   //
126   // Fill the grid,
127   // given a set of values of the input variable, 
128   // with weight (by default w=1)
129   //
130   fData->Fill(var,weight);
131 }
132
133 //___________________________________________________________________
134 TH1D *AliCFGridSparse::Project(Int_t ivar) const
135 {
136   //
137   // Make a 1D projection along variable ivar 
138   //
139
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)));
144   return hist;
145 }
146 //___________________________________________________________________
147 TH2D *AliCFGridSparse::Project(Int_t ivar1, Int_t ivar2) const
148 {
149   //
150   // Make a 2D projection along variables ivar1 & ivar2 
151   //
152
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)));
158   return hist;
159 }
160 //___________________________________________________________________
161 TH3D *AliCFGridSparse::Project(Int_t ivar1, Int_t ivar2, Int_t ivar3) const
162 {
163   //
164   // Make a 3D projection along variables ivar1 & ivar2 & ivar3 
165   //
166
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)));
173   return hist;
174 }
175
176 //___________________________________________________________________
177 AliCFGridSparse* AliCFGridSparse::Project(Int_t nVars, const Int_t* vars, const Double_t* varMin, const Double_t* varMax, Bool_t useBins) const
178 {
179   //
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
183   //
184
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]);
189   }
190   
191   // create new grid sparse
192   AliCFGridSparse* out = new AliCFGridSparse(fName,fTitle,nVars,bins);
193
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]);
200     }
201   }
202   else AliInfo("Keeping same axis ranges");
203
204   out->SetGrid(clone->Projection(nVars,vars));
205   delete [] bins;
206   return out;
207 }
208
209
210 //____________________________________________________________________
211 Float_t AliCFGridSparse::GetBinCenter(Int_t ivar, Int_t ibin) const
212 {
213   //
214   // Returns the center of specified bin for variable axis ivar
215   // 
216   
217   return (Float_t) fData->GetAxis(ivar)->GetBinCenter(ibin);
218 }
219
220 //____________________________________________________________________
221 Float_t AliCFGridSparse::GetBinSize(Int_t ivar, Int_t ibin) const
222 {
223   //
224   // Returns the size of specified bin for variable axis ivar
225   // 
226   
227   return (Float_t) fData->GetAxis(ivar)->GetBinUpEdge(ibin) - fData->GetAxis(ivar)->GetBinLowEdge(ibin);
228 }
229
230 //____________________________________________________________________
231 Float_t AliCFGridSparse::GetEntries() const
232 {
233   //
234   // total entries (including overflows and underflows)
235   //
236
237   return fData->GetEntries();
238 }
239
240 //____________________________________________________________________
241 Float_t AliCFGridSparse::GetElement(Long_t index) const
242 {
243   //
244   // Returns content of grid element index 
245   //
246   
247   return fData->GetBinContent(index);
248 }
249 //____________________________________________________________________
250 Float_t AliCFGridSparse::GetElement(const Int_t *bin) const
251 {
252   //
253   // Get the content in a bin corresponding to a set of bin indexes
254   //
255   return fData->GetBinContent(bin);
256
257 }  
258 //____________________________________________________________________
259 Float_t AliCFGridSparse::GetElement(const Double_t *var) const
260 {
261   //
262   // Get the content in a bin corresponding to a set of input variables
263   //
264
265   Long_t index = fData->GetBin(var,kFALSE);
266   if (index<0) return 0.;
267   return fData->GetBinContent(index);
268
269
270 //____________________________________________________________________
271 Float_t AliCFGridSparse::GetElementError(Long_t index) const
272 {
273   //
274   // Returns the error on the content 
275   //
276
277   return fData->GetBinContent(index);
278 }
279 //____________________________________________________________________
280 Float_t AliCFGridSparse::GetElementError(const Int_t *bin) const
281 {
282  //
283   // Get the error in a bin corresponding to a set of bin indexes
284   //
285   return fData->GetBinError(bin);
286
287 }  
288 //____________________________________________________________________
289 Float_t AliCFGridSparse::GetElementError(const Double_t *var) const
290 {
291   //
292   // Get the error in a bin corresponding to a set of input variables
293   //
294
295   Long_t index=fData->GetBin(var,kFALSE); //this is the THnSparse index (do not allocate new cells if content is empy)
296   if (index<0) return 0.;
297   return fData->GetBinError(index);
298
299
300 //____________________________________________________________________
301 void AliCFGridSparse::SetElement(Long_t index, Float_t val)
302 {
303   //
304   // Sets grid element value
305   //
306   Int_t* bin = new Int_t[GetNVar()];
307   fData->GetBinContent(index,bin); //affects the bin coordinates
308   SetElement(bin,val);
309   delete [] bin ;
310 }
311
312 //____________________________________________________________________
313 void AliCFGridSparse::SetElement(const Int_t *bin, Float_t val)
314 {
315   //
316   // Sets grid element of bin indeces bin to val
317   //
318   fData->SetBinContent(bin,val);
319 }
320 //____________________________________________________________________
321 void AliCFGridSparse::SetElement(const Double_t *var, Float_t val) 
322 {
323   //
324   // Set the content in a bin to value val corresponding to a set of input variables
325   //
326   Long_t index=fData->GetBin(var,kTRUE); //THnSparse index: allocate the cell
327   Int_t *bin = new Int_t[GetNVar()];
328   fData->GetBinContent(index,bin); //trick to access the array of bins
329   SetElement(bin,val);
330   delete [] bin;
331 }
332
333 //____________________________________________________________________
334 void AliCFGridSparse::SetElementError(Long_t index, Float_t val)
335 {
336   //
337   // Sets grid element iel error to val (linear indexing) in AliCFFrame
338   //
339   Int_t *bin = new Int_t[GetNVar()];
340   fData->GetBinContent(index,bin);
341   SetElementError(bin,val);
342   delete [] bin;
343 }
344
345 //____________________________________________________________________
346 void AliCFGridSparse::SetElementError(const Int_t *bin, Float_t val)
347 {
348   //
349   // Sets grid element error of bin indeces bin to val
350   //
351   fData->SetBinError(bin,val);
352 }
353 //____________________________________________________________________
354 void AliCFGridSparse::SetElementError(const Double_t *var, Float_t val) 
355 {
356   //
357   // Set the error in a bin to value val corresponding to a set of input variables
358   //
359   Long_t index=fData->GetBin(var); //THnSparse index
360   Int_t *bin = new Int_t[GetNVar()];
361   fData->GetBinContent(index,bin); //trick to access the array of bins
362   SetElementError(bin,val);
363   delete [] bin;
364 }
365
366 //____________________________________________________________________
367 void AliCFGridSparse::SumW2()
368 {
369   //
370   //set calculation of the squared sum of the weighted entries
371   //
372   if(!fSumW2){
373     fData->CalculateErrors(kTRUE); 
374   }
375   fSumW2=kTRUE;
376 }
377
378 //____________________________________________________________________
379 void AliCFGridSparse::Add(const AliCFGridSparse* aGrid, Double_t c)
380 {
381   //
382   //add aGrid to the current one
383   //
384
385   if (aGrid->GetNVar() != GetNVar()){
386     AliError("Different number of variables, cannot add the grids");
387     return;
388   } 
389   
390   if (!fSumW2  && aGrid->GetSumW2()) SumW2();
391   fData->Add(aGrid->GetGrid(),c);
392 }
393
394 //____________________________________________________________________
395 void AliCFGridSparse::Add(const AliCFGridSparse* aGrid1, const AliCFGridSparse* aGrid2, Double_t c1,Double_t c2)
396 {
397   //
398   //Add aGrid1 and aGrid2 and deposit the result into the current one
399   //
400
401   if (GetNVar() != aGrid1->GetNVar() || GetNVar() != aGrid2->GetNVar()) {
402     AliInfo("Different number of variables, cannot add the grids");
403     return;
404   } 
405   
406   if (!fSumW2  && (aGrid1->GetSumW2() || aGrid2->GetSumW2())) SumW2();
407
408   fData->Reset();
409   fData->Add(aGrid1->GetGrid(),c1);
410   fData->Add(aGrid2->GetGrid(),c2);
411 }
412
413 //____________________________________________________________________
414 void AliCFGridSparse::Multiply(const AliCFGridSparse* aGrid, Double_t c)
415 {
416   //
417   // Multiply aGrid to the current one
418   //
419
420   if (aGrid->GetNVar() != GetNVar()) {
421     AliError("Different number of variables, cannot multiply the grids");
422     return;
423   } 
424   
425   if(!fSumW2  && aGrid->GetSumW2()) SumW2();
426   THnSparse *h = aGrid->GetGrid();
427   fData->Multiply(h);
428   fData->Scale(c);
429 }
430
431 //____________________________________________________________________
432 void AliCFGridSparse::Multiply(const AliCFGridSparse* aGrid1, const AliCFGridSparse* aGrid2, Double_t c1,Double_t c2)
433 {
434   //
435   //Multiply aGrid1 and aGrid2 and deposit the result into the current one
436   //
437
438   if (GetNVar() != aGrid1->GetNVar() || GetNVar() != aGrid2->GetNVar()) {
439     AliError("Different number of variables, cannot multiply the grids");
440     return;
441   }
442   
443   if(!fSumW2  && (aGrid1->GetSumW2() || aGrid2->GetSumW2())) SumW2();
444
445   fData->Reset();
446   THnSparse *h1 = aGrid1->GetGrid();
447   THnSparse *h2 = aGrid2->GetGrid();
448   h2->Multiply(h1);
449   h2->Scale(c1*c2);
450   fData->Add(h2);
451 }
452
453 //____________________________________________________________________
454 void AliCFGridSparse::Divide(const AliCFGridSparse* aGrid, Double_t c)
455 {
456   //
457   // Divide aGrid to the current one
458   //
459
460   if (aGrid->GetNVar() != GetNVar()) {
461     AliError("Different number of variables, cannot divide the grids");
462     return;
463   } 
464   
465   if (!fSumW2  && aGrid->GetSumW2()) SumW2();
466
467   THnSparse *h1 = aGrid->GetGrid();
468   THnSparse *h2 = (THnSparse*)fData->Clone();
469   fData->Divide(h2,h1);
470   fData->Scale(c);
471 }
472
473 //____________________________________________________________________
474 void AliCFGridSparse::Divide(const AliCFGridSparse* aGrid1, const AliCFGridSparse* aGrid2, Double_t c1,Double_t c2, Option_t *option)
475 {
476   //
477   //Divide aGrid1 and aGrid2 and deposit the result into the current one
478   //binomial errors are supported
479   //
480
481   if (GetNVar() != aGrid1->GetNVar() || GetNVar() != aGrid2->GetNVar()) {
482     AliError("Different number of variables, cannot divide the grids");
483     return;
484   } 
485   
486   if (!fSumW2  && (aGrid1->GetSumW2() || aGrid2->GetSumW2())) SumW2();
487
488   THnSparse *h1= aGrid1->GetGrid();
489   THnSparse *h2= aGrid2->GetGrid();
490   fData->Divide(h1,h2,c1,c2,option);
491 }
492
493
494 //____________________________________________________________________
495 void AliCFGridSparse::Rebin(const Int_t* group)
496 {
497   //
498   // rebin the grid according to Rebin() as in THnSparse
499   // Please notice that the original number of bins on
500   // a given axis has to be divisible by the rebin group.
501   //
502
503   for(Int_t i=0;i<GetNVar();i++){
504     if (group[i]!=1) AliInfo(Form(" merging bins along dimension %i in groups of %i bins", i,group[i]));
505   }
506
507   THnSparse *rebinned =fData->Rebin(group);
508   fData->Reset();
509   fData = rebinned;
510 }
511 //____________________________________________________________________
512 void AliCFGridSparse::Scale(Long_t index, const Double_t *fact)
513 {
514   //
515   //scale content of a certain cell by (positive) fact (with error)
516   //
517
518   if (GetElement(index)==0 || fact[0]==0) return;
519
520   Double_t in[2], out[2];
521   in[0]=GetElement(index);
522   in[1]=GetElementError(index);
523   GetScaledValues(fact,in,out);
524   SetElement(index,out[0]);
525   if (fSumW2) SetElementError(index,out[1]);
526 }
527 //____________________________________________________________________
528 void AliCFGridSparse::Scale(const Int_t *bin, const Double_t *fact)
529 {
530   //
531   //scale content of a certain cell by (positive) fact (with error)
532   //
533   if(GetElement(bin)==0 || fact[0]==0)return;
534
535   Double_t in[2], out[2];
536   in[0]=GetElement(bin);
537   in[1]=GetElementError(bin);
538   GetScaledValues(fact,in,out);
539   SetElement(bin,out[0]);
540   if(fSumW2)SetElementError(bin,out[1]);
541   
542 }
543 //____________________________________________________________________
544 void AliCFGridSparse::Scale(const Double_t *var, const Double_t *fact) 
545 {
546   //
547   //scale content of a certain cell by (positive) fact (with error)
548   //
549   if(GetElement(var)==0 || fact[0]==0)return;
550
551   Double_t in[2], out[2];
552   in[0]=GetElement(var);
553   in[1]=GetElementError(var);
554   GetScaledValues(fact,in,out);
555   SetElement(var,out[0]);
556   if(fSumW2)SetElementError(var,out[1]);
557   
558 }
559 //____________________________________________________________________
560 void AliCFGridSparse::Scale(const Double_t* fact)
561 {
562   //
563   //scale contents of the whole grid by fact
564   //
565
566   for (Long_t iel=0; iel<GetNFilledBins(); iel++) {
567     Scale(iel,fact);
568   }
569 }
570 //____________________________________________________________________
571 Long_t AliCFGridSparse::GetEmptyBins() const {
572   //
573   // Get empty bins 
574   //
575
576   return (GetNBinsTotal() - GetNFilledBins()) ;
577
578
579 //_____________________________________________________________________
580 // Int_t AliCFGridSparse::GetEmptyBins( Double_t *varMin, Double_t* varMax ) const 
581 // {
582 //   //
583 //   // Get empty bins in a range specified by varMin and varMax
584 //   //
585
586 //   Int_t *indexMin = new Int_t[GetNVar()];
587 //   Int_t *indexMax = new Int_t[GetNVar()];
588
589 //   //Find out the min and max bins
590
591 //   for (Int_t i=0;i<GetNVar();i++) {
592 //     Double_t xmin=varMin[i]; // the min values  
593 //     Double_t xmax=varMax[i]; // the min values  
594 //     Int_t nbins = fData->GetAxis(i)->GetNbins()+1;
595 //     Double_t *bins=new Double_t[nbins];
596 //     for(Int_t ibin =0;ibin<nbins;ibin++){
597 //      bins[ibin] = fVarBinLimits[ibin+fOffset[i]];
598 //     }
599 //     indexMin[i] = TMath::BinarySearch(nbins,bins,xmin);
600 //     indexMax[i] = TMath::BinarySearch(nbins,bins,xmax);
601 //     delete [] bins;
602 //   }
603
604 //   Int_t val=0;
605 //   for(Int_t i=0;i<fNDim;i++){
606 //     for (Int_t j=0;j<GetNVar();j++)fIndex[j]=GetBinIndex(j,i);
607 //     Bool_t isIn=kTRUE;
608 //     for (Int_t j=0;j<GetNVar();j++){
609 //       if(!(fIndex[j]>=indexMin[j] && fIndex[j]<=indexMax[j]))isIn=kFALSE;   
610 //     }
611 //     if(isIn && GetElement(i)<=0)val++;     
612 //   }
613 //   AliInfo(Form(" the empty bins = %i ",val)); 
614
615 //   delete [] indexMin;
616 //   delete [] indexMax;
617 //   return val;
618 // } 
619
620 //____________________________________________________________________
621 Int_t AliCFGridSparse::CheckStats(Double_t thr) const
622 {
623   //
624   // Count the cells below a certain threshold
625   //
626   Int_t ncellsLow=0;
627   for (Int_t i=0; i<GetNBinsTotal(); i++) {
628     if (GetElement(i)<thr) ncellsLow++;
629   }
630   return ncellsLow;
631 }
632
633 //_____________________________________________________________________
634 Double_t AliCFGridSparse::GetIntegral() const 
635 {
636   //
637   // Get full Integral
638   //
639   return fData->ComputeIntegral();  
640
641
642 //_____________________________________________________________________
643 // Double_t AliCFGridSparse::GetIntegral(Int_t *binMin, Int_t* binMax ) const 
644 // {
645 //   //
646 //   // Get Integral in a range of bin indeces (extremes included)
647 //   //
648
649 //   Double_t val=0;
650
651 //   for(Int_t i=0;i<GetNVar();i++){
652 //     if(binMin[i]<1)binMin[i]=1;
653 //     if(binMax[i]>fNVarBins[i])binMax[i]=fNVarBins[i];
654 //     if((binMin[i]>binMax[i])){
655 //       AliInfo(Form(" Bin indices in variable %i in reverse order, please check!", i));
656 //       return val;
657 //     }
658 //   }
659 //   val=GetSum(0,binMin,binMax);
660
661 //   return val;
662 // } 
663
664 //_____________________________________________________________________
665 // Double_t AliCFGridSparse::GetIntegral(Double_t *varMin, Double_t* varMax ) const 
666 // {
667 //   //
668 //   // Get Integral in a range (extremes included)
669 //   //
670
671 //   Int_t *indexMin=new Int_t[GetNVar()];
672 //   Int_t *indexMax=new Int_t[GetNVar()];
673
674 //   //Find out the min and max bins
675
676 //   for(Int_t i=0;i<GetNVar();i++){
677 //     Double_t xmin=varMin[i]; // the min values  
678 //     Double_t xmax=varMax[i]; // the min values  
679 //     Int_t nbins=fNVarBins[i]+1;
680 //     Double_t *bins=new Double_t[nbins];
681 //     for(Int_t ibin =0;ibin<nbins;ibin++){
682 //      bins[ibin] = fVarBinLimits[ibin+fOffset[i]];
683 //     }
684 //     indexMin[i] = TMath::BinarySearch(nbins,bins,xmin);
685 //     indexMax[i] = TMath::BinarySearch(nbins,bins,xmax);
686 //     delete [] bins;
687 //   }
688
689 //   //move to the TH/THnSparse convention in N-dim bin numbering 
690 //   for(Int_t i=0;i<GetNVar();i++){
691 //     indexMin[i]+=1;
692 //     indexMax[i]+=1;
693 //   }
694
695 //   Double_t val=GetIntegral(indexMin,indexMax);
696
697 //   delete [] indexMin;
698 //   delete [] indexMax;
699
700 //   return val;
701 // } 
702
703 // //_____________________________________________________________________
704 // Double_t AliCFGridSparse::GetSum(Int_t ivar, Int_t *binMin, Int_t* binMax) const 
705 // {
706 //   //
707 //   // recursively add over nested loops.... 
708 //   //
709
710 //   static Double_t val;
711 //   if (ivar==0) val=0.;
712 //   for(Int_t ibin=binMin[ivar]-1 ; ibin<=binMax[ivar]-1 ; ibin++){
713 //     //-1 is to move from TH/ThnSparse N-dim bin convention to one in AliCFFrame
714 //     fIndex[ivar]=ibin;
715 //     if(ivar<GetNVar()-1) {
716 //       val=GetSum(ivar+1,binMin,binMax);
717 //     }
718 //     else {
719 //       Int_t iel=GetBinIndex(fIndex);
720 //       val+=GetElement(iel);
721 //     }
722 //   }
723   
724 //   return val;
725 // }
726
727 //____________________________________________________________________
728 Long64_t AliCFGridSparse::Merge(TCollection* list)
729 {
730   //
731   // Merge a list of AliCFGridSparse with this (needed for PROOF). 
732   // Returns the number of merged objects (including this).
733   //
734
735   if (!list)
736     return 0;
737   
738   if (list->IsEmpty())
739     return 1;
740
741   TIterator* iter = list->MakeIterator();
742   TObject* obj;
743   
744   Int_t count = 0;
745   while ((obj = iter->Next())) {
746     AliCFGridSparse* entry = dynamic_cast<AliCFGridSparse*> (obj);
747     if (entry == 0) 
748       continue;
749     this->Add(entry);
750     count++;
751   }
752
753   return count+1;
754 }
755
756 //____________________________________________________________________
757 void AliCFGridSparse::GetScaledValues(const Double_t *fact, const Double_t *in, Double_t *out) const{
758   //
759   // scale input *in and its error by (positive) fact (with error)
760   // and erite it to *out
761   //
762   out[0]=in[0]*fact[0];
763   out[1]=TMath::Sqrt(in[1]*in[1]/in[0]/in[0]
764                      +fact[1]*fact[1]/fact[0]/fact[0])*out[0];
765     
766 }
767
768 //____________________________________________________________________
769 void AliCFGridSparse::Copy(TObject& c) const
770 {
771   //
772   // copy function
773   //
774   AliCFFrame::Copy(c);
775   AliCFGridSparse& target = (AliCFGridSparse &) c;
776   target.fSumW2 = fSumW2 ;
777   if (fData) {
778     target.fData = (THnSparse*)fData->Clone();
779   }
780 }
781
782 //____________________________________________________________________
783 TH1D* AliCFGridSparse::Slice(Int_t iVar, const Double_t *varMin, const Double_t *varMax, Bool_t useBins) const
784 {
785   //
786   // return a slice (1D-projection) on variable iVar while axis ranges are defined with varMin,varMax
787   // arrays varMin and varMax contain the min and max values of each variable.
788   // therefore varMin and varMax must have their dimensions equal to GetNVar()
789   // If useBins=true, varMin and varMax are taken as bin numbers
790   
791   THnSparse* clone = (THnSparse*)fData->Clone();
792   for (Int_t iAxis=0; iAxis<GetNVar(); iAxis++) {
793     if (useBins)  clone->GetAxis(iAxis)->SetRange((Int_t)varMin[iAxis],(Int_t)varMax[iAxis]);
794     else          clone->GetAxis(iAxis)->SetRangeUser(varMin[iAxis],varMax[iAxis]);
795   }
796   return clone->Projection(iVar);
797 }
798
799 //____________________________________________________________________
800 TH2D* AliCFGridSparse::Slice(Int_t iVar1, Int_t iVar2, const Double_t *varMin, const Double_t *varMax, Bool_t useBins) const
801 {
802   //
803   // return a slice (2D-projection) on variables iVar1 and iVar2 while axis ranges are defined with varMin,varMax
804   // arrays varMin and varMax contain the min and max values of each variable.
805   // therefore varMin and varMax must have their dimensions equal to GetNVar()
806   // If useBins=true, varMin and varMax are taken as bin numbers
807   
808   THnSparse* clone = (THnSparse*)fData->Clone();
809   for (Int_t iAxis=0; iAxis<GetNVar(); iAxis++) {
810     if (useBins)  clone->GetAxis(iAxis)->SetRange((Int_t)varMin[iAxis],(Int_t)varMax[iAxis]);
811     else          clone->GetAxis(iAxis)->SetRangeUser(varMin[iAxis],varMax[iAxis]);
812   }
813   return clone->Projection(iVar1,iVar2);
814 }
815
816 //____________________________________________________________________
817 TH3D* AliCFGridSparse::Slice(Int_t iVar1, Int_t iVar2, Int_t iVar3, const Double_t *varMin, const Double_t *varMax, Bool_t useBins) const
818 {
819   //
820   // return a slice (3D-projection) on variables iVar1, iVar2 and iVar3 while axis ranges are defined with varMin,varMax
821   // arrays varMin and varMax contain the min and max values of each variable.
822   // therefore varMin and varMax must have their dimensions equal to GetNVar()
823   // If useBins=true, varMin and varMax are taken as bin numbers
824
825   THnSparse* clone = (THnSparse*)fData->Clone();
826   for (Int_t iAxis=0; iAxis<GetNVar(); iAxis++) {
827     if (useBins)  clone->GetAxis(iAxis)->SetRange((Int_t)varMin[iAxis],(Int_t)varMax[iAxis]);
828     else          clone->GetAxis(iAxis)->SetRangeUser(varMin[iAxis],varMax[iAxis]);
829   }
830   return clone->Projection(iVar1,iVar2,iVar3);
831 }
832
833 //____________________________________________________________________
834 void AliCFGridSparse::SetRangeUser(Int_t iVar, Double_t varMin, Double_t varMax) {
835   //
836   // set range of axis iVar. 
837   //
838   fData->GetAxis(iVar)->SetRangeUser(varMin,varMax);
839   AliWarning(Form("THnSparse axis %d range has been modified",iVar));
840 }
841
842 //____________________________________________________________________
843 void AliCFGridSparse::SetRangeUser(const Double_t *varMin, const Double_t *varMax) {
844   //
845   // set range of every axis. varMin and varMax must be of dimension GetNVar()
846   //
847   for (Int_t iAxis=0; iAxis<GetNVar() ; iAxis++) { // set new range for every axis
848     SetRangeUser(iAxis,varMin[iAxis],varMax[iAxis]);
849   }
850   AliWarning("THnSparse axes ranges have been modified");
851 }
852
853 //____________________________________________________________________
854 void AliCFGridSparse::UseAxisRange(Bool_t b) const {
855   for (Int_t iAxis=0; iAxis<GetNVar(); iAxis++) fData->GetAxis(iAxis)->SetBit(TAxis::kAxisRange,b);
856 }
857
858 //____________________________________________________________________
859 Float_t AliCFGridSparse::GetOverFlows(Int_t ivar, Bool_t exclusive) const
860 {
861   //
862   // Returns overflows in variable ivar
863   // Set 'exclusive' to true for an exclusive check on variable ivar
864   //
865   Int_t* bin = new Int_t[GetNVar()];
866   memset(bin, 0, sizeof(Int_t) * GetNVar());
867   Float_t ovfl=0.;
868   for (Long64_t i = 0; i < fData->GetNbins(); i++) {
869     Double_t v = fData->GetBinContent(i, bin);
870     Bool_t add=kTRUE;
871     if (exclusive) {
872       for(Int_t j=0;j<GetNVar();j++){
873         if(ivar==j)continue;
874         if((bin[j]==0) || (bin[j]==GetNBins(j)+1))add=kFALSE;
875       }
876     }
877     if(bin[ivar]==GetNBins(ivar)+1 && add) ovfl+=v;
878   }
879
880   delete[] bin;
881   return ovfl;
882 }
883
884 //____________________________________________________________________
885 Float_t AliCFGridSparse::GetUnderFlows(Int_t ivar, Bool_t exclusive) const
886 {
887   //
888   // Returns exclusive overflows in variable ivar
889   // Set 'exclusive' to true for an exclusive check on variable ivar
890   //
891   Int_t* bin = new Int_t[GetNVar()];
892   memset(bin, 0, sizeof(Int_t) * GetNVar());
893   Float_t unfl=0.;
894   for (Long64_t i = 0; i < fData->GetNbins(); i++) {
895     Double_t v = fData->GetBinContent(i, bin);
896     Bool_t add=kTRUE;
897     if (exclusive) {
898       for(Int_t j=0;j<GetNVar();j++){
899         if(ivar==j)continue;
900         if((bin[j]==0) || (bin[j]==GetNBins(j)+1))add=kFALSE;
901       }
902     }
903     if(bin[ivar]==0 && add) unfl+=v;
904   }
905
906   delete[] bin;
907   return unfl;
908 }
909