TGraph.cxx

Go to the documentation of this file.
00001 // @(#)root/hist:$Id: TGraph.cxx 36530 2010-11-08 10:24:42Z couet $
00002 // Author: Rene Brun, Olivier Couet   12/12/94
00003 
00004 /*************************************************************************
00005  * Copyright (C) 1995-2000, 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 "TEnv.h"
00017 #include "TGraph.h"
00018 #include "TGaxis.h"
00019 #include "TH1.h"
00020 #include "TF1.h"
00021 #include "TStyle.h"
00022 #include "TMath.h"
00023 #include "TFrame.h"
00024 #include "TVector.h"
00025 #include "TVectorD.h"
00026 #include "Foption.h"
00027 #include "TRandom.h"
00028 #include "TSpline.h"
00029 #include "TVirtualFitter.h"
00030 #include "TVirtualPad.h"
00031 #include "TVirtualGraphPainter.h"
00032 #include "TBrowser.h"
00033 #include "TClass.h"
00034 #include "TSystem.h"
00035 #include "TPluginManager.h"
00036 #include <stdlib.h>
00037 #include <string>
00038 #include <cassert>
00039 
00040 #include "HFitInterface.h"
00041 #include "Fit/DataRange.h"
00042 #include "Math/MinimizerOptions.h"
00043 
00044 extern void H1LeastSquareSeqnd(Int_t n, Double_t *a, Int_t idim, Int_t &ifail, Int_t k, Double_t *b);
00045 
00046 ClassImp(TGraph)
00047 
00048 
00049 //______________________________________________________________________________
00050 /* Begin_Html
00051 <center><h2>Graph class</h2></center>
00052 A Graph is a graphics object made of two arrays X and Y with npoints each.
00053 <p>
00054 The TGraph painting is permofed thanks to the
00055 <a href="http://root.cern.ch/root/html/TGraphPainter.html">TGraphPainter</a>
00056 class. All details about the various painting options are given in
00057 <a href="http://root.cern.ch/root/html/TGraphPainter.html">this class</a>.
00058 <p>
00059 The picture below gives an example:
00060 End_Html
00061 Begin_Macro(source)
00062 {
00063    TCanvas *c1 = new TCanvas("c1","A Simple Graph Example",200,10,700,500);
00064    Double_t x[100], y[100];
00065    Int_t n = 20;
00066    for (Int_t i=0;i<n;i++) {
00067      x[i] = i*0.1;
00068      y[i] = 10*sin(x[i]+0.2);
00069    }
00070    gr = new TGraph(n,x,y);
00071    gr->Draw("AC*");
00072    return c1;
00073 }
00074 End_Macro */
00075 
00076 
00077 //______________________________________________________________________________
00078 TGraph::TGraph(): TNamed(), TAttLine(), TAttFill(1,1001), TAttMarker()
00079 {
00080    // Graph default constructor.
00081 
00082    fNpoints = -1;  //will be reset to 0 in CtorAllocate
00083    if (!CtorAllocate()) return;
00084 }
00085 
00086 
00087 //______________________________________________________________________________
00088 TGraph::TGraph(Int_t n)
00089        : TNamed("Graph","Graph"), TAttLine(), TAttFill(1,1001), TAttMarker()
00090 {
00091    // Constructor with only the number of points set
00092    // the arrsys x and y will be set later
00093 
00094    fNpoints = n;
00095    if (!CtorAllocate()) return;
00096    FillZero(0, fNpoints);
00097 }
00098 
00099 
00100 //______________________________________________________________________________
00101 TGraph::TGraph(Int_t n, const Int_t *x, const Int_t *y)
00102        : TNamed("Graph","Graph"), TAttLine(), TAttFill(1,1001), TAttMarker()
00103 {
00104    // Graph normal constructor with ints.
00105 
00106    if (!x || !y) {
00107       fNpoints = 0;
00108    } else {
00109       fNpoints = n;
00110    }
00111    if (!CtorAllocate()) return;
00112    for (Int_t i=0;i<n;i++) {
00113       fX[i] = (Double_t)x[i];
00114       fY[i] = (Double_t)y[i];
00115    }
00116 }
00117 
00118 
00119 //______________________________________________________________________________
00120 TGraph::TGraph(Int_t n, const Float_t *x, const Float_t *y)
00121        : TNamed("Graph","Graph"), TAttLine(), TAttFill(1,1001), TAttMarker()
00122 {
00123    // Graph normal constructor with floats.
00124 
00125    if (!x || !y) {
00126       fNpoints = 0;
00127    } else {
00128       fNpoints = n;
00129    }
00130    if (!CtorAllocate()) return;
00131    for (Int_t i=0;i<n;i++) {
00132       fX[i] = x[i];
00133       fY[i] = y[i];
00134    }
00135 }
00136 
00137 
00138 //______________________________________________________________________________
00139 TGraph::TGraph(Int_t n, const Double_t *x, const Double_t *y)
00140        : TNamed("Graph","Graph"), TAttLine(), TAttFill(1,1001), TAttMarker()
00141 {
00142    // Graph normal constructor with doubles.
00143 
00144    if (!x || !y) {
00145       fNpoints = 0;
00146    } else {
00147       fNpoints = n;
00148    }
00149    if (!CtorAllocate()) return;
00150    n = fNpoints*sizeof(Double_t);
00151    memcpy(fX, x, n);
00152    memcpy(fY, y, n);
00153 }
00154 
00155 
00156 //______________________________________________________________________________
00157 TGraph::TGraph(const TGraph &gr)
00158        : TNamed(gr), TAttLine(gr), TAttFill(gr), TAttMarker(gr)
00159 {
00160    // Copy constructor for this graph
00161 
00162    fNpoints = gr.fNpoints;
00163    fMaxSize = gr.fMaxSize;
00164    if (gr.fFunctions) fFunctions = (TList*)gr.fFunctions->Clone();
00165    else fFunctions = new TList;
00166    fHistogram = 0;
00167    fMinimum = gr.fMinimum;
00168    fMaximum = gr.fMaximum;
00169    if (!fMaxSize) {
00170       fX = fY = 0;
00171       return;
00172    } else {
00173       fX = new Double_t[fMaxSize];
00174       fY = new Double_t[fMaxSize];
00175    }
00176 
00177    Int_t n = gr.GetN()*sizeof(Double_t);
00178    memcpy(fX, gr.fX, n);
00179    memcpy(fY, gr.fY, n);
00180 }
00181 
00182 
00183 //______________________________________________________________________________
00184 TGraph& TGraph::operator=(const TGraph &gr)
00185 {
00186    // Equal operator for this graph
00187 
00188    if(this!=&gr) {
00189       TNamed::operator=(gr);
00190       TAttLine::operator=(gr);
00191       TAttFill::operator=(gr);
00192       TAttMarker::operator=(gr);
00193 
00194       fNpoints = gr.fNpoints;
00195       fMaxSize = gr.fMaxSize;
00196       if (gr.fFunctions) fFunctions = (TList*)gr.fFunctions->Clone();
00197       else fFunctions = new TList;
00198       if (gr.fHistogram) fHistogram = new TH1F(*(gr.fHistogram));
00199       else fHistogram = 0;
00200       fMinimum = gr.fMinimum;
00201       fMaximum = gr.fMaximum;
00202       if (!fMaxSize) {
00203          fX = fY = 0;
00204          return *this;
00205       } else {
00206          fX = new Double_t[fMaxSize];
00207          fY = new Double_t[fMaxSize];
00208       }
00209 
00210       Int_t n = gr.GetN()*sizeof(Double_t);
00211       if (n>0) {
00212          memcpy(fX, gr.fX, n);
00213          memcpy(fY, gr.fY, n);
00214       }
00215    }
00216    return *this;
00217 }
00218 
00219 
00220 //______________________________________________________________________________
00221 TGraph::TGraph(const TVectorF &vx, const TVectorF &vy)
00222        : TNamed("Graph","Graph"), TAttLine(), TAttFill(1,1001), TAttMarker()
00223 {
00224    // Graph constructor with two vectors of floats in input
00225    // A graph is build with the X coordinates taken from vx and Y coord from vy
00226    // The number of points in the graph is the minimum of number of points
00227    // in vx and vy.
00228 
00229    fNpoints = TMath::Min(vx.GetNrows(), vy.GetNrows());
00230    if (!CtorAllocate()) return;
00231    Int_t ivxlow  = vx.GetLwb();
00232    Int_t ivylow  = vy.GetLwb();
00233    for (Int_t i=0;i<fNpoints;i++) {
00234       fX[i]  = vx(i+ivxlow);
00235       fY[i]  = vy(i+ivylow);
00236    }
00237 }
00238 
00239 
00240 //______________________________________________________________________________
00241 TGraph::TGraph(const TVectorD &vx, const TVectorD &vy)
00242        : TNamed("Graph","Graph"), TAttLine(), TAttFill(1,1001), TAttMarker()
00243 {
00244    // Graph constructor with two vectors of doubles in input
00245    // A graph is build with the X coordinates taken from vx and Y coord from vy
00246    // The number of points in the graph is the minimum of number of points
00247    // in vx and vy.
00248 
00249    fNpoints = TMath::Min(vx.GetNrows(), vy.GetNrows());
00250    if (!CtorAllocate()) return;
00251    Int_t ivxlow  = vx.GetLwb();
00252    Int_t ivylow  = vy.GetLwb();
00253    for (Int_t i=0;i<fNpoints;i++) {
00254       fX[i]  = vx(i+ivxlow);
00255       fY[i]  = vy(i+ivylow);
00256    }
00257 }
00258 
00259 
00260 //______________________________________________________________________________
00261 TGraph::TGraph(const TH1 *h)
00262        : TNamed("Graph","Graph"), TAttLine(), TAttFill(1,1001), TAttMarker()
00263 {
00264    // Graph constructor importing its parameters from the TH1 object passed as argument
00265 
00266    if (!h) {
00267       Error("TGraph", "Pointer to histogram is null");
00268       fNpoints = 0;
00269       return;
00270    }
00271    if (h->GetDimension() != 1) {
00272       Error("TGraph", "Histogram must be 1-D; h %s is %d-D",h->GetName(),h->GetDimension());
00273       fNpoints = 0;
00274    } else {
00275       fNpoints = ((TH1*)h)->GetXaxis()->GetNbins();
00276    }
00277 
00278    if (!CtorAllocate()) return;
00279 
00280    TAxis *xaxis = ((TH1*)h)->GetXaxis();
00281    for (Int_t i=0;i<fNpoints;i++) {
00282       fX[i] = xaxis->GetBinCenter(i+1);
00283       fY[i] = h->GetBinContent(i+1);
00284    }
00285    h->TAttLine::Copy(*this);
00286    h->TAttFill::Copy(*this);
00287    h->TAttMarker::Copy(*this);
00288 
00289    std::string gname = "Graph_from_" + std::string(h->GetName() );
00290    SetName(gname.c_str());
00291    SetTitle(h->GetTitle());
00292 }
00293 
00294 
00295 //______________________________________________________________________________
00296 TGraph::TGraph(const TF1 *f, Option_t *option)
00297        : TNamed("Graph","Graph"), TAttLine(), TAttFill(1,1001), TAttMarker()
00298 {
00299    // Graph constructor importing its parameters from the TF1 object passed as argument
00300    // if option =="" (default), a TGraph is created with points computed
00301    //                at the fNpx points of f.
00302    // if option =="d", a TGraph is created with points computed with the derivatives
00303    //                at the fNpx points of f.
00304    // if option =="i", a TGraph is created with points computed with the integral
00305    //                at the fNpx points of f.
00306    // if option =="I", a TGraph is created with points computed with the integral
00307    //                at the fNpx+1 points of f and the integral is normalized to 1.
00308 
00309    char coption = ' ';
00310    if (!f) {
00311       Error("TGraph", "Pointer to function is null");
00312       fNpoints = 0;
00313    } else {
00314       fNpoints   = f->GetNpx();
00315       if (option) coption = *option;
00316       if (coption == 'i' || coption == 'I') fNpoints++;
00317    }
00318    if (!CtorAllocate()) return;
00319 
00320    Double_t xmin = f->GetXmin();
00321    Double_t xmax = f->GetXmax();
00322    Double_t dx   = (xmax-xmin)/fNpoints;
00323    Double_t integ = 0;
00324    Int_t i;
00325    for (i=0;i<fNpoints;i++) {
00326       if (coption == 'i' || coption == 'I') {
00327          fX[i] = xmin +i*dx;
00328          if (i == 0) fY[i] = 0;
00329          else        fY[i] = integ + ((TF1*)f)->Integral(fX[i]-dx,fX[i]);
00330          integ = fY[i];
00331       } else if (coption == 'd' || coption == 'D') {
00332          fX[i] = xmin + (i+0.5)*dx;
00333          fY[i] = ((TF1*)f)->Derivative(fX[i]);
00334       } else {
00335          fX[i] = xmin + (i+0.5)*dx;
00336          fY[i] = ((TF1*)f)->Eval(fX[i]);
00337       }
00338    }
00339    if (integ != 0 && coption == 'I') {
00340       for (i=1;i<fNpoints;i++) fY[i] /= integ;
00341    }
00342 
00343    f->TAttLine::Copy(*this);
00344    f->TAttFill::Copy(*this);
00345    f->TAttMarker::Copy(*this);
00346 
00347    SetName(f->GetName());
00348    SetTitle(f->GetTitle());
00349 }
00350 
00351 
00352 //______________________________________________________________________________
00353 TGraph::TGraph(const char *filename, const char *format, Option_t *)
00354        : TNamed("Graph",filename), TAttLine(), TAttFill(1,1001), TAttMarker()
00355 {
00356    // Graph constructor reading input from filename
00357    // filename is assumed to contain at least two columns of numbers
00358    // the string format is by default "%lg %lg".
00359    // this is a standard c formatting for scanf. If columns of numbers should be skipped,
00360    // a "%*lg" for each column can be added, e.g. "%lg %*lg %lg" would read x-values from
00361    // the first and y-values from the third column.
00362 
00363    Double_t x,y;
00364    TString fname = filename;
00365    gSystem->ExpandPathName(fname);
00366 
00367    ifstream infile(fname.Data());
00368    if(!infile.good()){
00369       MakeZombie();
00370       Error("TGraph", "Cannot open file: %s, TGraph is Zombie",filename);
00371       fNpoints = 0;
00372    } else {
00373       fNpoints = 100;  //initial number of points
00374    }
00375    if (!CtorAllocate()) return;
00376    std::string line;
00377    Int_t np=0;
00378    while(std::getline(infile,line,'\n')){
00379       if(2 != sscanf(line.c_str(),format,&x,&y) ) {
00380          continue; // skip empty and ill-formed lines
00381       }
00382       SetPoint(np,x,y);
00383       np++;
00384    }
00385    Set(np);
00386 }
00387 
00388 
00389 //______________________________________________________________________________
00390 TGraph::~TGraph()
00391 {
00392    // Graph default destructor.
00393 
00394    delete [] fX;
00395    delete [] fY;
00396    if (fFunctions) {
00397       fFunctions->SetBit(kInvalidObject);
00398       //special logic to support the case where the same object is
00399       //added multiple times in fFunctions.
00400       //This case happens when the same object is added with different
00401       //drawing modes
00402       TObject *obj;
00403       while ((obj  = fFunctions->First())) {
00404          while(fFunctions->Remove(obj)) { }
00405          delete obj;
00406       }
00407       delete fFunctions;
00408       fFunctions = 0; //to avoid accessing a deleted object in RecursiveRemove
00409    }
00410    delete fHistogram;
00411 }
00412 
00413 
00414 //______________________________________________________________________________
00415 Double_t** TGraph::AllocateArrays(Int_t Narrays, Int_t arraySize)
00416 {
00417    // Allocate arrays.
00418 
00419    if (arraySize < 0) { arraySize = 0; }
00420    Double_t **newarrays = new Double_t*[Narrays];
00421    if (!arraySize) {
00422       for (Int_t i = 0; i < Narrays; ++i)
00423          newarrays[i] = 0;
00424    } else {
00425       for (Int_t i = 0; i < Narrays; ++i)
00426          newarrays[i] = new Double_t[arraySize];
00427    }
00428    fMaxSize = arraySize;
00429    return newarrays;
00430 }
00431 
00432 
00433 //______________________________________________________________________________
00434 void TGraph::Apply(TF1 *f)
00435 {
00436    // Apply function f to all the data points
00437    // f may be a 1-D function TF1 or 2-d function TF2
00438    // The Y values of the graph are replaced by the new values computed
00439    // using the function
00440 
00441    for (Int_t i=0;i<fNpoints;i++) {
00442       fY[i] = f->Eval(fX[i],fY[i]);
00443    }
00444 }
00445 
00446 
00447 //______________________________________________________________________________
00448 void TGraph::Browse(TBrowser *b)
00449 {
00450    // Browse
00451 
00452    TString opt = gEnv->GetValue("TGraph.BrowseOption","");
00453    if (opt.IsNull()) {
00454       opt = b ? b->GetDrawOption() : "alp";
00455    }
00456    Draw(opt.Data());
00457    gPad->Update();
00458 }
00459 
00460 
00461 //______________________________________________________________________________
00462 Double_t TGraph::Chisquare(const TF1 *f1) const
00463 {
00464    // Return the chisquare of this graph with respect to f1.
00465    // The chisquare is computed as the sum of the quantity below at each point:
00466    // Begin_Latex
00467    // #frac{(y-f1(x))^{2}}{ey^{2}+(#frac{1}{2}(exl+exh)f1'(x))^{2}}
00468    // End_latex
00469    // where x and y are the graph point coordinates and f1'(x) is the derivative of function f1(x).
00470    // This method to approximate the uncertainty in y because of the errors in x, is called
00471    // "effective variance" method.
00472    // In case of a pure TGraph, the denominator is 1.
00473    // In case of a TGraphErrors or TGraphAsymmErrors the errors are taken
00474    // into account.
00475 
00476    if (!f1) return 0;
00477    Double_t cu,eu,exh,exl,ey,eux,fu,fsum;
00478    Double_t x[1];
00479    Double_t chi2 = 0;
00480    TF1 *func = (TF1*)f1; //EvalPar is not const !
00481    for (Int_t i=0;i<fNpoints;i++) {
00482       func->InitArgs(x,0); //must be inside the loop because of TF1::Derivative calling InitArgs
00483       x[0] = fX[i];
00484       if (!func->IsInside(x)) continue;
00485       cu   = fY[i];
00486       TF1::RejectPoint(kFALSE);
00487       fu   = func->EvalPar(x);
00488       if (TF1::RejectedPoint()) continue;
00489       fsum = (cu-fu);
00490       //npfits++;
00491       exh = GetErrorXhigh(i);
00492       exl = GetErrorXlow(i);
00493       if (fsum < 0)
00494          ey = GetErrorYhigh(i);
00495       else
00496          ey = GetErrorYlow(i);
00497       if (exl < 0) exl = 0;
00498       if (exh < 0) exh = 0;
00499       if (ey < 0)  ey  = 0;
00500       if (exh > 0 || exl > 0) {
00501          //"Effective Variance" method introduced by Anna Kreshuk
00502          //a copy of the algorithm in GraphFitChisquare from TFitter
00503          eux = 0.5*(exl + exh)*func->Derivative(x[0]);
00504       } else
00505          eux = 0.;
00506       eu = ey*ey+eux*eux;
00507       if (eu <= 0) eu = 1;
00508       chi2 += fsum*fsum/eu;
00509    }
00510    return chi2;
00511 }
00512 
00513 
00514 //______________________________________________________________________________
00515 Bool_t TGraph::CompareArg(const TGraph* gr, Int_t left, Int_t right)
00516 {
00517    // Return kTRUE if point number "left"'s argument (angle with respect to positive
00518    // x-axis) is bigger than that of point number "right". Can be used by Sort.
00519 
00520    Double_t xl,yl,xr,yr;
00521    gr->GetPoint(left,xl,yl);
00522    gr->GetPoint(right,xr,yr);
00523    return (TMath::ATan2(yl, xl) > TMath::ATan2(yr, xr));
00524 }
00525 
00526 
00527 //______________________________________________________________________________
00528 Bool_t TGraph::CompareX(const TGraph* gr, Int_t left, Int_t right)
00529 {
00530    // Return kTRUE if fX[left] > fX[right]. Can be used by Sort.
00531 
00532    return gr->fX[left]>gr->fX[right];
00533 }
00534 
00535 
00536 //______________________________________________________________________________
00537 Bool_t TGraph::CompareY(const TGraph* gr, Int_t left, Int_t right)
00538 {
00539    // Return kTRUE if fY[left] > fY[right]. Can be used by Sort.
00540 
00541    return gr->fY[left]>gr->fY[right];
00542 }
00543 
00544 
00545 //______________________________________________________________________________
00546 Bool_t TGraph::CompareRadius(const TGraph* gr, Int_t left, Int_t right)
00547 {
00548    // Return kTRUE if point number "left"'s distance to origin is bigger than
00549    // that of point number "right". Can be used by Sort.
00550 
00551    return gr->fX[left]*gr->fX[left]+gr->fY[left]*gr->fY[left]
00552       >gr->fX[right]*gr->fX[right]+gr->fY[right]*gr->fY[right];
00553 }
00554 
00555 
00556 //______________________________________________________________________________
00557 void TGraph::ComputeRange(Double_t &xmin, Double_t &ymin, Double_t &xmax, Double_t &ymax) const
00558 {
00559    // Compute the x/y range of the points in this graph
00560    if (fNpoints <= 0) {
00561       xmin=xmax=ymin=ymax = 0;
00562       return;
00563    }
00564    xmin = xmax = fX[0];
00565    ymin = ymax = fY[0];
00566    for (Int_t i=1;i<fNpoints;i++) {
00567       if (fX[i] < xmin) xmin = fX[i];
00568       if (fX[i] > xmax) xmax = fX[i];
00569       if (fY[i] < ymin) ymin = fY[i];
00570       if (fY[i] > ymax) ymax = fY[i];
00571    }
00572 }
00573 
00574 
00575 //______________________________________________________________________________
00576 void TGraph::CopyAndRelease(Double_t **newarrays, Int_t ibegin, Int_t iend,
00577                            Int_t obegin)
00578 {
00579    // Copy points from fX and fY to arrays[0] and arrays[1]
00580    // or to fX and fY if arrays == 0 and ibegin != iend.
00581    // If newarrays is non null, replace fX, fY with pointers from newarrays[0,1].
00582    // Delete newarrays, old fX and fY
00583 
00584    CopyPoints(newarrays, ibegin, iend, obegin);
00585    if (newarrays) {
00586       delete[] fX;
00587       fX = newarrays[0];
00588       delete[] fY;
00589       fY = newarrays[1];
00590       delete[] newarrays;
00591    }
00592 }
00593 
00594 
00595 //______________________________________________________________________________
00596 Bool_t TGraph::CopyPoints(Double_t **arrays, Int_t ibegin, Int_t iend,
00597                         Int_t obegin)
00598 {
00599    // Copy points from fX and fY to arrays[0] and arrays[1]
00600    // or to fX and fY if arrays == 0 and ibegin != iend.
00601 
00602    if (ibegin < 0 || iend <= ibegin || obegin < 0) { // Error;
00603       return kFALSE;
00604    }
00605    if (!arrays && ibegin == obegin) { // No copying is needed
00606       return kFALSE;
00607    }
00608    Int_t n = (iend - ibegin)*sizeof(Double_t);
00609    if (arrays) {
00610       memmove(&arrays[0][obegin], &fX[ibegin], n);
00611       memmove(&arrays[1][obegin], &fY[ibegin], n);
00612    } else {
00613       memmove(&fX[obegin], &fX[ibegin], n);
00614       memmove(&fY[obegin], &fY[ibegin], n);
00615    }
00616    return kTRUE;
00617 }
00618 
00619 
00620 //______________________________________________________________________________
00621 Bool_t TGraph::CtorAllocate()
00622 {
00623    // In constructors set fNpoints than call this method.
00624    // Return kFALSE if the graph will contain no points.
00625 
00626    fHistogram = 0;
00627    fMaximum = -1111;
00628    fMinimum = -1111;
00629    SetBit(kClipFrame);
00630    fFunctions = new TList;
00631    if (fNpoints <= 0) {
00632       fNpoints = 0;
00633       fMaxSize   = 0;
00634       fX         = 0;
00635       fY         = 0;
00636       return kFALSE;
00637    } else {
00638       fMaxSize   = fNpoints;
00639       fX = new Double_t[fMaxSize];
00640       fY = new Double_t[fMaxSize];
00641    }
00642    return kTRUE;
00643 }
00644 
00645 
00646 //______________________________________________________________________________
00647 void TGraph::Draw(Option_t *option)
00648 {
00649    /* Begin_Html
00650    Draw this graph with its current attributes.
00651    <p>
00652    The options to draw a graph are described in
00653    <a href="http://root.cern.ch/root/html/TGraphPainter.html">TGraphPainter</a>
00654    class.
00655    End_Html */
00656 
00657    TString opt = option;
00658    opt.ToLower();
00659 
00660    if (opt.Contains("same")) {
00661       opt.ReplaceAll("same","");
00662    }
00663 
00664    // in case of option *, set marker style to 3 (star) and replace
00665    // * option by option P.
00666    Ssiz_t pos;
00667    if ((pos = opt.Index("*")) != kNPOS) {
00668       SetMarkerStyle(3);
00669       opt.Replace(pos, 1, "p");
00670    }
00671    if (gPad) {
00672       if (!gPad->IsEditable()) gROOT->MakeDefCanvas();
00673       if (opt.Contains("a")) gPad->Clear();
00674    }
00675    AppendPad(opt);
00676 }
00677 
00678 
00679 //______________________________________________________________________________
00680 Int_t TGraph::DistancetoPrimitive(Int_t px, Int_t py)
00681 {
00682    // Compute distance from point px,py to a graph.
00683    //
00684    //  Compute the closest distance of approach from point px,py to this line.
00685    //  The distance is computed in pixels units.
00686 
00687    return TVirtualGraphPainter::GetPainter()->DistancetoPrimitiveHelper(this, px,py);
00688 }
00689 
00690 
00691 //______________________________________________________________________________
00692 void TGraph::DrawGraph(Int_t n, const Int_t *x, const Int_t *y, Option_t *option)
00693 {
00694    // Draw this graph with new attributes.
00695 
00696    TGraph *newgraph = new TGraph(n, x, y);
00697    TAttLine::Copy(*newgraph);
00698    TAttFill::Copy(*newgraph);
00699    TAttMarker::Copy(*newgraph);
00700    newgraph->SetBit(kCanDelete);
00701    newgraph->AppendPad(option);
00702 }
00703 
00704 
00705 //______________________________________________________________________________
00706 void TGraph::DrawGraph(Int_t n, const Float_t *x, const Float_t *y, Option_t *option)
00707 {
00708    // Draw this graph with new attributes.
00709 
00710    TGraph *newgraph = new TGraph(n, x, y);
00711    TAttLine::Copy(*newgraph);
00712    TAttFill::Copy(*newgraph);
00713    TAttMarker::Copy(*newgraph);
00714    newgraph->SetBit(kCanDelete);
00715    newgraph->AppendPad(option);
00716 }
00717 
00718 
00719 //______________________________________________________________________________
00720 void TGraph::DrawGraph(Int_t n, const Double_t *x, const Double_t *y, Option_t *option)
00721 {
00722    // Draw this graph with new attributes.
00723 
00724    const Double_t *xx = x;
00725    const Double_t *yy = y;
00726    if (!xx) xx = fX;
00727    if (!yy) yy = fY;
00728    TGraph *newgraph = new TGraph(n, xx, yy);
00729    TAttLine::Copy(*newgraph);
00730    TAttFill::Copy(*newgraph);
00731    TAttMarker::Copy(*newgraph);
00732    newgraph->SetBit(kCanDelete);
00733    newgraph->AppendPad(option);
00734 }
00735 
00736 
00737 //______________________________________________________________________________
00738 void TGraph::DrawPanel()
00739 {
00740    // Display a panel with all graph drawing options.
00741 
00742    TVirtualGraphPainter::GetPainter()->DrawPanelHelper(this);
00743 }
00744 
00745 
00746 //______________________________________________________________________________
00747 Double_t TGraph::Eval(Double_t x, TSpline *spline, Option_t *option) const
00748 {
00749    // Interpolate points in this graph at x using a TSpline
00750    //  -if spline==0 and option="" a linear interpolation between the two points
00751    //   close to x is computed. If x is outside the graph range, a linear
00752    //   extrapolation is computed.
00753    //  -if spline==0 and option="S" a TSpline3 object is created using this graph
00754    //   and the interpolated value from the spline is returned.
00755    //   the internally created spline is deleted on return.
00756    //  -if spline is specified, it is used to return the interpolated value.
00757 
00758 
00759    if (!spline) {
00760 
00761       if (fNpoints == 0) return 0;
00762       if (fNpoints == 1) return fY[0];
00763 
00764 
00765       TString opt = option;
00766       opt.ToLower();
00767       if (opt.Contains("s")) {
00768 
00769          // points must be sorted before using a TSpline
00770          std::vector<Double_t> xsort(fNpoints);
00771          std::vector<Double_t> ysort(fNpoints);
00772          std::vector<Int_t> indxsort(fNpoints);
00773          TMath::Sort(fNpoints, fX, &indxsort[0], false );
00774          for (Int_t i = 0; i < fNpoints; ++i) {
00775             xsort[i] = fX[ indxsort[i] ];
00776             ysort[i] = fY[ indxsort[i] ];
00777          }
00778 
00779          // spline interpolation creating a new spline
00780          TSpline3 *s = new TSpline3("",&xsort[0], &ysort[0], fNpoints);
00781          Double_t result = s->Eval(x);
00782          delete s;
00783          return result;
00784       }
00785       //linear interpolation
00786       //In case x is < fX[0] or > fX[fNpoints-1] return the extrapolated point
00787 
00788       //find points in graph around x assuming points are not sorted
00789       // (if point are sorted could use binary search)
00790 
00791       // find neighbours simply looping  all points
00792       // and find also the 2 adjacent points: (low2 < low < x < up < up2 )
00793       // needed in case x is outside the graph ascissa interval
00794       Int_t low  = -1;  Int_t up  = -1;
00795       Int_t low2 = -1;  Int_t up2 = -1;
00796 
00797       for (Int_t i = 0; i < fNpoints; ++i) {
00798          if ( fX[i] < x ) {
00799             if  (low == -1 || fX[i] > fX[low] )  {  low2 = low;   low = i; }
00800             else if ( low2 == -1  ) low2 = i;
00801          }
00802          else if ( fX[i] > x) {
00803             if (up  == -1 || fX[i] < fX[up]  )  {  up2 = up;     up = i;  }
00804             else if (up2 == -1) up2 = i;
00805          }
00806          else // case x == fX[i]
00807             return fY[i]; // no interpolation needed
00808       }
00809 
00810       // treat cases when x is outside graph min max abscissa
00811       if (up == -1)  {up  = low; low = low2;}
00812       if (low == -1) {low = up;  up  = up2;  }
00813 
00814       assert( low != -1 && up != -1);
00815 
00816       if (fX[low] == fX[up]) return fY[low];
00817       Double_t yn = fY[up] + (x - fX[up] ) * (fY[low]-fY[up] ) / ( fX[low] - fX[up] );
00818       return yn;
00819    } else {
00820       //spline interpolation using the input spline
00821       return spline->Eval(x);
00822    }
00823 }
00824 
00825 
00826 //______________________________________________________________________________
00827 void TGraph::ExecuteEvent(Int_t event, Int_t px, Int_t py)
00828 {
00829    // Execute action corresponding to one event.
00830    //
00831    //  This member function is called when a graph is clicked with the locator
00832    //
00833    //  If Left button clicked on one of the line end points, this point
00834    //     follows the cursor until button is released.
00835    //
00836    //  if Middle button clicked, the line is moved parallel to itself
00837    //     until the button is released.
00838 
00839    TVirtualGraphPainter::GetPainter()->ExecuteEventHelper(this, event, px, py);
00840 }
00841 
00842 
00843 //______________________________________________________________________________
00844 void TGraph::Expand(Int_t newsize)
00845 {
00846    // If array sizes <= newsize, expand storage to 2*newsize.
00847 
00848    Double_t **ps = ExpandAndCopy(newsize, fNpoints);
00849    CopyAndRelease(ps, 0, 0, 0);
00850 }
00851 
00852 
00853 //______________________________________________________________________________
00854 void TGraph::Expand(Int_t newsize, Int_t step)
00855 {
00856    // If graph capacity is less than newsize points then make array sizes
00857    // equal to least multiple of step to contain newsize points.
00858    // Returns kTRUE if size was altered
00859 
00860    if (newsize <= fMaxSize) {
00861       return;
00862    }
00863    Double_t **ps = Allocate(step*(newsize/step + (newsize%step?1:0)));
00864    CopyAndRelease(ps, 0, fNpoints, 0);
00865 }
00866 
00867 
00868 //______________________________________________________________________________
00869 Double_t **TGraph::ExpandAndCopy(Int_t size, Int_t iend)
00870 {
00871    // if size > fMaxSize allocate new arrays of 2*size points
00872    //  and copy oend first points.
00873    // Return pointer to new arrays.
00874 
00875    if (size <= fMaxSize) { return 0; }
00876    Double_t **newarrays = Allocate(2*size);
00877    CopyPoints(newarrays, 0, iend, 0);
00878    return newarrays;
00879 }
00880 
00881 
00882 //______________________________________________________________________________
00883 void TGraph::FillZero(Int_t begin, Int_t end, Bool_t)
00884 {
00885    // Set zero values for point arrays in the range [begin, end)
00886    // Should be redefined in descendant classes
00887 
00888    memset(fX + begin, 0, (end - begin)*sizeof(Double_t));
00889    memset(fY + begin, 0, (end - begin)*sizeof(Double_t));
00890 }
00891 
00892 
00893 //______________________________________________________________________________
00894 TObject *TGraph::FindObject(const char *name) const
00895 {
00896    // Search object named name in the list of functions
00897 
00898    if (fFunctions) return fFunctions->FindObject(name);
00899    return 0;
00900 }
00901 
00902 
00903 //______________________________________________________________________________
00904 TObject *TGraph::FindObject(const TObject *obj) const
00905 {
00906    // Search object obj in the list of functions
00907 
00908    if (fFunctions) return fFunctions->FindObject(obj);
00909    return 0;
00910 }
00911 
00912 
00913 //______________________________________________________________________________
00914 TFitResultPtr TGraph::Fit(const char *fname, Option_t *option, Option_t *, Axis_t xmin, Axis_t xmax)
00915 {
00916    // Fit this graph with function with name fname.
00917    //
00918    //  interface to TGraph::Fit(TF1 *f1...
00919    //
00920    //      fname is the name of an already predefined function created by TF1 or TF2
00921    //      Predefined functions such as gaus, expo and poln are automatically
00922    //      created by ROOT.
00923    //      fname can also be a formula, accepted by the linear fitter (linear parts divided
00924    //      by "++" sign), for example "x++sin(x)" for fitting "[0]*x+[1]*sin(x)"
00925 
00926    char *linear;
00927    linear= (char*) strstr(fname, "++");
00928    TF1 *f1=0;
00929    if (linear)
00930       f1=new TF1(fname, fname, xmin, xmax);
00931    else {
00932       f1 = (TF1*)gROOT->GetFunction(fname);
00933       if (!f1) { Printf("Unknown function: %s",fname); return -1; }
00934    }
00935    return Fit(f1,option,"",xmin,xmax);
00936 }
00937 
00938 
00939 //______________________________________________________________________________
00940 TFitResultPtr TGraph::Fit(TF1 *f1, Option_t *option, Option_t *goption, Axis_t rxmin, Axis_t rxmax)
00941 {
00942    // Fit this graph with function f1.
00943    //
00944    //   f1 is an already predefined function created by TF1.
00945    //   Predefined functions such as gaus, expo and poln are automatically
00946    //   created by ROOT.
00947    //
00948    //   The list of fit options is given in parameter option.
00949    //      option = "W" Set all weights to 1; ignore error bars
00950    //             = "U" Use a User specified fitting algorithm (via SetFCN)
00951    //             = "Q" Quiet mode (minimum printing)
00952    //             = "V" Verbose mode (default is between Q and V)
00953    //             = "E"  Perform better Errors estimation using Minos technique
00954    //             = "B"  User defined parameter settings are used for predefined functions
00955    //                    like "gaus", "expo", "poln", "landau".
00956    //                    Use this option when you want to fix one or more parameters for these functions.
00957    //             = "M"  More. Improve fit results.
00958    //                    It uses the IMPROVE command of TMinuit (see TMinuit::mnimpr)
00959    //                    This algorithm attempts to improve the found local minimum by
00960    //                    searching for a better one.
00961    //             = "R" Use the Range specified in the function range
00962    //             = "N" Do not store the graphics function, do not draw
00963    //             = "0" Do not plot the result of the fit. By default the fitted function
00964    //                   is drawn unless the option "N" above is specified.
00965    //             = "+" Add this new fitted function to the list of fitted functions
00966    //                   (by default, any previous function is deleted)
00967    //             = "C" In case of linear fitting, do not calculate the chisquare
00968    //                    (saves time)
00969    //             = "F" If fitting a polN, use the minuit fitter
00970    //             = "EX0" When fitting a TGraphErrors do not consider errors in the coordinate
00971    //             = "ROB" In case of linear fitting, compute the LTS regression
00972    //                     coefficients (robust (resistant) regression), using
00973    //                     the default fraction of good points
00974    //               "ROB=0.x" - compute the LTS regression coefficients, using
00975    //                           0.x as a fraction of good points
00976    //             = "S"  The result of the fit is returned in the TFitResultPtr
00977    //                     (see below Access to the Fit Result)
00978    //
00979    //   When the fit is drawn (by default), the parameter goption may be used
00980    //   to specify a list of graphics options. See TGraphPainter for a complete
00981    //   list of these options.
00982    //
00983    //   In order to use the Range option, one must first create a function
00984    //   with the expression to be fitted. For example, if your graph
00985    //   has a defined range between -4 and 4 and you want to fit a gaussian
00986    //   only in the interval 1 to 3, you can do:
00987    //        TF1 *f1 = new TF1("f1","gaus",1,3);
00988    //        graph->Fit("f1","R");
00989    //
00990    //
00991    // Who is calling this function:
00992    //
00993    //   Note that this function is called when calling TGraphErrors::Fit
00994    //   or TGraphAsymmErrors::Fit ot TGraphBentErrors::Fit
00995    //   See the discussion below on error calulation.
00996    //
00997    // Linear fitting:
00998    // ===============
00999    //
01000    //   When the fitting function is linear (contains the "++" sign) or the fitting
01001    //   function is a polynomial, a linear fitter is initialised.
01002    //   To create a linear function, use the following syntax: linear parts
01003    //   separated by "++" sign.
01004    //   Example: to fit the parameters of "[0]*x + [1]*sin(x)", create a
01005    //    TF1 *f1=new TF1("f1", "x++sin(x)", xmin, xmax);
01006    //   For such a TF1 you don't have to set the initial conditions.
01007    //   Going via the linear fitter for functions, linear in parameters, gives a
01008    //   considerable advantage in speed.
01009    //
01010    // Setting initial conditions:
01011    // ===========================
01012    //
01013    //   Parameters must be initialized before invoking the Fit function.
01014    //   The setting of the parameter initial values is automatic for the
01015    //   predefined functions : poln, expo, gaus, landau. One can however disable
01016    //   this automatic computation by specifying the option "B".
01017    //   You can specify boundary limits for some or all parameters via
01018    //        f1->SetParLimits(p_number, parmin, parmax);
01019    //   If parmin>=parmax, the parameter is fixed
01020    //   Note that you are not forced to fix the limits for all parameters.
01021    //   For example, if you fit a function with 6 parameters, you can do:
01022    //     func->SetParameters(0,3.1,1.e-6,0.1,-8,100);
01023    //     func->SetParLimits(4,-10,-4);
01024    //     func->SetParLimits(5, 1,1);
01025    //   With this setup, parameters 0->3 can vary freely.
01026    //   Parameter 4 has boundaries [-10,-4] with initial value -8.
01027    //   Parameter 5 is fixed to 100.
01028    //
01029    // Fit range:
01030    // ==========
01031    //
01032    //   The fit range can be specified in two ways:
01033    //     - specify rxmax > rxmin (default is rxmin=rxmax=0)
01034    //     - specify the option "R". In this case, the function will be taken
01035    //       instead of the full graph range.
01036    //
01037    // Changing the fitting function:
01038    // ==============================
01039    //
01040    //   By default a chi2 fitting function is used for fitting a TGraph.
01041    //   The function is implemented in FitUtil::EvaluateChi2.
01042    //   In case of TGraphErrors an effective chi2 is used (see below TGraphErrors fit)
01043    //   To specify a User defined fitting function, specify option "U" and
01044    //   call the following functions:
01045    //     TVirtualFitter::Fitter(mygraph)->SetFCN(MyFittingFunction)
01046    //   where MyFittingFunction is of type:
01047    //   extern void MyFittingFunction(Int_t &npar, Double_t *gin, Double_t &f,
01048    //                                 Double_t *u, Int_t flag);
01049    //
01050    //
01051    // TGraphErrors fit:
01052    // =================
01053    //
01054    //   In case of a TGraphErrors object, when x errors are present, the error along x,
01055    //   is projected along the y-direction by calculating the function at the points x-exlow and
01056    //   x+exhigh. The chisquare is then computed as the sum of the quantity below at each point:
01057    //
01058    // Begin_Latex
01059    // #frac{(y-f(x))^{2}}{ey^{2}+(#frac{1}{2}(exl+exh)f'(x))^{2}}
01060    // End_Latex
01061    //
01062    //   where x and y are the point coordinates, and f'(x) is the derivative of the
01063    //   function f(x).
01064    //
01065    //   In case the function lies below (above) the data point, ey is ey_low (ey_high).
01066    //
01067    //   thanks to Andy Haas (haas@yahoo.com) for adding the case with TGraphAsymmErrors
01068    //             University of Washington
01069    //
01070    //   The approach used to approximate the uncertainty in y because of the
01071    //   errors in x is to make it equal the error in x times the slope of the line.
01072    //   The improvement, compared to the first method (f(x+ exhigh) - f(x-exlow))/2
01073    //   is of (error of x)**2 order. This approach is called "effective variance method".
01074    //   This improvement has been made in version 4.00/08 by Anna Kreshuk.
01075    //   The implementation is provided in the function FitUtil::EvaluateChi2Effective
01076    //
01077    // NOTE:
01078    //   1) By using the "effective variance" method a simple linear regression
01079    //      becomes a non-linear case, which takes several iterations
01080    //      instead of 0 as in the linear case.
01081    //
01082    //   2) The effective variance technique assumes that there is no correlation
01083    //      between the x and y coordinate.
01084    //
01085    //   3) The standard chi2 (least square) method without error in the coordinates (x) can
01086    //       be forced by using option "EX0"
01087    //
01088    //   4)  The linear fitter doesn't take into account the errors in x. When fitting a
01089    //       TGraphErrors with a linear functions the errors in x willnot be considere.
01090    //        If errors in x are important, go through minuit (use option "F" for polynomial fitting).
01091    //
01092    //   5) When fitting a TGraph (i.e. no errors associated with each point),
01093    //   a correction is applied to the errors on the parameters with the following
01094    //   formula:
01095    //      errorp *= sqrt(chisquare/(ndf-1))
01096    //
01097    //   Access to the fit result
01098    //   ========================
01099    //  The function returns a TFitResultPtr which can hold a  pointer to a TFitResult object.
01100    //  By default the TFitResultPtr contains only the status of the fit which is return by an
01101    //  automatic conversion of the TFitResultPtr to an integer. One can write in this case
01102    //  directly:
01103    //  Int_t fitStatus =  h->Fit(myFunc)
01104    //
01105    //  If the option "S" is instead used, TFitResultPtr contains the TFitResult and behaves
01106    //  as a smart pointer to it. For example one can do:
01107    //  TFitResultPtr r = h->Fit(myFunc,"S");
01108    //  TMatrixDSym cov = r->GetCovarianceMatrix();  //  to access the covariance matrix
01109    //  Double_t chi2   = r->Chi2(); // to retrieve the fit chi2
01110    //  Double_t par0   = r->Value(0); // retrieve the value for the parameter 0
01111    //  Double_t err0   = r->Error(0); // retrieve the error for the parameter 0
01112    //  r->Print("V");     // print full information of fit including covariance matrix
01113    //  r->Write();        // store the result in a file
01114    //
01115    //  The fit parameters, error and chi2 (but not covariance matrix) can be retrieved also
01116    //  from the fitted function.
01117    //  If the histogram is made persistent, the list of
01118    //  associated functions is also persistent. Given a pointer (see above)
01119    //  to an associated function myfunc, one can retrieve the function/fit
01120    //  parameters with calls such as:
01121    //    Double_t chi2 = myfunc->GetChisquare();
01122    //    Double_t par0 = myfunc->GetParameter(0); //value of 1st parameter
01123    //    Double_t err0 = myfunc->GetParError(0);  //error on first parameter
01124    //
01125    //
01126    //  Access to the fit status
01127    //  =====================
01128    //  The status of the fit can be obtained converting the TFitResultPtr to an integer
01129    //  indipendently if the fit option "S" is used or not:
01130    //  TFitResultPtr r = h=>Fit(myFunc,opt);
01131    //  Int_t fitStatus = r;
01132    //
01133    //  The fitStatus is 0 if the fit is OK (i.e. no error occurred).
01134    //  The value of the fit status code is negative in case of an error not connected with the
01135    //  minimization procedure, for example when a wrong function is used.
01136    //  Otherwise the return value is the one returned from the minimization procedure.
01137    //  When TMinuit (default case) or Minuit2 are used as minimizer the status returned is :
01138    //  fitStatus =  migradResult + 10*minosResult + 100*hesseResult + 1000*improveResult.
01139    //  TMinuit will return 0 (for migrad, minos, hesse or improve) in case of success and 4 in
01140    //  case of error (see the documentation of TMinuit::mnexcm). So for example, for an error
01141    //  only in Minos but not in Migrad a fitStatus of 40 will be returned.
01142    //  Minuit2 will return also 0 in case of success and different values in migrad, minos or
01143    //  hesse depending on the error.   See in this case the documentation of
01144    //  Minuit2Minimizer::Minimize for the migradResult, Minuit2Minimizer::GetMinosError for the
01145    //  minosResult and Minuit2Minimizer::Hesse for the hesseResult.
01146    //  If other minimizers are used see their specific documentation for the status code
01147    //  returned. For example in the case of Fumili, for the status returned see TFumili::Minimize.
01148    //
01149    // Associated functions:
01150    // =====================
01151    //
01152    //   One or more object (typically a TF1*) can be added to the list
01153    //   of functions (fFunctions) associated with each graph.
01154    //   When TGraph::Fit is invoked, the fitted function is added to this list.
01155    //   Given a graph gr, one can retrieve an associated function
01156    //   with:  TF1 *myfunc = gr->GetFunction("myfunc");
01157    //
01158    //   If the graph is made persistent, the list of associated functions is also
01159    //   persistent. Given a pointer (see above) to an associated function myfunc,
01160    //   one can retrieve the function/fit parameters with calls such as:
01161    //     Double_t chi2 = myfunc->GetChisquare();
01162    //     Double_t par0 = myfunc->GetParameter(0); //value of 1st parameter
01163    //     Double_t err0 = myfunc->GetParError(0);  //error on first parameter
01164    //
01165    // Fit Statistics
01166    // ==============
01167    //
01168    //   You can change the statistics box to display the fit parameters with
01169    //   the TStyle::SetOptFit(mode) method. This mode has four digits.
01170    //   mode = pcev  (default = 0111)
01171    //     v = 1;  print name/values of parameters
01172    //     e = 1;  print errors (if e=1, v must be 1)
01173    //     c = 1;  print Chisquare/Number of degress of freedom
01174    //     p = 1;  print Probability
01175    //
01176    //   For example: gStyle->SetOptFit(1011);
01177    //   prints the fit probability, parameter names/values, and errors.
01178    //   You can change the position of the statistics box with these lines
01179    //   (where g is a pointer to the TGraph):
01180    //
01181    //   Root > TPaveStats *st = (TPaveStats*)g->GetListOfFunctions()->FindObject("stats")
01182    //   Root > st->SetX1NDC(newx1); //new x start position
01183    //   Root > st->SetX2NDC(newx2); //new x end position
01184    //
01185 
01186    Foption_t fitOption;
01187    ROOT::Fit::FitOptionsMake(option,fitOption);
01188    // create range and minimizer options with default values
01189    ROOT::Fit::DataRange range(rxmin,rxmax);
01190    ROOT::Math::MinimizerOptions minOption;
01191    return ROOT::Fit::FitObject(this, f1 , fitOption , minOption, goption, range);
01192 }
01193 
01194 
01195 //______________________________________________________________________________
01196 void TGraph::FitPanel()
01197 {
01198    // Display a GUI panel with all graph fit options.
01199    //
01200    //   See class TFitEditor for example
01201 
01202    if (!gPad)
01203       gROOT->MakeDefCanvas();
01204 
01205    if (!gPad) {
01206       Error("FitPanel", "Unable to create a default canvas");
01207       return;
01208    }
01209 
01210    // use plugin manager to create instance of TFitEditor
01211    TPluginHandler *handler = gROOT->GetPluginManager()->FindHandler("TFitEditor");
01212    if (handler && handler->LoadPlugin() != -1) {
01213       if (handler->ExecPlugin(2, gPad, this) == 0)
01214          Error("FitPanel", "Unable to crate the FitPanel");
01215    }
01216    else
01217          Error("FitPanel", "Unable to find the FitPanel plug-in");
01218 }
01219 
01220 
01221 //______________________________________________________________________________
01222 Double_t TGraph::GetCorrelationFactor() const
01223 {
01224    // Return graph correlation factor
01225 
01226    Double_t rms1 = GetRMS(1);
01227    if (rms1 == 0) return 0;
01228    Double_t rms2 = GetRMS(2);
01229    if (rms2 == 0) return 0;
01230    return GetCovariance()/rms1/rms2;
01231 }
01232 
01233 
01234 //______________________________________________________________________________
01235 Double_t TGraph::GetCovariance() const
01236 {
01237    // Return covariance of vectors x,y
01238 
01239    if (fNpoints <= 0) return 0;
01240    Double_t sum = fNpoints, sumx = 0, sumy = 0, sumxy = 0;
01241 
01242    for (Int_t i=0;i<fNpoints;i++) {
01243       sumx  += fX[i];
01244       sumy  += fY[i];
01245       sumxy += fX[i]*fY[i];
01246    }
01247    return sumxy/sum - sumx/sum*sumy/sum;
01248 }
01249 
01250 
01251 //______________________________________________________________________________
01252 Double_t TGraph::GetMean(Int_t axis) const
01253 {
01254    // Return mean value of X (axis=1)  or Y (axis=2)
01255 
01256    if (axis < 1 || axis > 2) return 0;
01257    if (fNpoints <= 0) return 0;
01258    Double_t sumx = 0;
01259    for (Int_t i=0;i<fNpoints;i++) {
01260       if (axis == 1) sumx += fX[i];
01261       else           sumx += fY[i];
01262    }
01263    return sumx/fNpoints;
01264 }
01265 
01266 
01267 //______________________________________________________________________________
01268 Double_t TGraph::GetRMS(Int_t axis) const
01269 {
01270    // Return RMS of X (axis=1)  or Y (axis=2)
01271 
01272    if (axis < 1 || axis > 2) return 0;
01273    if (fNpoints <= 0) return 0;
01274    Double_t sumx = 0, sumx2 = 0;
01275    for (Int_t i=0;i<fNpoints;i++) {
01276       if (axis == 1) {sumx += fX[i]; sumx2 += fX[i]*fX[i];}
01277       else           {sumx += fY[i]; sumx2 += fY[i]*fY[i];}
01278    }
01279    Double_t x = sumx/fNpoints;
01280    Double_t rms2 = TMath::Abs(sumx2/fNpoints -x*x);
01281    return TMath::Sqrt(rms2);
01282 }
01283 
01284 
01285 //______________________________________________________________________________
01286 Double_t TGraph::GetErrorX(Int_t) const
01287 {
01288    // This function is called by GraphFitChisquare.
01289    // It always returns a negative value. Real implementation in TGraphErrors
01290 
01291    return -1;
01292 }
01293 
01294 
01295 //______________________________________________________________________________
01296 Double_t TGraph::GetErrorY(Int_t) const
01297 {
01298    // This function is called by GraphFitChisquare.
01299    // It always returns a negative value. Real implementation in TGraphErrors
01300 
01301    return -1;
01302 }
01303 
01304 
01305 //______________________________________________________________________________
01306 Double_t TGraph::GetErrorXhigh(Int_t) const
01307 {
01308    // This function is called by GraphFitChisquare.
01309    // It always returns a negative value. Real implementation in TGraphErrors
01310    // and TGraphAsymmErrors
01311 
01312    return -1;
01313 }
01314 
01315 
01316 //______________________________________________________________________________
01317 Double_t TGraph::GetErrorXlow(Int_t) const
01318 {
01319    // This function is called by GraphFitChisquare.
01320    // It always returns a negative value. Real implementation in TGraphErrors
01321    // and TGraphAsymmErrors
01322 
01323    return -1;
01324 }
01325 
01326 
01327 //______________________________________________________________________________
01328 Double_t TGraph::GetErrorYhigh(Int_t) const
01329 {
01330    // This function is called by GraphFitChisquare.
01331    // It always returns a negative value. Real implementation in TGraphErrors
01332    // and TGraphAsymmErrors
01333 
01334    return -1;
01335 }
01336 
01337 
01338 //______________________________________________________________________________
01339 Double_t TGraph::GetErrorYlow(Int_t) const
01340 {
01341    // This function is called by GraphFitChisquare.
01342    // It always returns a negative value. Real implementation in TGraphErrors
01343    // and TGraphAsymmErrors
01344 
01345    return -1;
01346 }
01347 
01348 
01349 //______________________________________________________________________________
01350 TF1 *TGraph::GetFunction(const char *name) const
01351 {
01352    // Return pointer to function with name.
01353    //
01354    // Functions such as TGraph::Fit store the fitted function in the list of
01355    // functions of this graph.
01356 
01357    if (!fFunctions) return 0;
01358    return (TF1*)fFunctions->FindObject(name);
01359 }
01360 
01361 
01362 //______________________________________________________________________________
01363 TH1F *TGraph::GetHistogram() const
01364 {
01365    // Returns a pointer to the histogram used to draw the axis
01366    // Takes into account the two following cases.
01367    //    1- option 'A' was specified in TGraph::Draw. Return fHistogram
01368    //    2- user had called TPad::DrawFrame. return pointer to hframe histogram
01369 
01370    Double_t rwxmin,rwxmax, rwymin, rwymax, maximum, minimum, dx, dy;
01371    Double_t uxmin, uxmax;
01372 
01373    ComputeRange(rwxmin, rwymin, rwxmax, rwymax);  //this is redefined in TGraphErrors
01374 
01375    // (if fHistogram exist) && (if the log scale is on) &&
01376    // (if the computed range minimum is > 0) && (if the fHistogram minimum is zero)
01377    // then it means fHistogram limits have been computed in linear scale
01378    // therefore they might be too strict and cut some points. In that case the
01379    // fHistogram limits should be recomputed ie: the existing fHistogram
01380    // should not be returned.
01381    TH1F *historg = 0;
01382    if (fHistogram) {
01383       if (gPad && gPad->GetLogx()) {
01384          if (rwxmin <= 0 || fHistogram->GetXaxis()->GetXmin() != 0) return fHistogram;
01385       } else if (gPad && gPad->GetLogy()) {
01386          if (rwymin <= 0 || fHistogram->GetMinimum() != 0) return fHistogram;
01387       } else {
01388          return fHistogram;
01389       }
01390       historg = fHistogram;
01391    }
01392 
01393    if (rwxmin == rwxmax) rwxmax += 1.;
01394    if (rwymin == rwymax) rwymax += 1.;
01395    dx = 0.1*(rwxmax-rwxmin);
01396    dy = 0.1*(rwymax-rwymin);
01397    uxmin    = rwxmin - dx;
01398    uxmax    = rwxmax + dx;
01399    minimum  = rwymin - dy;
01400    maximum  = rwymax + dy;
01401    if (fMinimum != -1111) minimum = fMinimum;
01402    if (fMaximum != -1111) maximum = fMaximum;
01403 
01404    // the graph is created with at least as many channels as there are points
01405    // to permit zooming on the full range
01406    if (uxmin < 0 && rwxmin >= 0) {
01407       if (gPad && gPad->GetLogx()) uxmin = 0.9*rwxmin;
01408       else                 uxmin = 0;
01409    }
01410    if (uxmax > 0 && rwxmax <= 0) {
01411       if (gPad && gPad->GetLogx()) uxmax = 1.1*rwxmax;
01412       else                 uxmax = 0;
01413    }
01414    if (minimum < 0 && rwymin >= 0) {
01415       if(gPad && gPad->GetLogy()) minimum = 0.9*rwymin;
01416       else                minimum = 0;
01417    }
01418    if (minimum <= 0 && gPad && gPad->GetLogy()) minimum = 0.001*maximum;
01419    if (uxmin <= 0 && gPad && gPad->GetLogx()) {
01420       if (uxmax > 1000) uxmin = 1;
01421       else              uxmin = 0.001*uxmax;
01422    }
01423 
01424    rwxmin = uxmin;
01425    rwxmax = uxmax;
01426    Int_t npt = 100;
01427    if (fNpoints > npt) npt = fNpoints;
01428    const char *gname = GetName();
01429    if (strlen(gname) == 0) gname = "Graph";
01430    ((TGraph*)this)->fHistogram = new TH1F(gname,GetTitle(),npt,rwxmin,rwxmax);
01431    if (!fHistogram) return 0;
01432    fHistogram->SetMinimum(minimum);
01433    fHistogram->SetBit(TH1::kNoStats);
01434    fHistogram->SetMaximum(maximum);
01435    fHistogram->GetYaxis()->SetLimits(minimum,maximum);
01436    fHistogram->SetDirectory(0);
01437    // Restore the axis attributes if needed
01438    if (historg) {
01439       fHistogram->GetXaxis()->SetTitle(historg->GetXaxis()->GetTitle());
01440       fHistogram->GetXaxis()->CenterTitle(historg->GetXaxis()->GetCenterTitle());
01441       fHistogram->GetXaxis()->RotateTitle(historg->GetXaxis()->GetRotateTitle());
01442       fHistogram->GetXaxis()->SetNoExponent(historg->GetXaxis()->GetNoExponent());
01443       fHistogram->GetXaxis()->SetNdivisions(historg->GetXaxis()->GetNdivisions());
01444       fHistogram->GetXaxis()->SetLabelFont(historg->GetXaxis()->GetLabelFont());
01445       fHistogram->GetXaxis()->SetLabelOffset(historg->GetXaxis()->GetLabelOffset());
01446       fHistogram->GetXaxis()->SetLabelSize(historg->GetXaxis()->GetLabelSize());
01447       fHistogram->GetXaxis()->SetTitleSize(historg->GetXaxis()->GetTitleSize());
01448       fHistogram->GetXaxis()->SetTitleOffset(historg->GetXaxis()->GetTitleOffset());
01449       fHistogram->GetXaxis()->SetTitleFont(historg->GetXaxis()->GetTitleFont());
01450 
01451       fHistogram->GetYaxis()->SetTitle(historg->GetYaxis()->GetTitle());
01452       fHistogram->GetYaxis()->CenterTitle(historg->GetYaxis()->GetCenterTitle());
01453       fHistogram->GetYaxis()->RotateTitle(historg->GetYaxis()->GetRotateTitle());
01454       fHistogram->GetYaxis()->SetNoExponent(historg->GetYaxis()->GetNoExponent());
01455       fHistogram->GetYaxis()->SetNdivisions(historg->GetYaxis()->GetNdivisions());
01456       fHistogram->GetYaxis()->SetLabelFont(historg->GetYaxis()->GetLabelFont());
01457       fHistogram->GetYaxis()->SetLabelOffset(historg->GetYaxis()->GetLabelOffset());
01458       fHistogram->GetYaxis()->SetLabelSize(historg->GetYaxis()->GetLabelSize());
01459       fHistogram->GetYaxis()->SetTitleSize(historg->GetYaxis()->GetTitleSize());
01460       fHistogram->GetYaxis()->SetTitleOffset(historg->GetYaxis()->GetTitleOffset());
01461       fHistogram->GetYaxis()->SetTitleFont(historg->GetYaxis()->GetTitleFont());
01462 
01463       delete historg;
01464    }
01465    return fHistogram;
01466 }
01467 
01468 
01469 //______________________________________________________________________________
01470 Int_t TGraph::GetPoint(Int_t i, Double_t &x, Double_t &y) const
01471 {
01472    // Get x and y values for point number i.
01473    // The function returns -1 in case of an invalid request or the point number otherwise
01474 
01475    if (i < 0 || i >= fNpoints) return -1;
01476    if (!fX || !fY) return -1;
01477    x = fX[i];
01478    y = fY[i];
01479    return i;
01480 }
01481 
01482 
01483 //______________________________________________________________________________
01484 TAxis *TGraph::GetXaxis() const
01485 {
01486    // Get x axis of the graph.
01487 
01488    TH1 *h = GetHistogram();
01489    if (!h) return 0;
01490    return h->GetXaxis();
01491 }
01492 
01493 
01494 //______________________________________________________________________________
01495 TAxis *TGraph::GetYaxis() const
01496 {
01497    // Get y axis of the graph.
01498 
01499    TH1 *h = GetHistogram();
01500    if (!h) return 0;
01501    return h->GetYaxis();
01502 }
01503 
01504 
01505 //______________________________________________________________________________
01506 void TGraph::InitGaus(Double_t xmin, Double_t xmax)
01507 {
01508    // Compute Initial values of parameters for a gaussian.
01509 
01510    Double_t allcha, sumx, sumx2, x, val, rms, mean;
01511    Int_t bin;
01512    const Double_t sqrtpi = 2.506628;
01513 
01514    // Compute mean value and RMS of the graph in the given range
01515    if (xmax <= xmin) {xmin = fX[0]; xmax = fX[fNpoints-1];}
01516    Int_t np = 0;
01517    allcha = sumx = sumx2 = 0;
01518    for (bin=0;bin<fNpoints;bin++) {
01519       x       = fX[bin];
01520       if (x < xmin || x > xmax) continue;
01521       np++;
01522       val     = fY[bin];
01523       sumx   += val*x;
01524       sumx2  += val*x*x;
01525       allcha += val;
01526    }
01527    if (np == 0 || allcha == 0) return;
01528    mean = sumx/allcha;
01529    rms  = TMath::Sqrt(sumx2/allcha - mean*mean);
01530    Double_t binwidx = TMath::Abs((xmax-xmin)/np);
01531    if (rms == 0) rms = 1;
01532    TVirtualFitter *grFitter = TVirtualFitter::GetFitter();
01533    TF1 *f1 = (TF1*)grFitter->GetUserFunc();
01534    f1->SetParameter(0,binwidx*allcha/(sqrtpi*rms));
01535    f1->SetParameter(1,mean);
01536    f1->SetParameter(2,rms);
01537    f1->SetParLimits(2,0,10*rms);
01538 }
01539 
01540 
01541 //______________________________________________________________________________
01542 void TGraph::InitExpo(Double_t xmin, Double_t xmax)
01543 {
01544    // Compute Initial values of parameters for an exponential.
01545 
01546    Double_t constant, slope;
01547    Int_t ifail;
01548    if (xmax <= xmin) {xmin = fX[0]; xmax = fX[fNpoints-1];}
01549    Int_t nchanx = fNpoints;
01550 
01551    LeastSquareLinearFit(-nchanx, constant, slope, ifail, xmin, xmax);
01552 
01553    TVirtualFitter *grFitter = TVirtualFitter::GetFitter();
01554    TF1 *f1 = (TF1*)grFitter->GetUserFunc();
01555    f1->SetParameter(0,constant);
01556    f1->SetParameter(1,slope);
01557 }
01558 
01559 
01560 //______________________________________________________________________________
01561 void TGraph::InitPolynom(Double_t xmin, Double_t xmax)
01562 {
01563    // Compute Initial values of parameters for a polynom.
01564 
01565    Double_t fitpar[25];
01566 
01567    TVirtualFitter *grFitter = TVirtualFitter::GetFitter();
01568    TF1 *f1 = (TF1*)grFitter->GetUserFunc();
01569    Int_t npar   = f1->GetNpar();
01570    if (xmax <= xmin) {xmin = fX[0]; xmax = fX[fNpoints-1];}
01571 
01572    LeastSquareFit(npar, fitpar, xmin, xmax);
01573 
01574    for (Int_t i=0;i<npar;i++) f1->SetParameter(i, fitpar[i]);
01575 }
01576 
01577 
01578 //______________________________________________________________________________
01579 Int_t TGraph::InsertPoint()
01580 {
01581    // Insert a new point at the mouse position
01582 
01583    Int_t px = gPad->GetEventX();
01584    Int_t py = gPad->GetEventY();
01585 
01586    //localize point where to insert
01587    Int_t ipoint = -2;
01588    Int_t i,d=0;
01589    // start with a small window (in case the mouse is very close to one point)
01590    for (i=0;i<fNpoints-1;i++) {
01591       d = DistancetoLine(px, py, gPad->XtoPad(fX[i]), gPad->YtoPad(fY[i]), gPad->XtoPad(fX[i+1]), gPad->YtoPad(fY[i+1]));
01592       if (d < 5) {ipoint = i+1; break;}
01593    }
01594    if (ipoint == -2) {
01595       //may be we are far from one point, try again with a larger window
01596       for (i=0;i<fNpoints-1;i++) {
01597          d = DistancetoLine(px, py, gPad->XtoPad(fX[i]), gPad->YtoPad(fY[i]), gPad->XtoPad(fX[i+1]), gPad->YtoPad(fY[i+1]));
01598          if (d < 10) {ipoint = i+1; break;}
01599       }
01600    }
01601    if (ipoint == -2) {
01602       //distinguish between first and last point
01603       Int_t dpx = px - gPad->XtoAbsPixel(gPad->XtoPad(fX[0]));
01604       Int_t dpy = py - gPad->YtoAbsPixel(gPad->XtoPad(fY[0]));
01605       if (dpx*dpx+dpy*dpy < 25) ipoint = 0;
01606       else                      ipoint = fNpoints;
01607    }
01608    Double_t **ps = ExpandAndCopy(fNpoints + 1, ipoint);
01609    CopyAndRelease(ps, ipoint, fNpoints++, ipoint + 1);
01610 
01611    // To avoid redefenitions in descendant classes
01612    FillZero(ipoint, ipoint + 1);
01613 
01614    fX[ipoint] = gPad->PadtoX(gPad->AbsPixeltoX(px));
01615    fY[ipoint] = gPad->PadtoY(gPad->AbsPixeltoY(py));
01616    gPad->Modified();
01617    return ipoint;
01618 }
01619 
01620 //______________________________________________________________________________
01621 Double_t TGraph::Integral(Int_t first, Int_t last) const
01622 {
01623    // Integrate the TGraph data within a given (index) range
01624    // NB: if last=-1 (default) last is set to the last point.
01625    //     if (first <0) the first point (0) is taken.
01626    //   : The graph segments should not intersect.
01627    //Method:
01628    // There are many ways to calculate the surface of a polygon. It all depends on what kind of data
01629    // you have to deal with. The most evident solution would be to divide the polygon in triangles and
01630    // calculate the surface of them. But this can quickly become complicated as you will have to test
01631    // every segments of every triangles and check if they are intersecting with a current polygon's
01632    // segment or if it goes outside the polygon. Many calculations that would lead to many problems...
01633    //      The solution (implemented by R.Brun)
01634    // Fortunately for us, there is a simple way to solve this problem, as long as the polygon's
01635    // segments don't intersect.
01636    // It takes the x coordinate of the current vertex and multiply it by the y coordinate of the next
01637    // vertex. Then it subtracts from it the result of the y coordinate of the current vertex multiplied
01638    // by the x coordinate of the next vertex. Then divide the result by 2 to get the surface/area.
01639    //      Sources
01640    //      http://forums.wolfram.com/mathgroup/archive/1998/Mar/msg00462.html
01641    //      http://stackoverflow.com/questions/451426/how-do-i-calculate-the-surface-area-of-a-2d-polygon
01642 
01643    if (first < 0) first = 0;
01644    if (last < 0) last = fNpoints-1;
01645    if(last >= fNpoints) last = fNpoints-1;
01646    if (first >= last) return 0;
01647    Int_t np = last-first+1;
01648    Double_t sum = 0.0;
01649    //for(Int_t i=first;i<=last;i++) {
01650    //   Int_t j = first + (i-first+1)%np;
01651    //   sum += TMath::Abs(fX[i]*fY[j]);
01652    //   sum -= TMath::Abs(fY[i]*fX[j]);
01653    //}
01654    for(Int_t i=first;i<=last;i++) {
01655       Int_t j = first + (i-first+1)%np;
01656       sum += (fY[i]+fY[j])*(fX[j]-fX[i]);
01657    }
01658    return 0.5*TMath::Abs(sum);
01659 }
01660 
01661 
01662 //______________________________________________________________________________
01663 Int_t TGraph::IsInside(Double_t x, Double_t y) const
01664 {
01665    // Return 1 if the point (x,y) is inside the polygon defined by
01666    // the graph vertices 0 otherwise.
01667    //
01668    // Algorithm:
01669    // The loop is executed with the end-point coordinates of a line segment
01670    // (X1,Y1)-(X2,Y2) and the Y-coordinate of a horizontal line.
01671    // The counter inter is incremented if the line (X1,Y1)-(X2,Y2) intersects
01672    // the horizontal line. In this case XINT is set to the X-coordinate of the
01673    // intersection point. If inter is an odd number, then the point x,y is within
01674    // the polygon.
01675 
01676    return (Int_t)TMath::IsInside(x, y, fNpoints, fX, fY);
01677 }
01678 
01679 
01680 //______________________________________________________________________________
01681 void TGraph::LeastSquareFit(Int_t m, Double_t *a, Double_t xmin, Double_t xmax)
01682 {
01683    // Least squares polynomial fitting without weights.
01684    //
01685    //  m     number of parameters
01686    //  a     array of parameters
01687    //  first 1st point number to fit (default =0)
01688    //  last  last point number to fit (default=fNpoints-1)
01689    //
01690    //   based on CERNLIB routine LSQ: Translated to C++ by Rene Brun
01691 
01692    const Double_t zero = 0.;
01693    const Double_t one = 1.;
01694    const Int_t idim = 20;
01695 
01696    Double_t  b[400]        /* was [20][20] */;
01697    Int_t i, k, l, ifail;
01698    Double_t power;
01699    Double_t da[20], xk, yk;
01700    Int_t n = fNpoints;
01701    if (xmax <= xmin) {xmin = fX[0]; xmax = fX[fNpoints-1];}
01702 
01703    if (m <= 2) {
01704       LeastSquareLinearFit(n, a[0], a[1], ifail, xmin, xmax);
01705       return;
01706    }
01707    if (m > idim || m > n) return;
01708    da[0] = zero;
01709    for (l = 2; l <= m; ++l) {
01710       b[l-1]           = zero;
01711       b[m + l*20 - 21] = zero;
01712       da[l-1]          = zero;
01713    }
01714    Int_t np = 0;
01715    for (k = 0; k < fNpoints; ++k) {
01716       xk     = fX[k];
01717       if (xk < xmin || xk > xmax) continue;
01718       np++;
01719       yk     = fY[k];
01720       power  = one;
01721       da[0] += yk;
01722       for (l = 2; l <= m; ++l) {
01723          power   *= xk;
01724          b[l-1]  += power;
01725          da[l-1] += power*yk;
01726       }
01727       for (l = 2; l <= m; ++l) {
01728          power            *= xk;
01729          b[m + l*20 - 21] += power;
01730       }
01731    }
01732    b[0]  = Double_t(np);
01733    for (i = 3; i <= m; ++i) {
01734       for (k = i; k <= m; ++k) {
01735          b[k - 1 + (i-1)*20 - 21] = b[k + (i-2)*20 - 21];
01736       }
01737    }
01738    H1LeastSquareSeqnd(m, b, idim, ifail, 1, da);
01739 
01740    if (ifail < 0) {
01741       a[0] = fY[0];
01742       for (i=1; i<m; ++i) a[i] = 0;
01743       return;
01744    }
01745    for (i=0; i<m; ++i) a[i] = da[i];
01746 }
01747 
01748 
01749 //______________________________________________________________________________
01750 void TGraph::LeastSquareLinearFit(Int_t ndata, Double_t &a0, Double_t &a1, Int_t &ifail, Double_t xmin, Double_t xmax)
01751 {
01752    // Least square linear fit without weights.
01753    //
01754    //  Fit a straight line (a0 + a1*x) to the data in this graph.
01755    //  ndata:  if ndata<0, fits the logarithm of the graph (used in InitExpo() to set
01756    //          the initial parameter values for a fit with exponential function.
01757    //  a0:     constant
01758    //  a1:     slope
01759    //  ifail:  return parameter indicating the status of the fit (ifail=0, fit is OK)
01760    //  xmin, xmax: fitting range
01761    //
01762    //  extracted from CERNLIB LLSQ: Translated to C++ by Rene Brun
01763 
01764    Double_t xbar, ybar, x2bar;
01765    Int_t i;
01766    Double_t xybar;
01767    Double_t fn, xk, yk;
01768    Double_t det;
01769    if (xmax <= xmin) {xmin = fX[0]; xmax = fX[fNpoints-1];}
01770 
01771    ifail = -2;
01772    xbar  = ybar = x2bar = xybar = 0;
01773    Int_t np = 0;
01774    for (i = 0; i < fNpoints; ++i) {
01775       xk = fX[i];
01776       if (xk < xmin || xk > xmax) continue;
01777       np++;
01778       yk = fY[i];
01779       if (ndata < 0) {
01780          if (yk <= 0) yk = 1e-9;
01781          yk = TMath::Log(yk);
01782       }
01783       xbar  += xk;
01784       ybar  += yk;
01785       x2bar += xk*xk;
01786       xybar += xk*yk;
01787    }
01788    fn    = Double_t(np);
01789    det   = fn*x2bar - xbar*xbar;
01790    ifail = -1;
01791    if (det <= 0) {
01792       if (fn > 0) a0 = ybar/fn;
01793       else        a0 = 0;
01794       a1 = 0;
01795       return;
01796    }
01797    ifail = 0;
01798    a0 = (x2bar*ybar - xbar*xybar) / det;
01799    a1 = (fn*xybar - xbar*ybar) / det;
01800 }
01801 
01802 
01803 //______________________________________________________________________________
01804 void TGraph::Paint(Option_t *option)
01805 {
01806    // Draw this graph with its current attributes.
01807 
01808    TVirtualGraphPainter::GetPainter()->PaintHelper(this, option);
01809 }
01810 
01811 
01812 //______________________________________________________________________________
01813 void TGraph::PaintGraph(Int_t npoints, const Double_t *x, const Double_t *y, Option_t *chopt)
01814 {
01815    // Draw the (x,y) as a graph.
01816 
01817    TVirtualGraphPainter::GetPainter()->PaintGraph(this, npoints, x, y, chopt);
01818 }
01819 
01820 
01821 //______________________________________________________________________________
01822 void TGraph::PaintGrapHist(Int_t npoints, const Double_t *x, const Double_t *y, Option_t *chopt)
01823 {
01824    // Draw the (x,y) as a histogram.
01825 
01826    TVirtualGraphPainter::GetPainter()->PaintGrapHist(this, npoints, x, y, chopt);
01827 }
01828 
01829 
01830 //______________________________________________________________________________
01831 void TGraph::PaintStats(TF1 *fit)
01832 {
01833    // Draw the stats
01834 
01835    TVirtualGraphPainter::GetPainter()->PaintStats(this, fit);
01836 }
01837 
01838 
01839 //______________________________________________________________________________
01840 void TGraph::Print(Option_t *) const
01841 {
01842    // Print graph values.
01843 
01844    for (Int_t i=0;i<fNpoints;i++) {
01845       printf("x[%d]=%g, y[%d]=%g\n",i,fX[i],i,fY[i]);
01846    }
01847 }
01848 
01849 
01850 //______________________________________________________________________________
01851 void TGraph::RecursiveRemove(TObject *obj)
01852 {
01853    // Recursively remove object from the list of functions
01854 
01855    if (fFunctions) {
01856       if (!fFunctions->TestBit(kInvalidObject)) fFunctions->RecursiveRemove(obj);
01857    }
01858    if (fHistogram == obj) fHistogram = 0;
01859 }
01860 
01861 
01862 //______________________________________________________________________________
01863 Int_t TGraph::RemovePoint()
01864 {
01865    // Delete point close to the mouse position
01866 
01867    Int_t px = gPad->GetEventX();
01868    Int_t py = gPad->GetEventY();
01869 
01870    //localize point to be deleted
01871    Int_t ipoint = -2;
01872    Int_t i;
01873    // start with a small window (in case the mouse is very close to one point)
01874    for (i=0;i<fNpoints;i++) {
01875       Int_t dpx = px - gPad->XtoAbsPixel(gPad->XtoPad(fX[i]));
01876       Int_t dpy = py - gPad->YtoAbsPixel(gPad->YtoPad(fY[i]));
01877       if (dpx*dpx+dpy*dpy < 100) {ipoint = i; break;}
01878    }
01879    return RemovePoint(ipoint);
01880 }
01881 
01882 
01883 //______________________________________________________________________________
01884 Int_t TGraph::RemovePoint(Int_t ipoint)
01885 {
01886    // Delete point number ipoint
01887 
01888    if (ipoint < 0) return -1;
01889    if (ipoint >= fNpoints) return -1;
01890 
01891    Double_t **ps = ShrinkAndCopy(fNpoints - 1, ipoint);
01892    CopyAndRelease(ps, ipoint+1, fNpoints--, ipoint);
01893    if (gPad) gPad->Modified();
01894    return ipoint;
01895 }
01896 
01897 
01898 //______________________________________________________________________________
01899 void TGraph::SavePrimitive(ostream &out, Option_t *option /*= ""*/)
01900 {
01901     // Save primitive as a C++ statement(s) on output stream out
01902 
01903    char quote = '"';
01904    out<<"   "<<endl;
01905    if (gROOT->ClassSaved(TGraph::Class())) {
01906       out<<"   ";
01907    } else {
01908       out<<"   TGraph *";
01909    }
01910    out<<"graph = new TGraph("<<fNpoints<<");"<<endl;
01911    out<<"   graph->SetName("<<quote<<GetName()<<quote<<");"<<endl;
01912    out<<"   graph->SetTitle("<<quote<<GetTitle()<<quote<<");"<<endl;
01913 
01914    SaveFillAttributes(out,"graph",0,1001);
01915    SaveLineAttributes(out,"graph",1,1,1);
01916    SaveMarkerAttributes(out,"graph",1,1,1);
01917 
01918    for (Int_t i=0;i<fNpoints;i++) {
01919       out<<"   graph->SetPoint("<<i<<","<<fX[i]<<","<<fY[i]<<");"<<endl;
01920    }
01921 
01922    static Int_t frameNumber = 0;
01923    if (fHistogram) {
01924       frameNumber++;
01925       TString hname = fHistogram->GetName();
01926       hname += frameNumber;
01927       fHistogram->SetName(Form("Graph_%s",hname.Data()));
01928       fHistogram->SavePrimitive(out,"nodraw");
01929       out<<"   graph->SetHistogram("<<fHistogram->GetName()<<");"<<endl;
01930       out<<"   "<<endl;
01931    }
01932 
01933    // save list of functions
01934    TIter next(fFunctions);
01935    TObject *obj;
01936    while ((obj=next())) {
01937       obj->SavePrimitive(out,"nodraw");
01938       if (obj->InheritsFrom("TPaveStats")) {
01939          out<<"   graph->GetListOfFunctions()->Add(ptstats);"<<endl;
01940          out<<"   ptstats->SetParent(graph->GetListOfFunctions());"<<endl;
01941       } else {
01942          out<<"   graph->GetListOfFunctions()->Add("<<obj->GetName()<<");"<<endl;
01943       }
01944    }
01945 
01946    const char *l;
01947    l = strstr(option,"multigraph");
01948    if (l) {
01949       out<<"   multigraph->Add(graph,"<<quote<<l+10<<quote<<");"<<endl;
01950       return;
01951    } 
01952    l = strstr(option,"th2poly");
01953    if (l) {
01954       out<<"   "<<l+7<<"->AddBin(graph);"<<endl;
01955       return;
01956    } 
01957    out<<"   graph->Draw("<<quote<<option<<quote<<");"<<endl;
01958 }
01959 
01960 
01961 //______________________________________________________________________________
01962 void TGraph::Set(Int_t n)
01963 {
01964    // Set number of points in the graph
01965    // Existing coordinates are preserved
01966    // New coordinates above fNpoints are preset to 0.
01967 
01968    if (n < 0) n = 0;
01969    if (n == fNpoints) return;
01970    Double_t **ps = Allocate(n);
01971    CopyAndRelease(ps, 0, TMath::Min(fNpoints,n), 0);
01972    if (n > fNpoints) {
01973       FillZero(fNpoints, n, kFALSE);
01974    }
01975    fNpoints = n;
01976 }
01977 
01978 
01979 //______________________________________________________________________________
01980 Bool_t TGraph::GetEditable() const
01981 {
01982    // Return kTRUE if kNotEditable bit is not set, kFALSE otherwise.
01983 
01984    return TestBit(kNotEditable) ? kFALSE : kTRUE;
01985 }
01986 
01987 
01988 //______________________________________________________________________________
01989 void TGraph::SetEditable(Bool_t editable)
01990 {
01991    // if editable=kFALSE, the graph cannot be modified with the mouse
01992    //  by default a TGraph is editable
01993 
01994    if (editable) ResetBit(kNotEditable);
01995    else          SetBit(kNotEditable);
01996 }
01997 
01998 
01999 //______________________________________________________________________________
02000 void TGraph::SetMaximum(Double_t maximum)
02001 {
02002    // Set the maximum of the graph.
02003 
02004    fMaximum = maximum;
02005    GetHistogram()->SetMaximum(maximum);
02006 }
02007 
02008 
02009 //______________________________________________________________________________
02010 void TGraph::SetMinimum(Double_t minimum)
02011 {
02012    // Set the minimum of the graph.
02013 
02014    fMinimum = minimum;
02015    GetHistogram()->SetMinimum(minimum);
02016 }
02017 
02018 
02019 //______________________________________________________________________________
02020 void TGraph::SetPoint(Int_t i, Double_t x, Double_t y)
02021 {
02022    // Set x and y values for point number i.
02023 
02024    if (i < 0) return;
02025    if (fHistogram) {
02026       delete fHistogram;
02027       fHistogram = 0;
02028    }
02029    if (i >= fMaxSize) {
02030       Double_t **ps = ExpandAndCopy(i+1, fNpoints);
02031       CopyAndRelease(ps, 0,0,0);
02032    }
02033    if (i >= fNpoints) {
02034       // points above i can be not initialized
02035       // set zero up to i-th point to avoid redefenition
02036       // of this method in descendant classes
02037       FillZero(fNpoints, i + 1);
02038       fNpoints = i+1;
02039    }
02040    fX[i] = x;
02041    fY[i] = y;
02042    if (gPad) gPad->Modified();
02043 }
02044 
02045 
02046 //______________________________________________________________________________
02047 void TGraph::SetTitle(const char* title)
02048 {
02049    // Set graph title.
02050 
02051    fTitle = title;
02052    if (fHistogram) fHistogram->SetTitle(title);
02053 }
02054 
02055 
02056 //______________________________________________________________________________
02057 Double_t **TGraph::ShrinkAndCopy(Int_t size, Int_t oend)
02058 {
02059    // if size*2 <= fMaxSize allocate new arrays of size points,
02060    // copy points [0,oend).
02061    // Return newarray (passed or new instance if it was zero
02062    // and allocations are needed)
02063    if (size*2 > fMaxSize || !fMaxSize) {
02064       return 0;
02065    }
02066    Double_t **newarrays = Allocate(size);
02067    CopyPoints(newarrays, 0, oend, 0);
02068    return newarrays;
02069 }
02070 
02071 
02072 //______________________________________________________________________________
02073 void TGraph::Sort(Bool_t (*greaterfunc)(const TGraph*, Int_t, Int_t) /*=TGraph::CompareX()*/,
02074                   Bool_t ascending /*=kTRUE*/, Int_t low /* =0 */, Int_t high /* =-1111 */)
02075 {
02076    // Sorts the points of this TGraph using in-place quicksort (see e.g. older glibc).
02077    // To compare two points the function parameter greaterfunc is used (see TGraph::CompareX for an
02078    // example of such a method, which is also the default comparison function for Sort). After
02079    // the sort, greaterfunc(this, i, j) will return kTRUE for all i>j if ascending == kTRUE, and
02080    // kFALSE otherwise.
02081    //
02082    // The last two parameters are used for the recursive quick sort, stating the range to be sorted
02083    //
02084    // Examples:
02085    //   // sort points along x axis
02086    //   graph->Sort();
02087    //   // sort points along their distance to origin
02088    //   graph->Sort(&TGraph::CompareRadius);
02089    //
02090    //   Bool_t CompareErrors(const TGraph* gr, Int_t i, Int_t j) {
02091    //     const TGraphErrors* ge=(const TGraphErrors*)gr;
02092    //     return (ge->GetEY()[i]>ge->GetEY()[j]); }
02093    //   // sort using the above comparison function, largest errors first
02094    //   graph->Sort(&CompareErrors, kFALSE);
02095 
02096    if (high == -1111) high = GetN()-1;
02097    //  Termination condition
02098    if (high <= low) return;
02099 
02100    int left, right;
02101    left = low; // low is the pivot element
02102    right = high;
02103    while (left < right) {
02104       // move left while item < pivot
02105       while(left <= high && greaterfunc(this, left, low) != ascending)
02106          left++;
02107       // move right while item > pivot
02108       while(right > low && greaterfunc(this, right, low) == ascending)
02109          right--;
02110       if (left < right && left < high && right > low)
02111          SwapPoints(left, right);
02112    }
02113    // right is final position for the pivot
02114    if (right > low)
02115       SwapPoints(low, right);
02116    Sort( greaterfunc, ascending, low, right-1 );
02117    Sort( greaterfunc, ascending, right+1, high );
02118 }
02119 
02120 
02121 //______________________________________________________________________________
02122 void TGraph::Streamer(TBuffer &b)
02123 {
02124    // Stream an object of class TGraph.
02125 
02126    if (b.IsReading()) {
02127       UInt_t R__s, R__c;
02128       Version_t R__v = b.ReadVersion(&R__s, &R__c);
02129       if (R__v > 2) {
02130          b.ReadClassBuffer(TGraph::Class(), this, R__v, R__s, R__c);
02131          if (fHistogram) fHistogram->SetDirectory(0);
02132          TIter next(fFunctions);
02133          TObject *obj;
02134          while ((obj = next())) {
02135             if (obj->InheritsFrom(TF1::Class())) {
02136                TF1 *f1 = (TF1*)obj;
02137                f1->SetParent(this);
02138             }
02139          }
02140          fMaxSize = fNpoints;
02141          return;
02142       }
02143       //====process old versions before automatic schema evolution
02144       TNamed::Streamer(b);
02145       TAttLine::Streamer(b);
02146       TAttFill::Streamer(b);
02147       TAttMarker::Streamer(b);
02148       b >> fNpoints;
02149       fMaxSize = fNpoints;
02150       fX = new Double_t[fNpoints];
02151       fY = new Double_t[fNpoints];
02152       if (R__v < 2) {
02153          Float_t *x = new Float_t[fNpoints];
02154          Float_t *y = new Float_t[fNpoints];
02155          b.ReadFastArray(x,fNpoints);
02156          b.ReadFastArray(y,fNpoints);
02157          for (Int_t i=0;i<fNpoints;i++) {
02158             fX[i] = x[i];
02159             fY[i] = y[i];
02160          }
02161          delete [] y;
02162          delete [] x;
02163       } else {
02164          b.ReadFastArray(fX,fNpoints);
02165          b.ReadFastArray(fY,fNpoints);
02166       }
02167       b >> fFunctions;
02168       b >> fHistogram;
02169       if (fHistogram) fHistogram->SetDirectory(0);
02170       if (R__v < 2) {
02171          Float_t mi,ma;
02172          b >> mi;
02173          b >> ma;
02174          fMinimum = mi;
02175          fMaximum = ma;
02176       } else {
02177          b >> fMinimum;
02178          b >> fMaximum;
02179       }
02180       b.CheckByteCount(R__s, R__c, TGraph::IsA());
02181       //====end of old versions
02182 
02183    } else {
02184       b.WriteClassBuffer(TGraph::Class(),this);
02185    }
02186 }
02187 
02188 
02189 //______________________________________________________________________________
02190 void TGraph::SwapPoints(Int_t pos1, Int_t pos2)
02191 {
02192    // Swap points.
02193 
02194    SwapValues(fX, pos1, pos2);
02195    SwapValues(fY, pos1, pos2);
02196 }
02197 
02198 
02199 //______________________________________________________________________________
02200 void TGraph::SwapValues(Double_t* arr, Int_t pos1, Int_t pos2)
02201 {
02202    // Swap values.
02203 
02204    Double_t tmp=arr[pos1];
02205    arr[pos1]=arr[pos2];
02206    arr[pos2]=tmp;
02207 }
02208 
02209 
02210 //______________________________________________________________________________
02211 void TGraph::UseCurrentStyle()
02212 {
02213    // Set current style settings in this graph
02214    // This function is called when either TCanvas::UseCurrentStyle
02215    // or TROOT::ForceStyle have been invoked.
02216 
02217    if (gStyle->IsReading()) {
02218       SetFillColor(gStyle->GetHistFillColor());
02219       SetFillStyle(gStyle->GetHistFillStyle());
02220       SetLineColor(gStyle->GetHistLineColor());
02221       SetLineStyle(gStyle->GetHistLineStyle());
02222       SetLineWidth(gStyle->GetHistLineWidth());
02223       SetMarkerColor(gStyle->GetMarkerColor());
02224       SetMarkerStyle(gStyle->GetMarkerStyle());
02225       SetMarkerSize(gStyle->GetMarkerSize());
02226    } else {
02227       gStyle->SetHistFillColor(GetFillColor());
02228       gStyle->SetHistFillStyle(GetFillStyle());
02229       gStyle->SetHistLineColor(GetLineColor());
02230       gStyle->SetHistLineStyle(GetLineStyle());
02231       gStyle->SetHistLineWidth(GetLineWidth());
02232       gStyle->SetMarkerColor(GetMarkerColor());
02233       gStyle->SetMarkerStyle(GetMarkerStyle());
02234       gStyle->SetMarkerSize(GetMarkerSize());
02235    }
02236    if (fHistogram) fHistogram->UseCurrentStyle();
02237 
02238    TIter next(GetListOfFunctions());
02239    TObject *obj;
02240 
02241    while ((obj = next())) {
02242       obj->UseCurrentStyle();
02243    }
02244 }
02245 
02246 
02247 //______________________________________________________________________________
02248 Int_t TGraph::Merge(TCollection* li)
02249 {
02250    // Adds all graphs from the collection to this graph.
02251    // Returns the total number of poins in the result or -1 in case of an error.
02252 
02253    TIter next(li);
02254    while (TObject* o = next()) {
02255       TGraph *g = dynamic_cast<TGraph*> (o);
02256       if (!g) {
02257          Error("Merge",
02258              "Cannot merge - an object which doesn't inherit from TGraph found in the list");
02259          return -1;
02260       }
02261       Double_t x, y;
02262       for (Int_t i = 0 ; i < g->GetN(); i++) {
02263          g->GetPoint(i, x, y);
02264          SetPoint(GetN(), x, y);
02265       }
02266    }
02267    return GetN();
02268 }
02269 
02270 
02271 //______________________________________________________________________________
02272 void TGraph::Zero(Int_t &k,Double_t AZ,Double_t BZ,Double_t E2,Double_t &X,Double_t &Y
02273                  ,Int_t maxiterations)
02274 {
02275    // Find zero of a continuous function.
02276    // This function finds a real zero of the continuous real
02277    // function Y(X) in a given interval (A,B). See accompanying
02278    // notes for details of the argument list and calling sequence
02279 
02280    static Double_t a, b, ya, ytest, y1, x1, h;
02281    static Int_t j1, it, j3, j2;
02282    Double_t yb, x2;
02283    yb = 0;
02284 
02285    //       Calculate Y(X) at X=AZ.
02286    if (k <= 0) {
02287       a  = AZ;
02288       b  = BZ;
02289       X  = a;
02290       j1 = 1;
02291       it = 1;
02292       k  = j1;
02293       return;
02294    }
02295 
02296    //       Test whether Y(X) is sufficiently small.
02297 
02298    if (TMath::Abs(Y) <= E2) { k = 2; return; }
02299 
02300    //       Calculate Y(X) at X=BZ.
02301 
02302    if (j1 == 1) {
02303       ya = Y;
02304       X  = b;
02305       j1 = 2;
02306       return;
02307    }
02308    //       Test whether the signs of Y(AZ) and Y(BZ) are different.
02309    //       if not, begin the binary subdivision.
02310 
02311    if (j1 != 2) goto L100;
02312    if (ya*Y < 0) goto L120;
02313    x1 = a;
02314    y1 = ya;
02315    j1 = 3;
02316    h  = b - a;
02317    j2 = 1;
02318    x2 = a + 0.5*h;
02319    j3 = 1;
02320    it++;      //*-*-   Check whether (maxiterations) function values have been calculated.
02321    if (it >= maxiterations) k = j1;
02322    else                     X = x2;
02323    return;
02324 
02325    //      Test whether a bracket has been found .
02326    //      If not,continue the search
02327 
02328 L100:
02329    if (j1 > 3) goto L170;
02330    if (ya*Y >= 0) {
02331       if (j3 >= j2) {
02332          h  = 0.5*h; j2 = 2*j2;
02333          a  = x1;  ya = y1;  x2 = a + 0.5*h; j3 = 1;
02334       }
02335       else {
02336          a  = X;   ya = Y;   x2 = X + h;     j3++;
02337       }
02338       it++;
02339       if (it >= maxiterations) k = j1;
02340       else                     X = x2;
02341       return;
02342    }
02343 
02344    //       The first bracket has been found.calculate the next X by the
02345    //       secant method based on the bracket.
02346 
02347 L120:
02348    b  = X;
02349    yb = Y;
02350    j1 = 4;
02351 L130:
02352    if (TMath::Abs(ya) > TMath::Abs(yb)) { x1 = a; y1 = ya; X  = b; Y  = yb; }
02353    else                                 { x1 = b; y1 = yb; X  = a; Y  = ya; }
02354 
02355    //       Use the secant method based on the function values y1 and Y.
02356    //       check that x2 is inside the interval (a,b).
02357 
02358 L150:
02359    x2    = X-Y*(X-x1)/(Y-y1);
02360    x1    = X;
02361    y1    = Y;
02362    ytest = 0.5*TMath::Min(TMath::Abs(ya),TMath::Abs(yb));
02363    if ((x2-a)*(x2-b) < 0) {
02364       it++;
02365       if (it >= maxiterations) k = j1;
02366       else                     X = x2;
02367       return;
02368    }
02369 
02370    //       Calculate the next value of X by bisection . Check whether
02371    //       the maximum accuracy has been achieved.
02372 
02373 L160:
02374    x2    = 0.5*(a+b);
02375    ytest = 0;
02376    if ((x2-a)*(x2-b) >= 0) { k = 2;  return; }
02377    it++;
02378    if (it >= maxiterations) k = j1;
02379    else                     X = x2;
02380    return;
02381 
02382 
02383    //       Revise the bracket (a,b).
02384 
02385 L170:
02386    if (j1 != 4) return;
02387    if (ya*Y < 0) { b  = X; yb = Y; }
02388    else          { a  = X; ya = Y; }
02389 
02390    //       Use ytest to decide the method for the next value of X.
02391 
02392    if (ytest <= 0) goto L130;
02393    if (TMath::Abs(Y)-ytest <= 0) goto L150;
02394    goto L160;
02395 }

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