]> git.uio.no Git - u/mrichter/AliRoot.git/blobdiff - CORRFW/AliCFGridSparse.cxx
excluding D0Trigger from the inputs of the HLTGlobalTrigger because of a rare segfaul...
[u/mrichter/AliRoot.git] / CORRFW / AliCFGridSparse.cxx
index 3641f663cee1911dfc422f3e01a79abe5a2eceb6..1c4f31232032362146a26e7e4fd4fe5a2d299679 100755 (executable)
@@ -97,6 +97,19 @@ AliCFGridSparse& AliCFGridSparse::operator=(const AliCFGridSparse &c)
   return *this;
 } 
 
+//____________________________________________________________________
+void AliCFGridSparse::SetBinLimits(Int_t ivar, Double_t min, Double_t max)
+{
+  //
+  // set a uniform binning for variable ivar
+  //
+  Int_t nBins = GetNBins(ivar);
+  Double_t * array = new Double_t[nBins+1];
+  for (Int_t iEdge=0; iEdge<=nBins; iEdge++) array[iEdge] = min + iEdge * (max-min)/nBins ;
+  fData->SetBinEdges(ivar, array);
+  delete [] array ;
+} 
+
 //____________________________________________________________________
 void AliCFGridSparse::SetBinLimits(Int_t ivar, const Double_t *array)
 {
@@ -128,6 +141,10 @@ TH1D *AliCFGridSparse::Project(Int_t ivar) const
   hist->SetXTitle(GetVarTitle(ivar));
   hist->SetName(Form("%s_proj-%s",GetName(),GetVarTitle(ivar)));
   hist->SetTitle(Form("%s: projection on %s",GetTitle(),GetVarTitle(ivar)));
+  for (Int_t iBin=1; iBin<=GetNBins(ivar); iBin++) {
+    TString binLabel = GetAxis(ivar)->GetBinLabel(iBin) ;
+    if (binLabel.CompareTo("") != 0) hist->GetXaxis()->SetBinLabel(iBin,binLabel);
+  }
   return hist;
 }
 //___________________________________________________________________
@@ -142,6 +159,14 @@ TH2D *AliCFGridSparse::Project(Int_t ivar1, Int_t ivar2) const
   hist->SetYTitle(GetVarTitle(ivar2));
   hist->SetName(Form("%s_proj-%s,%s",GetName(),GetVarTitle(ivar1),GetVarTitle(ivar2)));
   hist->SetTitle(Form("%s: projection on %s-%s",GetTitle(),GetVarTitle(ivar1),GetVarTitle(ivar2)));
+  for (Int_t iBin=1; iBin<=GetNBins(ivar1); iBin++) {
+    TString binLabel = GetAxis(ivar1)->GetBinLabel(iBin) ;
+    if (binLabel.CompareTo("") != 0) hist->GetXaxis()->SetBinLabel(iBin,binLabel);
+  }
+  for (Int_t iBin=1; iBin<=GetNBins(ivar2); iBin++) {
+    TString binLabel = GetAxis(ivar2)->GetBinLabel(iBin) ;
+    if (binLabel.CompareTo("") != 0) hist->GetYaxis()->SetBinLabel(iBin,binLabel);
+  }
   return hist;
 }
 //___________________________________________________________________
@@ -157,15 +182,28 @@ TH3D *AliCFGridSparse::Project(Int_t ivar1, Int_t ivar2, Int_t ivar3) const
   hist->SetZTitle(GetVarTitle(ivar3));
   hist->SetName(Form("%s_proj-%s,%s,%s",GetName(),GetVarTitle(ivar1),GetVarTitle(ivar2),GetVarTitle(ivar3)));
   hist->SetTitle(Form("%s: projection on %s-%s-%s",GetTitle(),GetVarTitle(ivar1),GetVarTitle(ivar2),GetVarTitle(ivar3)));
+  for (Int_t iBin=1; iBin<=GetNBins(ivar1); iBin++) {
+    TString binLabel = GetAxis(ivar1)->GetBinLabel(iBin) ;
+    if (binLabel.CompareTo("") != 0) hist->GetXaxis()->SetBinLabel(iBin,binLabel);
+  }
+  for (Int_t iBin=1; iBin<=GetNBins(ivar2); iBin++) {
+    TString binLabel = GetAxis(ivar2)->GetBinLabel(iBin) ;
+    if (binLabel.CompareTo("") != 0) hist->GetYaxis()->SetBinLabel(iBin,binLabel);
+  }
+  for (Int_t iBin=1; iBin<=GetNBins(ivar3); iBin++) {
+    TString binLabel = GetAxis(ivar3)->GetBinLabel(iBin) ;
+    if (binLabel.CompareTo("") != 0) hist->GetZaxis()->SetBinLabel(iBin,binLabel);
+  }
   return hist;
 }
 
 //___________________________________________________________________
