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