TGraphBentErrors.cxx

Go to the documentation of this file.
00001 // @(#)root/hist:$Id: TGraphBentErrors.cxx 33818 2010-06-10 13:37:21Z couet $
00002 // Author: Dave Morrison  30/06/2003
00003 
00004 /*************************************************************************
00005  * Copyright (C) 1995-2003, Rene Brun and Fons Rademakers.               *
00006  * All rights reserved.                                                  *
00007  *                                                                       *
00008  * For the licensing terms see $ROOTSYS/LICENSE.                         *
00009  * For the list of contributors see $ROOTSYS/README/CREDITS.             *
00010  *************************************************************************/
00011 
00012 #include <string.h>
00013 
00014 #include "Riostream.h"
00015 #include "TROOT.h"
00016 #include "TGraphBentErrors.h"
00017 #include "TStyle.h"
00018 #include "TMath.h"
00019 #include "TArrow.h"
00020 #include "TBox.h"
00021 #include "TVirtualPad.h"
00022 #include "TH1.h"
00023 #include "TF1.h"
00024 
00025 ClassImp(TGraphBentErrors)
00026 
00027 
00028 //______________________________________________________________________________
00029 /* Begin_Html
00030 <center><h2>TGraphBentErrors class</h2></center>
00031 A TGraphBentErrors is a TGraph with bent, assymetric error bars.
00032 <p>
00033 The TGraphBentErrors painting is permofed thanks to the
00034 <a href="http://root.cern.ch/root/html/TGraphPainter.html">TGraphPainter</a>
00035 class. All details about the various painting options are given in
00036 <a href="http://root.cern.ch/root/html/TGraphPainter.html">this class</a>.
00037 <p>
00038 The picture below gives an example:
00039 End_Html
00040 Begin_Macro(source)
00041 {
00042    c1 = new TCanvas("c1","A Simple Graph with bent error bars",200,10,700,500);
00043    Int_t n = 10;
00044    Double_t x[n]    = {-0.22, 0.05, 0.25, 0.35, 0.5, 0.61,0.7,0.85,0.89,0.95};
00045    Double_t y[n]    = {1,2.9,5.6,7.4,9,9.6,8.7,6.3,4.5,1};
00046    Double_t exl[n]  = {.05,.1,.07,.07,.04,.05,.06,.07,.08,.05};
00047    Double_t eyl[n]  = {.8,.7,.6,.5,.4,.4,.5,.6,.7,.8};
00048    Double_t exh[n]  = {.02,.08,.05,.05,.03,.03,.04,.05,.06,.03};
00049    Double_t eyh[n]  = {.6,.5,.4,.3,.2,.2,.3,.4,.5,.6};
00050    Double_t exld[n] = {.0,.0,.0,.0,.0,.0,.0,.0,.0,.0};
00051    Double_t eyld[n] = {.0,.0,.05,.0,.0,.0,.0,.0,.0,.0};
00052    Double_t exhd[n] = {.0,.0,.0,.0,.0,.0,.0,.0,.0,.0};
00053    Double_t eyhd[n] = {.0,.0,.0,.0,.0,.0,.0,.0,.05,.0};
00054    gr = new TGraphBentErrors(n,x,y,exl,exh,eyl,eyh,exld,exhd,eyld,eyhd);
00055    gr->SetTitle("TGraphBentErrors Example");
00056    gr->SetMarkerColor(4);
00057    gr->SetMarkerStyle(21);
00058    gr->Draw("ALP");
00059    return c1;
00060 }
00061 End_Macro */
00062 
00063 
00064 //_____________________________________________________________________________
00065 TGraphBentErrors::TGraphBentErrors(): TGraph()
00066 {
00067    // TGraphBentErrors default constructor.
00068 
00069    if (!CtorAllocate()) return;
00070 }
00071 
00072 
00073 //______________________________________________________________________________
00074 TGraphBentErrors::TGraphBentErrors(const TGraphBentErrors &gr)
00075        : TGraph(gr)
00076 {
00077    // TGraphBentErrors copy constructor
00078 
00079    if (!CtorAllocate()) return;
00080    Int_t n = fNpoints*sizeof(Double_t);
00081    memcpy(fEXlow, gr.fEXlow, n);
00082    memcpy(fEYlow, gr.fEYlow, n);
00083    memcpy(fEXhigh, gr.fEXhigh, n);
00084    memcpy(fEYhigh, gr.fEYhigh, n);
00085    memcpy(fEXlowd, gr.fEXlowd, n);
00086    memcpy(fEYlowd, gr.fEYlowd, n);
00087    memcpy(fEXhighd, gr.fEXhighd, n);
00088    memcpy(fEYhighd, gr.fEYhighd, n);
00089 }
00090 
00091 
00092 //_____________________________________________________________________________
00093 TGraphBentErrors::TGraphBentErrors(Int_t n)
00094        : TGraph(n)
00095 {
00096    // TGraphBentErrors normal constructor.
00097    //
00098    //  the arrays are preset to zero
00099 
00100    if (!CtorAllocate()) return;
00101    FillZero(0, fNpoints);
00102 }
00103 
00104 
00105 //_____________________________________________________________________________
00106 TGraphBentErrors::TGraphBentErrors(Int_t n,
00107                                    const Float_t *x, const Float_t *y,
00108                                    const Float_t *exl, const Float_t *exh,
00109                                    const Float_t *eyl, const Float_t *eyh,
00110                                    const Float_t *exld, const Float_t *exhd,
00111                                    const Float_t *eyld, const Float_t *eyhd)
00112   : TGraph(n,x,y)
00113 {
00114    // TGraphBentErrors normal constructor.
00115    //
00116    // if exl,h or eyl,h are null, the corresponding arrays are preset to zero
00117 
00118    if (!CtorAllocate()) return;
00119 
00120    for (Int_t i=0;i<n;i++) {
00121       if (exl) fEXlow[i]  = exl[i];
00122       else     fEXlow[i]  = 0;
00123       if (exh) fEXhigh[i] = exh[i];
00124       else     fEXhigh[i] = 0;
00125       if (eyl) fEYlow[i]  = eyl[i];
00126       else     fEYlow[i]  = 0;
00127       if (eyh) fEYhigh[i] = eyh[i];
00128       else     fEYhigh[i] = 0;
00129 
00130       if (exld) fEXlowd[i]  = exld[i];
00131       else     fEXlowd[i]  = 0;
00132       if (exhd) fEXhighd[i] = exhd[i];
00133       else     fEXhighd[i] = 0;
00134       if (eyld) fEYlowd[i]  = eyld[i];
00135       else     fEYlowd[i]  = 0;
00136       if (eyhd) fEYhighd[i] = eyhd[i];
00137       else     fEYhighd[i] = 0;
00138    }
00139 }
00140 
00141 
00142 //_____________________________________________________________________________
00143 TGraphBentErrors::TGraphBentErrors(Int_t n,
00144                                    const Double_t *x, const Double_t *y,
00145                                    const Double_t *exl, const Double_t *exh,
00146                                    const Double_t *eyl, const Double_t *eyh,
00147                                    const Double_t *exld, const Double_t *exhd,
00148                                    const Double_t *eyld, const Double_t *eyhd)
00149   : TGraph(n,x,y)
00150 {
00151    // TGraphBentErrors normal constructor.
00152    //
00153    // if exl,h or eyl,h are null, the corresponding arrays are preset to zero
00154 
00155    if (!CtorAllocate()) return;
00156    n = sizeof(Double_t)*fNpoints;
00157 
00158       if (exl) memcpy(fEXlow, exl, n);
00159       else     memset(fEXlow, 0, n);
00160       if (exh) memcpy(fEXhigh, exh, n);
00161       else     memset(fEXhigh, 0, n);
00162       if (eyl) memcpy(fEYlow, eyl, n);
00163       else     memset(fEYlow, 0, n);
00164       if (eyh) memcpy(fEYhigh, eyh, n);
00165       else     memset(fEYhigh, 0, n);
00166 
00167       if (exld) memcpy(fEXlowd, exld, n);
00168       else      memset(fEXlowd, 0, n);
00169       if (exhd) memcpy(fEXhighd, exhd, n);
00170       else      memset(fEXhighd, 0, n);
00171       if (eyld) memcpy(fEYlowd,  eyld, n);
00172       else      memset(fEYlowd, 0, n);
00173       if (eyhd) memcpy(fEYhighd, eyhd, n);
00174       else      memset(fEYhighd, 0, n);
00175 }
00176 
00177 
00178 //_____________________________________________________________________________
00179 TGraphBentErrors::~TGraphBentErrors()
00180 {
00181    // TGraphBentErrors default destructor.
00182 
00183    delete [] fEXlow;
00184    delete [] fEXhigh;
00185    delete [] fEYlow;
00186    delete [] fEYhigh;
00187 
00188    delete [] fEXlowd;
00189    delete [] fEXhighd;
00190    delete [] fEYlowd;
00191    delete [] fEYhighd;
00192 }
00193 
00194 
00195 //_____________________________________________________________________________
00196 void TGraphBentErrors::Apply(TF1 *f)
00197 {
00198    // apply a function to all data points
00199    // y = f(x,y)
00200    //
00201    // Errors are calculated as eyh = f(x,y+eyh)-f(x,y) and
00202    // eyl = f(x,y)-f(x,y-eyl)
00203    //
00204    // Special treatment has to be applied for the functions where the
00205    // role of "up" and "down" is reversed.
00206    // function suggested/implemented by Miroslav Helbich <helbich@mail.desy.de>
00207 
00208    Double_t x,y,exl,exh,eyl,eyh,eyl_new,eyh_new,fxy;
00209 
00210    for (Int_t i=0;i<GetN();i++) {
00211       GetPoint(i,x,y);
00212       exl=GetErrorXlow(i);
00213       exh=GetErrorXhigh(i);
00214       eyl=GetErrorYlow(i);
00215       eyh=GetErrorYhigh(i);
00216 
00217       fxy = f->Eval(x,y);
00218       SetPoint(i,x,fxy);
00219 
00220       // in the case of the functions like y-> -1*y the roles of the
00221       // upper and lower error bars is reversed
00222       if (f->Eval(x,y-eyl)<f->Eval(x,y+eyh)) {
00223          eyl_new = TMath::Abs(fxy - f->Eval(x,y-eyl));
00224          eyh_new = TMath::Abs(f->Eval(x,y+eyh) - fxy);
00225       }
00226       else {
00227          eyh_new = TMath::Abs(fxy - f->Eval(x,y-eyl));
00228          eyl_new = TMath::Abs(f->Eval(x,y+eyh) - fxy);
00229       }
00230 
00231       //error on x doesn't change
00232       SetPointError(i,exl,exh,eyl_new,eyh_new);
00233    }
00234 }
00235 
00236 
00237 //_____________________________________________________________________________
00238 void TGraphBentErrors::ComputeRange(Double_t &xmin, Double_t &ymin, Double_t &xmax, Double_t &ymax) const
00239 {
00240    // Compute range.
00241 
00242    TGraph::ComputeRange(xmin,ymin,xmax,ymax);
00243 
00244    for (Int_t i=0;i<fNpoints;i++) {
00245       if (fX[i] -fEXlow[i] < xmin) {
00246          if (gPad && gPad->GetLogx()) {
00247             if (fEXlow[i] < fX[i]) xmin = fX[i]-fEXlow[i];
00248             else                   xmin = TMath::Min(xmin,fX[i]/3);
00249          } else {
00250             xmin = fX[i]-fEXlow[i];
00251          }
00252       }
00253       if (fX[i] +fEXhigh[i] > xmax) xmax = fX[i]+fEXhigh[i];
00254       if (fY[i] -fEYlow[i] < ymin) {
00255          if (gPad && gPad->GetLogy()) {
00256             if (fEYlow[i] < fY[i]) ymin = fY[i]-fEYlow[i];
00257             else                   ymin = TMath::Min(ymin,fY[i]/3);
00258          } else {
00259             ymin = fY[i]-fEYlow[i];
00260          }
00261       }
00262       if (fY[i] +fEYhigh[i] > ymax) ymax = fY[i]+fEYhigh[i];
00263    }
00264 }
00265 
00266 
00267 //______________________________________________________________________________
00268 void TGraphBentErrors::CopyAndRelease(Double_t **newarrays,
00269                                       Int_t ibegin, Int_t iend, Int_t obegin)
00270 {
00271    // Copy and release.
00272 
00273    CopyPoints(newarrays, ibegin, iend, obegin);
00274    if (newarrays) {
00275       delete[] fEXlow;
00276       fEXlow = newarrays[0];
00277       delete[] fEXhigh;
00278       fEXhigh = newarrays[1];
00279       delete[] fEYlow;
00280       fEYlow = newarrays[2];
00281       delete[] fEYhigh;
00282       fEYhigh = newarrays[3];
00283       delete[] fEXlowd;
00284       fEXlowd = newarrays[4];
00285       delete[] fEXhighd;
00286       fEXhighd = newarrays[5];
00287       delete[] fEYlowd;
00288       fEYlowd = newarrays[6];
00289       delete[] fEYhighd;
00290       fEYhighd = newarrays[7];
00291       delete[] fX;
00292       fX = newarrays[8];
00293       delete[] fY;
00294       fY = newarrays[9];
00295       delete[] newarrays;
00296    }
00297 }
00298 
00299 
00300 //______________________________________________________________________________
00301 Bool_t TGraphBentErrors::CopyPoints(Double_t **arrays,
00302                                     Int_t ibegin, Int_t iend, Int_t obegin)
00303 {
00304    // Copy errors from fE*** to arrays[***]
00305    // or to f*** Copy points.
00306 
00307    if (TGraph::CopyPoints(arrays ? arrays+8 : 0, ibegin, iend, obegin)) {
00308       Int_t n = (iend - ibegin)*sizeof(Double_t);
00309       if (arrays) {
00310          memmove(&arrays[0][obegin], &fEXlow[ibegin], n);
00311          memmove(&arrays[1][obegin], &fEXhigh[ibegin], n);
00312          memmove(&arrays[2][obegin], &fEYlow[ibegin], n);
00313          memmove(&arrays[3][obegin], &fEYhigh[ibegin], n);
00314          memmove(&arrays[4][obegin], &fEXlowd[ibegin], n);
00315          memmove(&arrays[5][obegin], &fEXhighd[ibegin], n);
00316          memmove(&arrays[6][obegin], &fEYlowd[ibegin], n);
00317          memmove(&arrays[7][obegin], &fEYhighd[ibegin], n);
00318       } else {
00319          memmove(&fEXlow[obegin], &fEXlow[ibegin], n);
00320          memmove(&fEXhigh[obegin], &fEXhigh[ibegin], n);
00321          memmove(&fEYlow[obegin], &fEYlow[ibegin], n);
00322          memmove(&fEYhigh[obegin], &fEYhigh[ibegin], n);
00323          memmove(&fEXlowd[obegin], &fEXlowd[ibegin], n);
00324          memmove(&fEXhighd[obegin], &fEXhighd[ibegin], n);
00325          memmove(&fEYlowd[obegin], &fEYlowd[ibegin], n);
00326          memmove(&fEYhighd[obegin], &fEYhighd[ibegin], n);
00327       }
00328       return kTRUE;
00329    } else {
00330       return kFALSE;
00331    }
00332 }
00333 
00334 
00335 //______________________________________________________________________________
00336 Bool_t TGraphBentErrors::CtorAllocate(void)
00337 {
00338    // Should be called from ctors after fNpoints has been set
00339 
00340    if (!fNpoints) {
00341       fEXlow = fEYlow = fEXhigh = fEYhigh = 0;
00342       fEXlowd = fEYlowd = fEXhighd = fEYhighd = 0;
00343       return kFALSE;
00344    }
00345    fEXlow = new Double_t[fMaxSize];
00346    fEYlow = new Double_t[fMaxSize];
00347    fEXhigh = new Double_t[fMaxSize];
00348    fEYhigh = new Double_t[fMaxSize];
00349    fEXlowd = new Double_t[fMaxSize];
00350    fEYlowd = new Double_t[fMaxSize];
00351    fEXhighd = new Double_t[fMaxSize];
00352    fEYhighd = new Double_t[fMaxSize];
00353    return kTRUE;
00354 }
00355 
00356 
00357 //_____________________________________________________________________________
00358 Double_t TGraphBentErrors::GetErrorX(Int_t i) const
00359 {
00360    // This function is called by GraphFitChisquare.
00361    // It returns the error along X at point i.
00362 
00363    if (i < 0 || i >= fNpoints) return -1;
00364    if (!fEXlow && !fEXhigh) return -1;
00365    Double_t elow=0, ehigh=0;
00366    if (fEXlow)  elow  = fEXlow[i];
00367    if (fEXhigh) ehigh = fEXhigh[i];
00368    return TMath::Sqrt(0.5*(elow*elow + ehigh*ehigh));
00369 }
00370 
00371 
00372 //_____________________________________________________________________________
00373 Double_t TGraphBentErrors::GetErrorY(Int_t i) const
00374 {
00375    // This function is called by GraphFitChisquare.
00376    // It returns the error along Y at point i.
00377 
00378    if (i < 0 || i >= fNpoints) return -1;
00379    if (!fEYlow && !fEYhigh) return -1;
00380    Double_t elow=0, ehigh=0;
00381    if (fEYlow)  elow  = fEYlow[i];
00382    if (fEYhigh) ehigh = fEYhigh[i];
00383    return TMath::Sqrt(0.5*(elow*elow + ehigh*ehigh));
00384 }
00385 
00386 
00387 //______________________________________________________________________________
00388 Double_t TGraphBentErrors::GetErrorXhigh(Int_t i) const
00389 {
00390    // Get high error on X[i].
00391 
00392    if (i<0 || i>fNpoints) return -1;
00393    if (fEXhigh) return fEXhigh[i];
00394    return -1;
00395 }
00396 
00397 
00398 //______________________________________________________________________________
00399 Double_t TGraphBentErrors::GetErrorXlow(Int_t i) const
00400 {
00401    // Get low error on X[i].
00402 
00403    if (i<0 || i>fNpoints) return -1;
00404    if (fEXlow) return fEXlow[i];
00405    return -1;
00406 }
00407 
00408 
00409 //______________________________________________________________________________
00410 Double_t TGraphBentErrors::GetErrorYhigh(Int_t i) const
00411 {
00412    // Get high error on Y[i].
00413 
00414    if (i<0 || i>fNpoints) return -1;
00415    if (fEYhigh) return fEYhigh[i];
00416    return -1;
00417 }
00418 
00419 
00420 //______________________________________________________________________________
00421 Double_t TGraphBentErrors::GetErrorYlow(Int_t i) const
00422 {
00423    // Get low error on Y[i].
00424 
00425    if (i<0 || i>fNpoints) return -1;
00426    if (fEYlow) return fEYlow[i];
00427    return -1;
00428 }
00429 
00430 
00431 //______________________________________________________________________________
00432 void TGraphBentErrors::FillZero(Int_t begin, Int_t end,
00433                                  Bool_t from_ctor)
00434 {
00435    // Set zero values for point arrays in the range [begin, end)
00436 
00437    if (!from_ctor) {
00438       TGraph::FillZero(begin, end, from_ctor);
00439    }
00440    Int_t n = (end - begin)*sizeof(Double_t);
00441    memset(fEXlow + begin, 0, n);
00442    memset(fEXhigh + begin, 0, n);
00443    memset(fEYlow + begin, 0, n);
00444    memset(fEYhigh + begin, 0, n);
00445    memset(fEXlowd + begin, 0, n);
00446    memset(fEXhighd + begin, 0, n);
00447    memset(fEYlowd + begin, 0, n);
00448    memset(fEYhighd + begin, 0, n);
00449 }
00450 
00451 
00452 //______________________________________________________________________________
00453 void TGraphBentErrors::Print(Option_t *) const
00454 {
00455    // Print graph and errors values.
00456 
00457    for (Int_t i=0;i<fNpoints;i++) {
00458       printf("x[%d]=%g, y[%d]=%g, exl[%d]=%g, exh[%d]=%g, eyl[%d]=%g, eyh[%d]=%g\n"
00459          ,i,fX[i],i,fY[i],i,fEXlow[i],i,fEXhigh[i],i,fEYlow[i],i,fEYhigh[i]);
00460    }
00461 }
00462 
00463 
00464 //______________________________________________________________________________
00465 void TGraphBentErrors::SavePrimitive(ostream &out, Option_t *option /*= ""*/)
00466 {
00467    // Save primitive as a C++ statement(s) on output stream out
00468 
00469    char quote = '"';
00470    out<<"   "<<endl;
00471    if (gROOT->ClassSaved(TGraphBentErrors::Class())) {
00472       out<<"   ";
00473    } else {
00474       out<<"   TGraphBentErrors *";
00475    }
00476    out<<"grbe = new TGraphBentErrors("<<fNpoints<<");"<<endl;
00477    out<<"   grbe->SetName("<<quote<<GetName()<<quote<<");"<<endl;
00478    out<<"   grbe->SetTitle("<<quote<<GetTitle()<<quote<<");"<<endl;
00479 
00480    SaveFillAttributes(out,"grbe",0,1001);
00481    SaveLineAttributes(out,"grbe",1,1,1);
00482    SaveMarkerAttributes(out,"grbe",1,1,1);
00483 
00484    for (Int_t i=0;i<fNpoints;i++) {
00485       out<<"   grbe->SetPoint("<<i<<","<<fX[i]<<","<<fY[i]<<");"<<endl;
00486       out<<"   grbe->SetPointError("<<i<<","<<fEXlow[i]<<","<<fEXhigh[i]
00487                                        <<","<<fEYlow[i]<<","<<fEYhigh[i]
00488                                        <<","<<fEXlowd[i]<<","<<fEXhighd[i]
00489                                        <<","<<fEYlowd[i]<<","<<fEYhighd[i]
00490                                        <<");"<<endl;
00491    }
00492 
00493    static Int_t frameNumber = 0;
00494    if (fHistogram) {
00495       frameNumber++;
00496       TString hname = fHistogram->GetName();
00497       hname += frameNumber;
00498       fHistogram->SetName(hname.Data());
00499       fHistogram->SavePrimitive(out,"nodraw");
00500       out<<"   grbe->SetHistogram("<<fHistogram->GetName()<<");"<<endl;
00501       out<<"   "<<endl;
00502    }
00503 
00504    // save list of functions
00505    TIter next(fFunctions);
00506    TObject *obj;
00507    while ((obj=next())) {
00508       obj->SavePrimitive(out,"nodraw");
00509       out<<"   grbe->GetListOfFunctions()->Add("<<obj->GetName()<<");"<<endl;
00510       if (obj->InheritsFrom("TPaveStats")) {
00511          out<<"   ptstats->SetParent(grbe->GetListOfFunctions());"<<endl;
00512       }
00513    }
00514 
00515    const char *l = strstr(option,"multigraph");
00516    if (l) {
00517       out<<"   multigraph->Add(grbe,"<<quote<<l+10<<quote<<");"<<endl;
00518    } else {
00519       out<<"   grbe->Draw("<<quote<<option<<quote<<");"<<endl;
00520    }
00521 }
00522 
00523 
00524 //______________________________________________________________________________
00525 void TGraphBentErrors::SetPointError(Double_t exl, Double_t exh, Double_t eyl, Double_t eyh,
00526                                      Double_t exld, Double_t exhd, Double_t eyld, Double_t eyhd)
00527 {
00528    // Set ex and ey values for point pointed by the mouse.
00529 
00530    Int_t px = gPad->GetEventX();
00531    Int_t py = gPad->GetEventY();
00532 
00533    //localize point to be deleted
00534    Int_t ipoint = -2;
00535    Int_t i;
00536    // start with a small window (in case the mouse is very close to one point)
00537    for (i=0;i<fNpoints;i++) {
00538       Int_t dpx = px - gPad->XtoAbsPixel(gPad->XtoPad(fX[i]));
00539       Int_t dpy = py - gPad->YtoAbsPixel(gPad->YtoPad(fY[i]));
00540       if (dpx*dpx+dpy*dpy < 25) {ipoint = i; break;}
00541    }
00542    if (ipoint == -2) return;
00543 
00544    fEXlow[ipoint]   = exl;
00545    fEYlow[ipoint]   = eyl;
00546    fEXhigh[ipoint]  = exh;
00547    fEYhigh[ipoint]  = eyh;
00548    fEXlowd[ipoint]  = exld;
00549    fEXhighd[ipoint] = exhd;
00550    fEYlowd[ipoint]  = eyld;
00551    fEYhighd[ipoint] = eyhd;
00552    gPad->Modified();
00553 }
00554 
00555 
00556 //______________________________________________________________________________
00557 void TGraphBentErrors::SetPointError(Int_t i, Double_t exl, Double_t exh, Double_t eyl, Double_t eyh,
00558                                      Double_t exld, Double_t exhd, Double_t eyld, Double_t eyhd)
00559 {
00560    // Set ex and ey values for point number i.
00561 
00562    if (i < 0) return;
00563    if (i >= fNpoints) {
00564    // re-allocate the object
00565       TGraphBentErrors::SetPoint(i,0,0);
00566    }
00567    fEXlow[i]   = exl;
00568    fEYlow[i]   = eyl;
00569    fEXhigh[i]  = exh;
00570    fEYhigh[i]  = eyh;
00571    fEXlowd[i]  = exld;
00572    fEXhighd[i] = exhd;
00573    fEYlowd[i]  = eyld;
00574    fEYhighd[i] = eyhd;
00575 }
00576 
00577 
00578 //_____________________________________________________________________________
00579 void TGraphBentErrors::SwapPoints(Int_t pos1, Int_t pos2)
00580 {
00581    // Swap points.
00582 
00583    SwapValues(fEXlow,  pos1, pos2);
00584    SwapValues(fEXhigh, pos1, pos2);
00585    SwapValues(fEYlow,  pos1, pos2);
00586    SwapValues(fEYhigh, pos1, pos2);
00587 
00588    SwapValues(fEXlowd,  pos1, pos2);
00589    SwapValues(fEXhighd, pos1, pos2);
00590    SwapValues(fEYlowd,  pos1, pos2);
00591    SwapValues(fEYhighd, pos1, pos2);
00592 
00593    TGraph::SwapPoints(pos1, pos2);
00594 }

Generated on Tue Jul 5 14:24:18 2011 for ROOT_528-00b_version by  doxygen 1.5.1