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() :
45 // default constructor
47 //____________________________________________________________________
48 AliCFGridSparse::AliCFGridSparse(const Char_t* name, const Char_t* title) :
49 AliCFVGrid(name,title),
52 // default constructor
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),
63 fData=new THnSparseF(name,title,fNVar,fNVarBins);
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];
72 fData->SetBinEdges(ivar, array);
77 //____________________________________________________________________
78 AliCFGridSparse::AliCFGridSparse(const AliCFGridSparse& c) :
85 ((AliCFGridSparse &)c).Copy(*this);
88 //____________________________________________________________________
89 AliCFGridSparse::~AliCFGridSparse()
94 if(fData) delete fData;
97 //____________________________________________________________________
98 AliCFGridSparse &AliCFGridSparse::operator=(const AliCFGridSparse &c)
101 // assigment operator
104 ((AliCFGridSparse &) c).Copy(*this);
108 //____________________________________________________________________
109 void AliCFGridSparse::SetBinLimits(Int_t ivar, Double_t *array)
112 // setting the arrays containing the bin limits
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];
123 //____________________________________________________________________
124 void AliCFGridSparse::Fill(Double_t *var, Double_t weight)
128 // given a set of values of the input variable,
129 // with weight (by default w=1)
131 fData->Fill(var,weight);
134 //___________________________________________________________________
135 TH1D *AliCFGridSparse::Project(Int_t ivar) const
138 // Make a 1D projection along variable ivar
141 TH1D *hist=fData->Projection(ivar);
144 //___________________________________________________________________
145 TH2D *AliCFGridSparse::Project(Int_t ivar1, Int_t ivar2) const
148 // Make a 2D projection along variables ivar1 & ivar2
151 TH2D *hist=fData->Projection(ivar2,ivar1); //notice inverted axis (THnSparse uses TH3 2d-projection convention...)
155 //___________________________________________________________________
156 TH3D *AliCFGridSparse::Project(Int_t ivar1, Int_t ivar2, Int_t ivar3) const
159 // Make a 3D projection along variables ivar1 & ivar2 & ivar3
162 TH3D *hist=fData->Projection(ivar1,ivar2,ivar3);
167 //___________________________________________________________________
168 AliCFGridSparse* AliCFGridSparse::Project(Int_t nVars, Int_t* vars, Double_t* varMin, Double_t* varMax) const
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]];
177 // create new grid sparse
178 AliCFGridSparse* out = new AliCFGridSparse(fName,fTitle,nVars,bins);
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]);
185 out->SetGrid(clone->Projection(nVars,vars));
189 //____________________________________________________________________
190 Float_t AliCFGridSparse::GetOverFlows(Int_t ivar) const
193 // Returns exclusive overflows in variable ivar
195 Int_t* bin = new Int_t[fNVar];
196 memset(bin, 0, sizeof(Int_t) * fNVar);
198 for (Long64_t i = 0; i < fData->GetNbins(); ++i) {
199 Double_t v = fData->GetBinContent(i, bin);
201 for(Int_t j=0;j<fNVar;j++){
203 if((bin[j]==0) || (bin[j]==fNVarBins[j]+1))add=kFALSE;
205 if(bin[ivar]==fNVarBins[ivar]+1 && add) ovfl+=v;
212 //____________________________________________________________________
213 Float_t AliCFGridSparse::GetUnderFlows(Int_t ivar) const
216 // Returns exclusive overflows in variable ivar
218 Int_t* bin = new Int_t[fNVar];
219 memset(bin, 0, sizeof(Int_t) * fNVar);
221 for (Long64_t i = 0; i < fData->GetNbins(); ++i) {
222 Double_t v = fData->GetBinContent(i, bin);
224 for(Int_t j=0;j<fNVar;j++){
226 if((bin[j]==0) || (bin[j]==fNVarBins[j]+1))add=kFALSE;
228 if(bin[ivar]==0 && add) unfl+=v;
236 //____________________________________________________________________
237 Float_t AliCFGridSparse::GetEntries() const
240 // total entries (including overflows and underflows)
243 return fData->GetEntries();
246 //____________________________________________________________________
247 Float_t AliCFGridSparse::GetElement(Int_t index) const
250 // Returns content of grid element index according to the
251 // linear indexing in AliCFFrame
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);
261 //____________________________________________________________________
262 Float_t AliCFGridSparse::GetElement(Int_t *bin) const
265 // Get the content in a bin corresponding to a set of bin indexes
267 return fData->GetBinContent(bin);
270 //____________________________________________________________________
271 Float_t AliCFGridSparse::GetElement(Double_t *var) const
274 // Get the content in a bin corresponding to a set of input variables
277 Long_t index=fData->GetBin(var,kFALSE); //this is the THnSparse index (do not allocate new cells if content is empty)
281 return fData->GetBinContent(index);
285 //____________________________________________________________________
286 Float_t AliCFGridSparse::GetElementError(Int_t index) const
289 // Returns the error on the content of a bin according to a linear
290 // indexing in AliCFFrame
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);
301 //____________________________________________________________________
302 Float_t AliCFGridSparse::GetElementError(Int_t *bin) const
305 // Get the error in a bin corresponding to a set of bin indexes
307 return fData->GetBinError(bin);
310 //____________________________________________________________________
311 Float_t AliCFGridSparse::GetElementError(Double_t *var) const
314 // Get the error in a bin corresponding to a set of input variables
317 Long_t index=fData->GetBin(var,kFALSE); //this is the THnSparse index (do not allocate new cells if content is empy)
321 return fData->GetBinError(index);
326 //____________________________________________________________________
327 void AliCFGridSparse::SetElement(Int_t index, Float_t val)
330 // Sets grid element iel to val (linear indexing) in AliCFFrame
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);
338 //____________________________________________________________________
339 void AliCFGridSparse::SetElement(Int_t *bin, Float_t val)
342 // Sets grid element of bin indeces bin to val
344 fData->SetBinContent(bin,val);
346 //____________________________________________________________________
347 void AliCFGridSparse::SetElement(Double_t *var, Float_t val)
350 // Set the content in a bin to value val corresponding to a set of input variables
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);
360 //____________________________________________________________________
361 void AliCFGridSparse::SetElementError(Int_t index, Float_t val)
364 // Sets grid element iel error to val (linear indexing) in AliCFFrame
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);
372 //____________________________________________________________________
373 void AliCFGridSparse::SetElementError(Int_t *bin, Float_t val)
376 // Sets grid element error of bin indeces bin to val
378 fData->SetBinError(bin,val);
380 //____________________________________________________________________
381 void AliCFGridSparse::SetElementError(Double_t *var, Float_t val)
384 // Set the error in a bin to value val corresponding to a set of input variables
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);
393 //____________________________________________________________________
394 void AliCFGridSparse::SumW2()
397 //set calculation of the squared sum of the weighted entries
400 fData->CalculateErrors(kTRUE);
406 //____________________________________________________________________
407 void AliCFGridSparse::Add(AliCFVGrid* aGrid, Double_t c)
410 //add aGrid to the current one
414 if(aGrid->GetNVar()!=fNVar){
415 AliInfo("Different number of variables, cannot add the grids");
418 if(aGrid->GetNDim()!=fNDim){
419 AliInfo("Different number of dimensions, cannot add the grids!");
423 if(!fSumW2 && aGrid->GetSumW2())SumW2();
426 fData->Add(((AliCFGridSparse*)aGrid)->GetGrid(),c);
430 //____________________________________________________________________
431 void AliCFGridSparse::Add(AliCFVGrid* aGrid1, AliCFVGrid* aGrid2, Double_t c1,Double_t c2)
434 //Add aGrid1 and aGrid2 and deposit the result into the current one
437 if(fNVar!=aGrid1->GetNVar()|| fNVar!=aGrid2->GetNVar()){
438 AliInfo("Different number of variables, cannot add the grids");
441 if(fNDim!=aGrid1->GetNDim()|| fNDim!=aGrid2->GetNDim()){
442 AliInfo("Different number of dimensions, cannot add the grids!");
446 if(!fSumW2 && (aGrid1->GetSumW2() || aGrid2->GetSumW2()))SumW2();
450 fData->Add(((AliCFGridSparse*)aGrid1)->GetGrid(),c1);
451 fData->Add(((AliCFGridSparse*)aGrid2)->GetGrid(),c2);
455 //____________________________________________________________________
456 void AliCFGridSparse::Multiply(AliCFVGrid* aGrid, Double_t c)
459 // Multiply aGrid to the current one
463 if(aGrid->GetNVar()!=fNVar){
464 AliInfo("Different number of variables, cannot multiply the grids");
467 if(aGrid->GetNDim()!=fNDim){
468 AliInfo("Different number of dimensions, cannot multiply the grids!");
472 if(!fSumW2 && aGrid->GetSumW2())SumW2();
474 THnSparse *h= ((AliCFGridSparse*)aGrid)->GetGrid();
481 //____________________________________________________________________
482 void AliCFGridSparse::Multiply(AliCFVGrid* aGrid1, AliCFVGrid* aGrid2, Double_t c1,Double_t c2)
485 //Multiply aGrid1 and aGrid2 and deposit the result into the current one
488 if(fNVar!=aGrid1->GetNVar()|| fNVar!=aGrid2->GetNVar()){
489 AliInfo("Different number of variables, cannot multiply the grids");
492 if(fNDim!=aGrid1->GetNDim()|| fNDim!=aGrid2->GetNDim()){
493 AliInfo("Different number of dimensions, cannot multiply the grids!");
497 if(!fSumW2 && (aGrid1->GetSumW2() || aGrid2->GetSumW2()))SumW2();
501 THnSparse *h1= ((AliCFGridSparse*)aGrid1)->GetGrid();
502 THnSparse *h2= ((AliCFGridSparse*)aGrid2)->GetGrid();
510 //____________________________________________________________________
511 void AliCFGridSparse::Divide(AliCFVGrid* aGrid, Double_t c)
514 // Divide aGrid to the current one
518 if(aGrid->GetNVar()!=fNVar){
519 AliInfo("Different number of variables, cannot divide the grids");
522 if(aGrid->GetNDim()!=fNDim){
523 AliInfo("Different number of dimensions, cannot divide the grids!");
527 if(!fSumW2 && aGrid->GetSumW2())SumW2();
529 THnSparse *h= ((AliCFGridSparse*)aGrid)->GetGrid();
536 //____________________________________________________________________
537 void AliCFGridSparse::Divide(AliCFVGrid* aGrid1, AliCFVGrid* aGrid2, Double_t c1,Double_t c2, Option_t *option)
540 //Divide aGrid1 and aGrid2 and deposit the result into the current one
541 //bynomial errors are supported
544 if(fNVar!=aGrid1->GetNVar()|| fNVar!=aGrid2->GetNVar()){
545 AliInfo("Different number of variables, cannot divide the grids");
548 if(fNDim!=aGrid1->GetNDim()|| fNDim!=aGrid2->GetNDim()){
549 AliInfo("Different number of dimensions, cannot divide the grids!");
553 if(!fSumW2 && (aGrid1->GetSumW2() || aGrid2->GetSumW2()))SumW2();
556 THnSparse *h1= ((AliCFGridSparse*)aGrid1)->GetGrid();
557 THnSparse *h2= ((AliCFGridSparse*)aGrid2)->GetGrid();
558 fData->Divide(h1,h2,c1,c2,option);
562 //____________________________________________________________________
563 void AliCFGridSparse::Rebin(const Int_t* group)
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.
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]));
575 THnSparse *rebinned =fData->Rebin(group);
579 //redefine the needed stuff
584 //number of bins in each dimension, auxiliary variables
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);
592 for(Int_t i =0;i<ivar;i++)offset+=(fNVarBins[i]+1);
593 fOffset[ivar]=offset;
595 for(Int_t i=0;i<ivar;i++)prod*=fNVarBins[i];
601 //now the array of bin limits
603 delete fVarBinLimits;
604 fNVarBinLimits=nbinTot;
605 fVarBinLimits=new Double_t[fNVarBinLimits];
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();
612 for(Int_t ibin=0;ibin<=fNVarBins[ivar];ibin++){
613 fVarBinLimits[ibin+fOffset[ivar]] = low + ibin*(high-low)/((Double_t) fNVarBins[ivar]);
618 for(Int_t ibin=0;ibin<=fNVarBins[ivar];ibin++) {
619 fVarBinLimits[ibin+fOffset[ivar]] = xbins->At(ibin);
625 //____________________________________________________________________
626 void AliCFGridSparse::Copy(TObject& c) const
631 AliCFGridSparse& target = (AliCFGridSparse &) c;
633 if(fData)target.fData = fData;
636 //____________________________________________________________________
637 TH1D* AliCFGridSparse::Slice(Int_t iVar, Double_t *varMin, Double_t *varMax) const
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
645 THnSparse* clone = (THnSparse*)fData->Clone();
646 for (Int_t iAxis=0; iAxis<fNVar; iAxis++) {
647 clone->GetAxis(iAxis)->SetRangeUser(varMin[iAxis],varMax[iAxis]);
649 return clone->Projection(iVar);
652 //____________________________________________________________________
653 TH2D* AliCFGridSparse::Slice(Int_t iVar1, Int_t iVar2, Double_t *varMin, Double_t *varMax) const
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
661 THnSparse* clone = (THnSparse*)fData->Clone();
662 for (Int_t iAxis=0; iAxis<fNVar; iAxis++) {
663 clone->GetAxis(iAxis)->SetRangeUser(varMin[iAxis],varMax[iAxis]);
665 return clone->Projection(iVar1,iVar2);
668 //____________________________________________________________________
669 TH3D* AliCFGridSparse::Slice(Int_t iVar1, Int_t iVar2, Int_t iVar3, Double_t *varMin, Double_t *varMax) const
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
677 THnSparse* clone = (THnSparse*)fData->Clone();
678 for (Int_t iAxis=0; iAxis<fNVar; iAxis++) {
679 clone->GetAxis(iAxis)->SetRangeUser(varMin[iAxis],varMax[iAxis]);
681 return clone->Projection(iVar1,iVar2,iVar3);
684 //____________________________________________________________________
685 void AliCFGridSparse::SetRangeUser(Int_t iVar, Double_t varMin, Double_t varMax) {
687 // set range of axis iVar.
689 fData->GetAxis(iVar)->SetRangeUser(varMin,varMax);
690 AliWarning(Form("THnSparse axis %d range has been modified",iVar));
693 //____________________________________________________________________
694 void AliCFGridSparse::SetRangeUser(Double_t *varMin, Double_t *varMax) {
696 // set range of every axis. varMin and varMax must be of dimension fNVar
698 for (Int_t iAxis=0; iAxis<fNVar ; iAxis++) { // set new range for every axis
699 SetRangeUser(iAxis,varMin[iAxis],varMax[iAxis]);
701 AliWarning("THnSparse axes ranges have been modified");
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);