-AliCFGridSparse* AliCFGridSparse::Project(Int_t nVars, const Int_t* vars, const Double_t* varMin, const Double_t* varMax) const
+AliCFGridSparse* AliCFGridSparse::Project(Int_t nVars, const Int_t* vars, const Double_t* varMin, const Double_t* varMax, Bool_t useBins) const
 {
   //
   // projects the grid on the nVars dimensions defined in vars.
   // axis ranges can be defined in arrays varMin, varMax
+  // If useBins=true, varMin and varMax are taken as bin numbers
   //
 
   // binning for new grid
@@ -181,12 +219,14 @@ AliCFGridSparse* AliCFGridSparse::Project(Int_t nVars, const Int_t* vars, const
   THnSparse* clone = ((THnSparse*)fData->Clone());
   if (varMin && varMax) {
     for (Int_t iAxis=0; iAxis<GetNVar(); iAxis++) {
-      clone->GetAxis(iAxis)->SetRangeUser(varMin[iAxis],varMax[iAxis]);
+      if (useBins)  clone->GetAxis(iAxis)->SetRange((Int_t)varMin[iAxis],(Int_t)varMax[iAxis]);
+      else          clone->GetAxis(iAxis)->SetRangeUser(varMin[iAxis],varMax[iAxis]);
     }
   }
   else AliInfo("Keeping same axis ranges");
 
   out->SetGrid(clone->Projection(nVars,vars));
+  delete [] bins;
   return out;
 }
 
@@ -764,51 +804,111 @@ void AliCFGridSparse::Copy(TObject& c) const
 }
 
 //____________________________________________________________________
-TH1D* AliCFGridSparse::Slice(Int_t iVar, const Double_t *varMin, const Double_t *varMax) const
+TH1D* AliCFGridSparse::Slice(Int_t iVar, const Double_t *varMin, const Double_t *varMax, Bool_t useBins) const
 {
   //
   // return a slice (1D-projection) on variable iVar while axis ranges are defined with varMin,varMax
   // arrays varMin and varMax contain the min and max values of each variable.
   // therefore varMin and varMax must have their dimensions equal to GetNVar()
-  //
+  // If useBins=true, varMin and varMax are taken as bin numbers
   
   THnSparse* clone = (THnSparse*)fData->Clone();
   for (Int_t iAxis=0; iAxis<GetNVar(); iAxis++) {
-    clone->GetAxis(iAxis)->SetRangeUser(varMin[iAxis],varMax[iAxis]);
+    if (useBins)  clone->GetAxis(iAxis)->SetRange((Int_t)varMin[iAxis],(Int_t)varMax[iAxis]);
+    else          clone->GetAxis(iAxis)->SetRangeUser(varMin[iAxis],varMax[iAxis]);
   }
-  return clone->Projection(iVar);
+
+  TH1D* projection = 0x0 ;
+
+  TObject* objInMem = gROOT->FindObject(Form("%s_proj_%d",clone->GetName(),iVar)) ;
+  if (objInMem) {
+    TString name(objInMem->ClassName());
+    if (name.CompareTo("TH1D")==0) projection = (TH1D*) objInMem ;
+    else projection = clone->Projection(iVar);
+  }
+  else projection = clone->Projection(iVar); 
+  delete clone;
+  for (Int_t iBin=1; iBin<=GetNBins(iVar); iBin++) {
+    TString binLabel = GetAxis(iVar)->GetBinLabel(iBin) ;
+    if (binLabel.CompareTo("") != 0) projection->GetXaxis()->SetBinLabel(iBin,binLabel);
+  }
+  return projection ;
 }
 
 //____________________________________________________________________
