#include <string.h>
#include "Riostream.h"
#include "TROOT.h"
#include "TGraphBentErrors.h"
#include "TStyle.h"
#include "TMath.h"
#include "TArrow.h"
#include "TBox.h"
#include "TVirtualPad.h"
#include "TH1.h"
#include "TF1.h"
ClassImp(TGraphBentErrors)
TGraphBentErrors::TGraphBentErrors(): TGraph()
{
CtorAllocate();
}
TGraphBentErrors::TGraphBentErrors(const TGraphBentErrors &gr)
: TGraph(gr)
{
if (!CtorAllocate()) return;
Int_t n = fNpoints*sizeof(Double_t);
memcpy(fEXlow, gr.fEXlow, n);
memcpy(fEYlow, gr.fEYlow, n);
memcpy(fEXhigh, gr.fEXhigh, n);
memcpy(fEYhigh, gr.fEYhigh, n);
memcpy(fEXlowd, gr.fEXlowd, n);
memcpy(fEYlowd, gr.fEYlowd, n);
memcpy(fEXhighd, gr.fEXhighd, n);
memcpy(fEYhighd, gr.fEYhighd, n);
}
TGraphBentErrors::TGraphBentErrors(Int_t n)
: TGraph(n)
{
if (!CtorAllocate()) return;
FillZero(0, fNpoints);
}
TGraphBentErrors::TGraphBentErrors(Int_t n,
const Float_t *x, const Float_t *y,
const Float_t *exl, const Float_t *exh,
const Float_t *eyl, const Float_t *eyh,
const Float_t *exld, const Float_t *exhd,
const Float_t *eyld, const Float_t *eyhd)
: TGraph(n,x,y)
{
if (!CtorAllocate()) return;
for (Int_t i=0;i<n;i++) {
if (exl) fEXlow[i] = exl[i];
else fEXlow[i] = 0;
if (exh) fEXhigh[i] = exh[i];
else fEXhigh[i] = 0;
if (eyl) fEYlow[i] = eyl[i];
else fEYlow[i] = 0;
if (eyh) fEYhigh[i] = eyh[i];
else fEYhigh[i] = 0;
if (exld) fEXlowd[i] = exld[i];
else fEXlowd[i] = 0;
if (exhd) fEXhighd[i] = exhd[i];
else fEXhighd[i] = 0;
if (eyld) fEYlowd[i] = eyld[i];
else fEYlowd[i] = 0;
if (eyhd) fEYhighd[i] = eyhd[i];
else fEYhighd[i] = 0;
}
}
TGraphBentErrors::TGraphBentErrors(Int_t n,
const Double_t *x, const Double_t *y,
const Double_t *exl, const Double_t *exh,
const Double_t *eyl, const Double_t *eyh,
const Double_t *exld, const Double_t *exhd,
const Double_t *eyld, const Double_t *eyhd)
: TGraph(n,x,y)
{
if (!CtorAllocate()) return;
n = sizeof(Double_t)*fNpoints;
if (exl) memcpy(fEXlow, exl, n);
else memset(fEXlow, 0, n);
if (exh) memcpy(fEXhigh, exh, n);
else memset(fEXhigh, 0, n);
if (eyl) memcpy(fEYlow, eyl, n);
else memset(fEYlow, 0, n);
if (eyh) memcpy(fEYhigh, eyh, n);
else memset(fEYhigh, 0, n);
if (exld) memcpy(fEXlowd, exld, n);
else memset(fEXlowd, 0, n);
if (exhd) memcpy(fEXhighd, exhd, n);
else memset(fEXhighd, 0, n);
if (eyld) memcpy(fEYlowd, eyld, n);
else memset(fEYlowd, 0, n);
if (eyhd) memcpy(fEYhighd, eyhd, n);
else memset(fEYhighd, 0, n);
}
TGraphBentErrors::~TGraphBentErrors()
{
delete [] fEXlow;
delete [] fEXhigh;
delete [] fEYlow;
delete [] fEYhigh;
delete [] fEXlowd;
delete [] fEXhighd;
delete [] fEYlowd;
delete [] fEYhighd;
}
void TGraphBentErrors::Apply(TF1 *f)
{
Double_t x,y,exl,exh,eyl,eyh,eyl_new,eyh_new,fxy;
for (Int_t i=0;i<GetN();i++) {
GetPoint(i,x,y);
exl=GetErrorXlow(i);
exh=GetErrorXhigh(i);
eyl=GetErrorYlow(i);
eyh=GetErrorYhigh(i);
fxy = f->Eval(x,y);
SetPoint(i,x,fxy);
if (f->Eval(x,y-eyl)<f->Eval(x,y+eyh)) {
eyl_new = TMath::Abs(fxy - f->Eval(x,y-eyl));
eyh_new = TMath::Abs(f->Eval(x,y+eyh) - fxy);
}
else {
eyh_new = TMath::Abs(fxy - f->Eval(x,y-eyl));
eyl_new = TMath::Abs(f->Eval(x,y+eyh) - fxy);
}
SetPointError(i,exl,exh,eyl_new,eyh_new);
}
}
void TGraphBentErrors::ComputeRange(Double_t &xmin, Double_t &ymin, Double_t &xmax, Double_t &ymax) const
{
TGraph::ComputeRange(xmin,ymin,xmax,ymax);
for (Int_t i=0;i<fNpoints;i++) {
if (fX[i] -fEXlow[i] < xmin) {
if (gPad && gPad->GetLogx()) {
if (fEXlow[i] < fX[i]) xmin = fX[i]-fEXlow[i];
else xmin = TMath::Min(xmin,fX[i]/3);
} else {
xmin = fX[i]-fEXlow[i];
}
}
if (fX[i] +fEXhigh[i] > xmax) xmax = fX[i]+fEXhigh[i];
if (fY[i] -fEYlow[i] < ymin) {
if (gPad && gPad->GetLogy()) {
if (fEYlow[i] < fY[i]) ymin = fY[i]-fEYlow[i];
else ymin = TMath::Min(ymin,fY[i]/3);
} else {
ymin = fY[i]-fEYlow[i];
}
}
if (fY[i] +fEYhigh[i] > ymax) ymax = fY[i]+fEYhigh[i];
}
}
void TGraphBentErrors::CopyAndRelease(Double_t **newarrays,
Int_t ibegin, Int_t iend, Int_t obegin)
{
CopyPoints(newarrays, ibegin, iend, obegin);
if (newarrays) {
delete[] fEXlow;
fEXlow = newarrays[0];
delete[] fEXhigh;
fEXhigh = newarrays[1];
delete[] fEYlow;
fEYlow = newarrays[2];
delete[] fEYhigh;
fEYhigh = newarrays[3];
delete[] fEXlowd;
fEXlowd = newarrays[4];
delete[] fEXhighd;
fEXhighd = newarrays[5];
delete[] fEYlowd;
fEYlowd = newarrays[6];
delete[] fEYhighd;
fEYhighd = newarrays[7];
delete[] fX;
fX = newarrays[8];
delete[] fY;
fY = newarrays[9];
delete[] newarrays;
}
}
Bool_t TGraphBentErrors::CopyPoints(Double_t **arrays,
Int_t ibegin, Int_t iend, Int_t obegin)
{
if (TGraph::CopyPoints(arrays ? arrays+8 : 0, ibegin, iend, obegin)) {
Int_t n = (iend - ibegin)*sizeof(Double_t);
if (arrays) {
memmove(&arrays[0][obegin], &fEXlow[ibegin], n);
memmove(&arrays[1][obegin], &fEXhigh[ibegin], n);
memmove(&arrays[2][obegin], &fEYlow[ibegin], n);
memmove(&arrays[3][obegin], &fEYhigh[ibegin], n);
memmove(&arrays[4][obegin], &fEXlowd[ibegin], n);
memmove(&arrays[5][obegin], &fEXhighd[ibegin], n);
memmove(&arrays[6][obegin], &fEYlowd[ibegin], n);
memmove(&arrays[7][obegin], &fEYhighd[ibegin], n);
} else {
memmove(&fEXlow[obegin], &fEXlow[ibegin], n);
memmove(&fEXhigh[obegin], &fEXhigh[ibegin], n);
memmove(&fEYlow[obegin], &fEYlow[ibegin], n);
memmove(&fEYhigh[obegin], &fEYhigh[ibegin], n);
memmove(&fEXlowd[obegin], &fEXlowd[ibegin], n);
memmove(&fEXhighd[obegin], &fEXhighd[ibegin], n);
memmove(&fEYlowd[obegin], &fEYlowd[ibegin], n);
memmove(&fEYhighd[obegin], &fEYhighd[ibegin], n);
}
return kTRUE;
} else {
return kFALSE;
}
}
Bool_t TGraphBentErrors::CtorAllocate(void)
{
if (!fNpoints) {
fEXlow = fEYlow = fEXhigh = fEYhigh = 0;
fEXlowd = fEYlowd = fEXhighd = fEYhighd = 0;
return kFALSE;
}
fEXlow = new Double_t[fMaxSize];
fEYlow = new Double_t[fMaxSize];
fEXhigh = new Double_t[fMaxSize];
fEYhigh = new Double_t[fMaxSize];
fEXlowd = new Double_t[fMaxSize];
fEYlowd = new Double_t[fMaxSize];
fEXhighd = new Double_t[fMaxSize];
fEYhighd = new Double_t[fMaxSize];
return kTRUE;
}
Double_t TGraphBentErrors::GetErrorX(Int_t i) const
{
if (i < 0 || i >= fNpoints) return -1;
if (!fEXlow && !fEXhigh) return -1;
Double_t elow=0, ehigh=0;
if (fEXlow) elow = fEXlow[i];
if (fEXhigh) ehigh = fEXhigh[i];
return TMath::Sqrt(0.5*(elow*elow + ehigh*ehigh));
}
Double_t TGraphBentErrors::GetErrorY(Int_t i) const
{
if (i < 0 || i >= fNpoints) return -1;
if (!fEYlow && !fEYhigh) return -1;
Double_t elow=0, ehigh=0;
if (fEYlow) elow = fEYlow[i];
if (fEYhigh) ehigh = fEYhigh[i];
return TMath::Sqrt(0.5*(elow*elow + ehigh*ehigh));
}
Double_t TGraphBentErrors::GetErrorXhigh(Int_t i) const
{
if (i<0 || i>fNpoints) return -1;
if (fEXhigh) return fEXhigh[i];
return -1;
}
Double_t TGraphBentErrors::GetErrorXlow(Int_t i) const
{
if (i<0 || i>fNpoints) return -1;
if (fEXlow) return fEXlow[i];
return -1;
}
Double_t TGraphBentErrors::GetErrorYhigh(Int_t i) const
{
if (i<0 || i>fNpoints) return -1;
if (fEYhigh) return fEYhigh[i];
return -1;
}
Double_t TGraphBentErrors::GetErrorYlow(Int_t i) const
{
if (i<0 || i>fNpoints) return -1;
if (fEYlow) return fEYlow[i];
return -1;
}
void TGraphBentErrors::FillZero(Int_t begin, Int_t end,
Bool_t from_ctor)
{
if (!from_ctor) {
TGraph::FillZero(begin, end, from_ctor);
}
Int_t n = (end - begin)*sizeof(Double_t);
memset(fEXlow + begin, 0, n);
memset(fEXhigh + begin, 0, n);
memset(fEYlow + begin, 0, n);
memset(fEYhigh + begin, 0, n);
memset(fEXlowd + begin, 0, n);
memset(fEXhighd + begin, 0, n);
memset(fEYlowd + begin, 0, n);
memset(fEYhighd + begin, 0, n);
}
void TGraphBentErrors::Print(Option_t *) const
{
for (Int_t i=0;i<fNpoints;i++) {
printf("x[%d]=%g, y[%d]=%g, exl[%d]=%g, exh[%d]=%g, eyl[%d]=%g, eyh[%d]=%g\n"
,i,fX[i],i,fY[i],i,fEXlow[i],i,fEXhigh[i],i,fEYlow[i],i,fEYhigh[i]);
}
}
void TGraphBentErrors::SavePrimitive(ostream &out, Option_t *option )
{
char quote = '"';
out<<" "<<endl;
if (gROOT->ClassSaved(TGraphBentErrors::Class())) {
out<<" ";
} else {
out<<" TGraphBentErrors *";
}
out<<"grbe = new TGraphBentErrors("<<fNpoints<<");"<<endl;
out<<" grbe->SetName("<<quote<<GetName()<<quote<<");"<<endl;
out<<" grbe->SetTitle("<<quote<<GetTitle()<<quote<<");"<<endl;
SaveFillAttributes(out,"grbe",0,1001);
SaveLineAttributes(out,"grbe",1,1,1);
SaveMarkerAttributes(out,"grbe",1,1,1);
for (Int_t i=0;i<fNpoints;i++) {
out<<" grbe->SetPoint("<<i<<","<<fX[i]<<","<<fY[i]<<");"<<endl;
out<<" grbe->SetPointError("<<i<<","<<fEXlow[i]<<","<<fEXhigh[i]
<<","<<fEYlow[i]<<","<<fEYhigh[i]
<<","<<fEXlowd[i]<<","<<fEXhighd[i]
<<","<<fEYlowd[i]<<","<<fEYhighd[i]
<<");"<<endl;
}
static Int_t frameNumber = 0;
if (fHistogram) {
frameNumber++;
TString hname = fHistogram->GetName();
hname += frameNumber;
fHistogram->SetName(hname.Data());
fHistogram->SavePrimitive(out,"nodraw");
out<<" grbe->SetHistogram("<<fHistogram->GetName()<<");"<<endl;
out<<" "<<endl;
}
TIter next(fFunctions);
TObject *obj;
while ((obj=next())) {
obj->SavePrimitive(out,"nodraw");
out<<" grbe->GetListOfFunctions()->Add("<<obj->GetName()<<");"<<endl;
if (obj->InheritsFrom("TPaveStats")) {
out<<" ptstats->SetParent(grbe->GetListOfFunctions());"<<endl;
}
}
const char *l = strstr(option,"multigraph");
if (l) {
out<<" multigraph->Add(grbe,"<<quote<<l+10<<quote<<");"<<endl;
} else {
out<<" grbe->Draw("<<quote<<option<<quote<<");"<<endl;
}
}
void TGraphBentErrors::SetPointError(Double_t exl, Double_t exh, Double_t eyl, Double_t eyh,
Double_t exld, Double_t exhd, Double_t eyld, Double_t eyhd)
{
Int_t px = gPad->GetEventX();
Int_t py = gPad->GetEventY();
Int_t ipoint = -2;
Int_t i;
for (i=0;i<fNpoints;i++) {
Int_t dpx = px - gPad->XtoAbsPixel(gPad->XtoPad(fX[i]));
Int_t dpy = py - gPad->YtoAbsPixel(gPad->YtoPad(fY[i]));
if (dpx*dpx+dpy*dpy < 25) {ipoint = i; break;}
}
if (ipoint == -2) return;
fEXlow[ipoint] = exl;
fEYlow[ipoint] = eyl;
fEXhigh[ipoint] = exh;
fEYhigh[ipoint] = eyh;
fEXlowd[ipoint] = exld;
fEXhighd[ipoint] = exhd;
fEYlowd[ipoint] = eyld;
fEYhighd[ipoint] = eyhd;
gPad->Modified();
}
void TGraphBentErrors::SetPointError(Int_t i, Double_t exl, Double_t exh, Double_t eyl, Double_t eyh,
Double_t exld, Double_t exhd, Double_t eyld, Double_t eyhd)
{
if (i < 0) return;
if (i >= fNpoints) {
TGraphBentErrors::SetPoint(i,0,0);
}
fEXlow[i] = exl;
fEYlow[i] = eyl;
fEXhigh[i] = exh;
fEYhigh[i] = eyh;
fEXlowd[i] = exld;
fEXhighd[i] = exhd;
fEYlowd[i] = eyld;
fEYhighd[i] = eyhd;
}
void TGraphBentErrors::SwapPoints(Int_t pos1, Int_t pos2)
{
SwapValues(fEXlow, pos1, pos2);
SwapValues(fEXhigh, pos1, pos2);
SwapValues(fEYlow, pos1, pos2);
SwapValues(fEYhigh, pos1, pos2);
SwapValues(fEXlowd, pos1, pos2);
SwapValues(fEXhighd, pos1, pos2);
SwapValues(fEYlowd, pos1, pos2);
SwapValues(fEYhighd, pos1, pos2);
TGraph::SwapPoints(pos1, pos2);
}