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 // projects the grid on the nVars dimensions defined in vars.
172 // axis ranges can be defined in arrays varMin, varMax
175 // binning for new grid
176 Int_t* bins = new Int_t[nVars];
177 for (Int_t iVar=0; iVar<nVars; iVar++) {
178 bins[iVar] = fNVarBins[vars[iVar]];
181 // create new grid sparse
182 AliCFGridSparse* out = new AliCFGridSparse(fName,fTitle,nVars,bins);
184 //set the range in the THnSparse to project
185 THnSparse* clone = ((THnSparse*)fData->Clone());
186 if (varMin && varMax) {
187 for (Int_t iAxis=0; iAxis<fNVar; iAxis++) {
188 clone->GetAxis(iAxis)->SetRangeUser(varMin[iAxis],varMax[iAxis]);
191 else AliInfo("Keeping same axis ranges");
193 out->SetGrid(clone->Projection(nVars,vars));
197 //____________________________________________________________________
198 Float_t AliCFGridSparse::GetOverFlows(Int_t ivar) const
201 // Returns exclusive overflows in variable ivar
203 Int_t* bin = new Int_t[fNVar];
204 memset(bin, 0, sizeof(Int_t) * fNVar);
206 for (Long64_t i = 0; i < fData->GetNbins(); ++i) {
207 Double_t v = fData->GetBinContent(i, bin);
209 for(Int_t j=0;j<fNVar;j++){
211 if((bin[j]==0) || (bin[j]==fNVarBins[j]+1))add=kFALSE;
213 if(bin[ivar]==fNVarBins[ivar]+1 && add) ovfl+=v;
220 //____________________________________________________________________
221 Float_t AliCFGridSparse::GetUnderFlows(Int_t ivar) const
224 // Returns exclusive overflows in variable ivar
226 Int_t* bin = new Int_t[fNVar];
227 memset(bin, 0, sizeof(Int_t) * fNVar);
229 for (Long64_t i = 0; i < fData->GetNbins(); ++i) {
230 Double_t v = fData->GetBinContent(i, bin);
232 for(Int_t j=0;j<fNVar;j++){
234 if((bin[j]==0) || (bin[j]==fNVarBins[j]+1))add=kFALSE;
236 if(bin[ivar]==0 && add) unfl+=v;
244 //____________________________________________________________________
245 Float_t AliCFGridSparse::GetEntries() const
248 // total entries (including overflows and underflows)
251 return fData->GetEntries();
254 //____________________________________________________________________
255 Float_t AliCFGridSparse::GetElement(Int_t index) const
258 // Returns content of grid element index according to the
259 // linear indexing in AliCFFrame
261 Int_t *bin = new Int_t[fNVar];
262 GetBinIndex(index, bin);
263 for(Int_t i=0;i<fNVar;i++)fIndex[i]=bin[i]+1; //consistency with AliCFGrid
264 Float_t val=GetElement(fIndex);
269 //____________________________________________________________________
270 Float_t AliCFGridSparse::GetElement(Int_t *bin) const
273 // Get the content in a bin corresponding to a set of bin indexes
275 return fData->GetBinContent(bin);
278 //____________________________________________________________________
279 Float_t AliCFGridSparse::GetElement(Double_t *var) const
282 // Get the content in a bin corresponding to a set of input variables
285 Long_t index=fData->GetBin(var,kFALSE); //this is the THnSparse index (do not allocate new cells if content is empty)
289 return fData->GetBinContent(index);
293 //____________________________________________________________________
294 Float_t AliCFGridSparse::GetElementError(Int_t index) const
297 // Returns the error on the content of a bin according to a linear
298 // indexing in AliCFFrame
301 Int_t *bin = new Int_t[fNVar];
302 GetBinIndex(index, bin);
303 for(Int_t i=0;i<fNVar;i++)fIndex[i]=bin[i]+1; //consistency with AliCFGrid
304 Float_t val=GetElementError(fIndex);
309 //____________________________________________________________________
310 Float_t AliCFGridSparse::GetElementError(Int_t *bin) const
313 // Get the error in a bin corresponding to a set of bin indexes
315 return fData->GetBinError(bin);
318 //____________________________________________________________________
319 Float_t AliCFGridSparse::GetElementError(Double_t *var) const
322 // Get the error in a bin corresponding to a set of input variables
325 Long_t index=fData->GetBin(var,kFALSE); //this is the THnSparse index (do not allocate new cells if content is empy)
329 return fData->GetBinError(index);
334 //____________________________________________________________________
335 void AliCFGridSparse::SetElement(Int_t index, Float_t val)
338 // Sets grid element iel to val (linear indexing) in AliCFFrame
340 Int_t *bin = new Int_t[fNVar];
341 GetBinIndex(index, bin);
342 for(Int_t i=0;i<fNVar;i++)fIndex[i]=bin[i]+1;
343 SetElement(fIndex,val);
346 //____________________________________________________________________
347 void AliCFGridSparse::SetElement(Int_t *bin, Float_t val)
350 // Sets grid element of bin indeces bin to val
352 fData->SetBinContent(bin,val);
354 //____________________________________________________________________
355 void AliCFGridSparse::SetElement(Double_t *var, Float_t val)
358 // Set the content in a bin to value val corresponding to a set of input variables
360 Long_t index=fData->GetBin(var); //THnSparse index: allocate the cell
361 Int_t *bin = new Int_t[fNVar];
362 fData->GetBinContent(index,bin); //trick to access the array of bins
363 fData->SetBinContent(bin,val);
368 //____________________________________________________________________
369 void AliCFGridSparse::SetElementError(Int_t index, Float_t val)
372 // Sets grid element iel error to val (linear indexing) in AliCFFrame
374 Int_t *bin = new Int_t[fNVar];
375 GetBinIndex(index, bin);
376 for(Int_t i=0;i<fNVar;i++)fIndex[i]=bin[i]+1;
377 SetElementError(fIndex,val);
380 //____________________________________________________________________
381 void AliCFGridSparse::SetElementError(Int_t *bin, Float_t val)
384 // Sets grid element error of bin indeces bin to val
386 fData->SetBinError(bin,val);
388 //____________________________________________________________________
389 void AliCFGridSparse::SetElementError(Double_t *var, Float_t val)
392 // Set the error in a bin to value val corresponding to a set of input variables
394 Long_t index=fData->GetBin(var); //THnSparse index
395 Int_t *bin = new Int_t[fNVar];
396 fData->GetBinContent(index,bin); //trick to access the array of bins
397 fData->SetBinError(bin,val);
401 //____________________________________________________________________
402 void AliCFGridSparse::SumW2()
405 //set calculation of the squared sum of the weighted entries
408 fData->CalculateErrors(kTRUE);
414 //____________________________________________________________________
415 void AliCFGridSparse::Add(AliCFVGrid* aGrid, Double_t c)
418 //add aGrid to the current one
422 if(aGrid->GetNVar()!=fNVar){
423 AliInfo("Different number of variables, cannot add the grids");
426 if(aGrid->GetNDim()!=fNDim){
427 AliInfo("Different number of dimensions, cannot add the grids!");
431 if(!fSumW2 && aGrid->GetSumW2())SumW2();
434 fData->Add(((AliCFGridSparse*)aGrid)->GetGrid(),c);
438 //____________________________________________________________________
439 void AliCFGridSparse::Add(AliCFVGrid* aGrid1, AliCFVGrid* aGrid2, Double_t c1,Double_t c2)
442 //Add aGrid1 and aGrid2 and deposit the result into the current one
445 if(fNVar!=aGrid1->GetNVar()|| fNVar!=aGrid2->GetNVar()){
446 AliInfo("Different number of variables, cannot add the grids");
449 if(fNDim!=aGrid1->GetNDim()|| fNDim!=aGrid2->GetNDim()){
450 AliInfo("Different number of dimensions, cannot add the grids!");
454 if(!fSumW2 && (aGrid1->GetSumW2() || aGrid2->GetSumW2()))SumW2();
458 fData->Add(((AliCFGridSparse*)aGrid1)->GetGrid(),c1);
459 fData->Add(((AliCFGridSparse*)aGrid2)->GetGrid(),c2);
463 //____________________________________________________________________
464 void AliCFGridSparse::Multiply(AliCFVGrid* aGrid, Double_t c)
467 // Multiply aGrid to the current one
471 if(aGrid->GetNVar()!=fNVar){
472 AliInfo("Different number of variables, cannot multiply the grids");
475 if(aGrid->GetNDim()!=fNDim){
476 AliInfo("Different number of dimensions, cannot multiply the grids!");
480 if(!fSumW2 && aGrid->GetSumW2())SumW2();
482 THnSparse *h= ((AliCFGridSparse*)aGrid)->GetGrid();
489 //____________________________________________________________________
490 void AliCFGridSparse::Multiply(AliCFVGrid* aGrid1, AliCFVGrid* aGrid2, Double_t c1,Double_t c2)
493 //Multiply aGrid1 and aGrid2 and deposit the result into the current one
496 if(fNVar!=aGrid1->GetNVar()|| fNVar!=aGrid2->GetNVar()){
497 AliInfo("Different number of variables, cannot multiply the grids");
500 if(fNDim!=aGrid1->GetNDim()|| fNDim!=aGrid2->GetNDim()){
501 AliInfo("Different number of dimensions, cannot multiply the grids!");
505 if(!fSumW2 && (aGrid1->GetSumW2() || aGrid2->GetSumW2()))SumW2();
509 THnSparse *h1= ((AliCFGridSparse*)aGrid1)->GetGrid();
510 THnSparse *h2= ((AliCFGridSparse*)aGrid2)->GetGrid();
518 //____________________________________________________________________
519 void AliCFGridSparse::Divide(AliCFVGrid* aGrid, Double_t c)
522 // Divide aGrid to the current one
526 if(aGrid->GetNVar()!=fNVar){
527 AliInfo("Different number of variables, cannot divide the grids");
530 if(aGrid->GetNDim()!=fNDim){
531 AliInfo("Different number of dimensions, cannot divide the grids!");
535 if(!fSumW2 && aGrid->GetSumW2())SumW2();
537 THnSparse *h= ((AliCFGridSparse*)aGrid)->GetGrid();
544 //____________________________________________________________________
545 void AliCFGridSparse::Divide(AliCFVGrid* aGrid1, AliCFVGrid* aGrid2, Double_t c1,Double_t c2, Option_t *option)
548 //Divide aGrid1 and aGrid2 and deposit the result into the current one
549 //bynomial errors are supported
552 if(fNVar!=aGrid1->GetNVar()|| fNVar!=aGrid2->GetNVar()){
553 AliInfo("Different number of variables, cannot divide the grids");
556 if(fNDim!=aGrid1->GetNDim()|| fNDim!=aGrid2->GetNDim()){
557 AliInfo("Different number of dimensions, cannot divide the grids!");
561 if(!fSumW2 && (aGrid1->GetSumW2() || aGrid2->GetSumW2()))SumW2();
564 THnSparse *h1= ((AliCFGridSparse*)aGrid1)->GetGrid();
565 THnSparse *h2= ((AliCFGridSparse*)aGrid2)->GetGrid();
566 fData->Divide(h1,h2,c1,c2,option);
570 //____________________________________________________________________
571 void AliCFGridSparse::Rebin(const Int_t* group)
574 // rebin the grid according to Rebin() as in THnSparse
575 // Please notice that the original number of bins on
576 // a given axis has to be divisible by the rebin group.
579 for(Int_t i=0;i<fNVar;i++){
580 if(group[i]!=1)AliInfo(Form(" merging bins along dimension %i in groups of %i bins", i,group[i]));
583 THnSparse *rebinned =fData->Rebin(group);
587 //redefine the needed stuff
592 //number of bins in each dimension, auxiliary variables
594 for(Int_t ivar=0;ivar<fNVar;ivar++){
595 Int_t nbins = fData->GetAxis(ivar)->GetNbins();
596 fNVarBins[ivar]=nbins;
597 ndimTot*=fNVarBins[ivar];
598 nbinTot+=(fNVarBins[ivar]+1);
600 for(Int_t i =0;i<ivar;i++)offset+=(fNVarBins[i]+1);
601 fOffset[ivar]=offset;
603 for(Int_t i=0;i<ivar;i++)prod*=fNVarBins[i];
609 //now the array of bin limits
611 delete fVarBinLimits;
612 fNVarBinLimits=nbinTot;
613 fVarBinLimits=new Double_t[fNVarBinLimits];
615 for(Int_t ivar=0;ivar<fNVar;ivar++){
616 Double_t low = fData->GetAxis(ivar)->GetXmin();
617 Double_t high = fData->GetAxis(ivar)->GetXmax();
618 const TArrayD *xbins = fData->GetAxis(ivar)->GetXbins();
620 for(Int_t ibin=0;ibin<=fNVarBins[ivar];ibin++){
621 fVarBinLimits[ibin+fOffset[ivar]] = low + ibin*(high-low)/((Double_t) fNVarBins[ivar]);
626 for(Int_t ibin=0;ibin<=fNVarBins[ivar];ibin++) {
627 fVarBinLimits[ibin+fOffset[ivar]] = xbins->At(ibin);
633 //____________________________________________________________________
634 void AliCFGridSparse::Copy(TObject& c) const
639 AliCFGridSparse& target = (AliCFGridSparse &) c;
641 if(fData)target.fData = fData;
644 //____________________________________________________________________
645 TH1D* AliCFGridSparse::Slice(Int_t iVar, Double_t *varMin, Double_t *varMax) const
648 // return a slice (1D-projection) on variable iVar while axis ranges are defined with varMin,varMax
649 // arrays varMin and varMax contain the min and max values of each variable.
650 // therefore varMin and varMax must have their dimensions equal to fNVar
653 THnSparse* clone = (THnSparse*)fData->Clone();
654 for (Int_t iAxis=0; iAxis<fNVar; iAxis++) {
655 clone->GetAxis(iAxis)->SetRangeUser(varMin[iAxis],varMax[iAxis]);
657 return clone->Projection(iVar);
660 //____________________________________________________________________
661 TH2D* AliCFGridSparse::Slice(Int_t iVar1, Int_t iVar2, Double_t *varMin, Double_t *varMax) const
664 // return a slice (2D-projection) on variables iVar1 and iVar2 while axis ranges are defined with varMin,varMax
665 // arrays varMin and varMax contain the min and max values of each variable.
666 // therefore varMin and varMax must have their dimensions equal to fNVar
669 THnSparse* clone = (THnSparse*)fData->Clone();
670 for (Int_t iAxis=0; iAxis<fNVar; iAxis++) {
671 clone->GetAxis(iAxis)->SetRangeUser(varMin[iAxis],varMax[iAxis]);
673 return clone->Projection(iVar1,iVar2);
676 //____________________________________________________________________
677 TH3D* AliCFGridSparse::Slice(Int_t iVar1, Int_t iVar2, Int_t iVar3, Double_t *varMin, Double_t *varMax) const
680 // return a slice (3D-projection) on variables iVar1, iVar2 and iVar3 while axis ranges are defined with varMin,varMax
681 // arrays varMin and varMax contain the min and max values of each variable.
682 // therefore varMin and varMax must have their dimensions equal to fNVar
685 THnSparse* clone = (THnSparse*)fData->Clone();
686 for (Int_t iAxis=0; iAxis<fNVar; iAxis++) {
687 clone->GetAxis(iAxis)->SetRangeUser(varMin[iAxis],varMax[iAxis]);
689 return clone->Projection(iVar1,iVar2,iVar3);
692 //____________________________________________________________________
693 void AliCFGridSparse::SetRangeUser(Int_t iVar, Double_t varMin, Double_t varMax) {
695 // set range of axis iVar.
697 fData->GetAxis(iVar)->SetRangeUser(varMin,varMax);
698 AliWarning(Form("THnSparse axis %d range has been modified",iVar));
701 //____________________________________________________________________
702 void AliCFGridSparse::SetRangeUser(Double_t *varMin, Double_t *varMax) {
704 // set range of every axis. varMin and varMax must be of dimension fNVar
706 for (Int_t iAxis=0; iAxis<fNVar ; iAxis++) { // set new range for every axis
707 SetRangeUser(iAxis,varMin[iAxis],varMax[iAxis]);
709 AliWarning("THnSparse axes ranges have been modified");
712 //____________________________________________________________________
713 void AliCFGridSparse::UseAxisRange(Bool_t b) const {
714 for (Int_t iAxis=0; iAxis<fNVar; iAxis++) fData->GetAxis(iAxis)->SetBit(TAxis::kAxisRange,b);