-TH2D* AliCFGridSparse::Slice(Int_t iVar1, Int_t iVar2, const Double_t *varMin, const Double_t *varMax) const
+TH2D* AliCFGridSparse::Slice(Int_t iVar1, Int_t iVar2, const Double_t *varMin, const Double_t *varMax, Bool_t useBins) const
 {
   //
   // return a slice (2D-projection) on variables iVar1 and iVar2 while axis ranges are defined with varMin,varMax
   // arrays varMin and varMax contain the min and max values of each variable.
   // therefore varMin and varMax must have their dimensions equal to GetNVar()
-  //
+  // If useBins=true, varMin and varMax are taken as bin numbers
   
   THnSparse* clone = (THnSparse*)fData->Clone();
   for (Int_t iAxis=0; iAxis<GetNVar(); iAxis++) {
-    clone->GetAxis(iAxis)->SetRangeUser(varMin[iAxis],varMax[iAxis]);
+    if (useBins)  clone->GetAxis(iAxis)->SetRange((Int_t)varMin[iAxis],(Int_t)varMax[iAxis]);
+    else          clone->GetAxis(iAxis)->SetRangeUser(varMin[iAxis],varMax[iAxis]);
+  }
+
+  TH2D* projection = 0x0 ;
+
+  TObject* objInMem = gROOT->FindObject(Form("%s_proj_%d_%d",clone->GetName(),iVar2,iVar1)) ;
+  if (objInMem) {
+    TString name(objInMem->ClassName());
+    if (name.CompareTo("TH2D")==0) projection = (TH2D*) objInMem ;
+    else projection = clone->Projection(iVar1,iVar2);
+  }
+  else projection = clone->Projection(iVar1,iVar2); 
+  delete clone;
+  for (Int_t iBin=1; iBin<=GetNBins(iVar1); iBin++) {
+    TString binLabel = GetAxis(iVar1)->GetBinLabel(iBin) ;
+    if (binLabel.CompareTo("") != 0) projection->GetXaxis()->SetBinLabel(iBin,binLabel);
+  }
+  for (Int_t iBin=1; iBin<=GetNBins(iVar2); iBin++) {
+    TString binLabel = GetAxis(iVar2)->GetBinLabel(iBin) ;
+    if (binLabel.CompareTo("") != 0) projection->GetYaxis()->SetBinLabel(iBin,binLabel);
   }
-  return clone->Projection(iVar1,iVar2);
+  return projection ;
 }
 
 //____________________________________________________________________
-TH3D* AliCFGridSparse::Slice(Int_t iVar1, Int_t iVar2, Int_t iVar3, const Double_t *varMin, const Double_t *varMax) const
+TH3D* AliCFGridSparse::Slice(Int_t iVar1, Int_t iVar2, Int_t iVar3, const Double_t *varMin, const Double_t *varMax, Bool_t useBins) const
 {
   //
   // return a slice (3D-projection) on variables iVar1, iVar2 and iVar3 while axis ranges are defined with varMin,varMax
   // arrays varMin and varMax contain the min and max values of each variable.
   // therefore varMin and varMax must have their dimensions equal to GetNVar()
-  //
+  // If useBins=true, varMin and varMax are taken as bin numbers
 
   THnSparse* clone = (THnSparse*)fData->Clone();
   for (Int_t iAxis=0; iAxis<GetNVar(); iAxis++) {
-    clone->GetAxis(iAxis)->SetRangeUser(varMin[iAxis],varMax[iAxis]);
+    if (useBins)  clone->GetAxis(iAxis)->SetRange((Int_t)varMin[iAxis],(Int_t)varMax[iAxis]);
+    else          clone->GetAxis(iAxis)->SetRangeUser(varMin[iAxis],varMax[iAxis]);
+  }
+
+  TH3D* projection = 0x0 ;
+
+  TObject* objInMem = gROOT->FindObject(Form("%s_proj_%d_%d_%d",clone->GetName(),iVar1,iVar2,iVar3)) ;
+  if (objInMem) {
+    TString name(objInMem->ClassName());
+    if (name.CompareTo("TH3D")==0) projection = (TH3D*) objInMem ;
+    else projection = clone->Projection(iVar1,iVar2,iVar3);
+  }
+  else projection = clone->Projection(iVar1,iVar2,iVar3); 
+  delete clone;
+  for (Int_t iBin=1; iBin<=GetNBins(iVar1); iBin++) {
+    TString binLabel = GetAxis(iVar1)->GetBinLabel(iBin) ;
+    if (binLabel.CompareTo("") != 0) projection->GetXaxis()->SetBinLabel(iBin,binLabel);
+  }
+  for (Int_t iBin=1; iBin<=GetNBins(iVar2); iBin++) {
+    TString binLabel = GetAxis(iVar2)->GetBinLabel(iBin) ;
+    if (binLabel.CompareTo("") != 0) projection->GetYaxis()->SetBinLabel(iBin,binLabel);
+  }
+  for (Int_t iBin=1; iBin<=GetNBins(iVar3); iBin++) {
+    TString binLabel = GetAxis(iVar3)->GetBinLabel(iBin) ;
+    if (binLabel.CompareTo("") != 0) projection->GetZaxis()->SetBinLabel(iBin,binLabel);
   }
-  return clone->Projection(iVar1,iVar2,iVar3);
+  return projection ;
 }
 
 //____________________________________________________________________