library: libGraf
#include "TMultiGraph.h"

TMultiGraph


class description - source file - inheritance tree (.pdf)

class TMultiGraph : public TNamed

Inheritance Chart:
TObject
<-
TNamed
<-
TMultiGraph

    public:
TMultiGraph() TMultiGraph(const char* name, const char* title) TMultiGraph(const TMultiGraph&) virtual ~TMultiGraph() virtual void Add(TGraph* graph, Option_t* chopt = "") virtual void Browse(TBrowser* b) static TClass* Class() virtual Int_t DistancetoPrimitive(Int_t px, Int_t py) virtual void Draw(Option_t* chopt = "") virtual Int_t Fit(const char* formula, Option_t* option = "", Option_t* goption = "", Axis_t xmin = 0, Axis_t xmax = 0) virtual Int_t Fit(TF1* f1, Option_t* option = "", Option_t* goption = "", Axis_t rxmin = 0, Axis_t rxmax = 0) TF1* GetFunction(const char* name) const virtual Option_t* GetGraphDrawOption(const TGraph* gr) const TH1F* GetHistogram() const TList* GetListOfFunctions() const TList* GetListOfGraphs() const TAxis* GetXaxis() const TAxis* GetYaxis() const virtual void InitExpo(Double_t xmin, Double_t xmax) virtual void InitGaus(Double_t xmin, Double_t xmax) virtual void InitPolynom(Double_t xmin, Double_t xmax) virtual TClass* IsA() const virtual void LeastSquareFit(Int_t m, Double_t* a, Double_t xmin, Double_t xmax) virtual void LeastSquareLinearFit(Int_t ndata, Double_t& a0, Double_t& a1, Int_t& ifail, Double_t xmin, Double_t xmax) TMultiGraph& operator=(const TMultiGraph&) virtual void Paint(Option_t* chopt = "") virtual void Print(Option_t* chopt = "") const virtual void RecursiveRemove(TObject* obj) virtual void SavePrimitive(ofstream& out, Option_t* option) virtual void SetMaximum(Double_t maximum = -1111) virtual void SetMinimum(Double_t minimum = -1111) virtual void ShowMembers(TMemberInspector& insp, char* parent) virtual void Streamer(TBuffer& b) void StreamerNVirtual(TBuffer& b)

Data Members


    protected:
TList* fGraphs Pointer to list of TGraphs TList* fFunctions Pointer to list of functions (fits and user) TH1F* fHistogram Pointer to histogram used for drawing axis Double_t fMaximum Maximum value for plotting along y Double_t fMinimum Minimum value for plotting along y

Class Description

   A TMultiGraph is a collection of TGraph (or derived) objects
   Use TMultiGraph::Add to add a new graph to the list.
   The TMultiGraph owns the objects in the list.
   Drawing options are the same as for TGraph
   Example;
     TGraph *gr1 = new TGraph(...
     TGraphErrors *gr2 = new TGraphErrors(...
     TMultiGraph *mg = new TMultiGraph();
     mg->Add(gr1,"lp");
     mg->Add(gr2,"cp");
     mg->Draw("a");

  The drawing option for each TGraph may be specified as an optional
  second argument of the Add function.
  If a draw option is specified, it will be used to draw the graph,
  otherwise the graph will be drawn with the option specified in
  TMultiGraph::Draw

TMultiGraph(): TNamed()
 TMultiGraph default constructor

TMultiGraph(const char *name, const char *title) : TNamed(name,title)
 constructor with name and title

~TMultiGraph()
 TMultiGraph destructor

void Add(TGraph *graph, Option_t *chopt)
 add a new graph to the list of graphs
 note that the graph is now owned by the TMultigraph.
 Deleting the TMultiGraph object will automatically delete the graphs.
 You should not delete the graphs when the TMultigraph is still active.

void Browse(TBrowser *)

Int_t DistancetoPrimitive(Int_t px, Int_t py)
 Compute distance from point px,py to each graph


void Draw(Option_t *option)
*-*-*-*-*-*-*-*-*-*-*Draw this multigraph with its current attributes*-*-*-*-*-*-*
*-*                  ==========================================

   Options to draw a graph are described in TGraph::PainGraph

  The drawing option for each TGraph may be specified as an optional
  second argument of the Add function. You can use GetGraphDrawOption
  to return this option.
  If a draw option is specified, it will be used to draw the graph,
  otherwise the graph will be drawn with the option specified in
  TMultiGraph::Draw. Use GetDrawOption to return the option specified
  when drawin the TMultiGraph.

Int_t Fit(const char *fname, Option_t *option, Option_t *, Axis_t xmin, Axis_t xmax)
*-*-*-*-*-*Fit this graph with function with name fname*-*-*-*-*-*-*-*-*-*
*-*        ============================================
  interface to TF1::Fit(TF1 *f1...

Int_t Fit(TF1 *f1, Option_t *option, Option_t *, Axis_t rxmin, Axis_t rxmax)
*-*-*-*-*-*-*-*-*-*-*Fit this multigraph with function f1*-*-*-*-*-*-*-*-*-*
*-*                  ==================================

   In this function all graphs of the multigraph are fitted simultaneously

   f1 is an already predefined function created by TF1.
   Predefined functions such as gaus, expo and poln are automatically
   created by ROOT.

   The list of fit options is given in parameter option.
      option = "W"  Set all errors to 1
             = "U" Use a User specified fitting algorithm (via SetFCN)
             = "Q" Quiet mode (minimum printing)
             = "V" Verbose mode (default is between Q and V)
             = "B" Use this option when you want to fix one or more parameters
                   and the fitting function is like "gaus","expo","poln","landau".
             = "R" Use the Range specified in the function range
             = "N" Do not store the graphics function, do not draw
             = "0" Do not plot the result of the fit. By default the fitted function
                   is drawn unless the option"N" above is specified.
             = "+" Add this new fitted function to the list of fitted functions
                   (by default, any previous function is deleted)
             = "C" In case of linear fitting, not calculate the chisquare
                    (saves time)
             = "F" If fitting a polN, switch to minuit fitter
             = "ROB" In case of linear fitting, compute the LTS regression
                     coefficients (robust(resistant) regression), using
                     the default fraction of good points
               "ROB=0.x" - compute the LTS regression coefficients, using
                           0.x as a fraction of good points

   When the fit is drawn (by default), the parameter goption may be used
   to specify a list of graphics options. See TGraph::Paint for a complete
   list of these options.

   In order to use the Range option, one must first create a function
   with the expression to be fitted. For example, if your graph
   has a defined range between -4 and 4 and you want to fit a gaussian
   only in the interval 1 to 3, you can do:
        TF1 *f1 = new TF1("f1","gaus",1,3);
        graph->Fit("f1","R");


   who is calling this function
   ============================
   Note that this function is called when calling TGraphErrors::Fit
   or TGraphAsymmErrors::Fit ot TGraphBentErrors::Fit
   see the discussion below on the errors calulation.

   Setting initial conditions
   ==========================
   Parameters must be initialized before invoking the Fit function.
   The setting of the parameter initial values is automatic for the
   predefined functions : poln, expo, gaus, landau. One can however disable
   this automatic computation by specifying the option "B".
   You can specify boundary limits for some or all parameters via
        f1->SetParLimits(p_number, parmin, parmax);
   if parmin>=parmax, the parameter is fixed
   Note that you are not forced to fix the limits for all parameters.
   For example, if you fit a function with 6 parameters, you can do:
     func->SetParameters(0,3.1,1.e-6,0.1,-8,100);
     func->SetParLimits(4,-10,-4);
     func->SetParLimits(5, 1,1);
   With this setup, parameters 0->3 can vary freely
   Parameter 4 has boundaries [-10,-4] with initial value -8
   Parameter 5 is fixed to 100.

  Fit range
  =========
  The fit range can be specified in two ways:
    - specify rxmax > rxmin (default is rxmin=rxmax=0)
    - specify the option "R". In this case, the function will be taken
      instead of the full graph range.

   Changing the fitting function
   =============================
  By default the fitting function GraphFitChisquare is used.
  To specify a User defined fitting function, specify option "U" and
  call the following functions:
    TVirtualFitter::Fitter(mygraph)->SetFCN(MyFittingFunction)
  where MyFittingFunction is of type:
  extern void MyFittingFunction(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag);

  How errors are used in the chisquare function (see TFitter GraphFitChisquare)//   Access to the fit results
   ============================================
 In case of a TGraphErrors object, ex, the error along x,  is projected
 along the y-direction by calculating the function at the points x-exlow and
 x+exhigh.

 The chisquare is computed as the sum of the quantity below at each point:

                     (y - f(x))**2
         -----------------------------------
         ey**2 + ((f(x+exhigh) - f(x-exlow))/2)**2

 where x and y are the point coordinates.

 In case the function lies below (above) the data point, ey is ey_low (ey_high).

  thanks to Andy Haas (haas@yahoo.com) for adding the case with TGraphasymmerrors
            University of Washington

 a little different approach to approximating the uncertainty in y because of the
 errors in x, is to make it equal the error in x times the slope of the line.
 The improvement, compared to the first method (f(x+ exhigh) - f(x-exlow))/2
 is of (error of x)**2 order. This approach is called "effective variance method".
 This improvement has been made in version 4.00/08 by Anna Kreshuk.

   Associated functions
   ====================
  One or more object (typically a TF1*) can be added to the list
  of functions (fFunctions) associated to each graph.
  When TGraph::Fit is invoked, the fitted function is added to this list.
  Given a graph gr, one can retrieve an associated function
  with:  TF1 *myfunc = gr->GetFunction("myfunc");

  If the graph is made persistent, the list of
  associated functions is also persistent. Given a pointer (see above)
  to an associated function myfunc, one can retrieve the function/fit
  parameters with calls such as:
    Double_t chi2 = myfunc->GetChisquare();
    Double_t par0 = myfunc->GetParameter(0); //value of 1st parameter
    Double_t err0 = myfunc->GetParError(0);  //error on first parameter

   Fit Statistics
   ==============
  You can change the statistics box to display the fit parameters with
  the TStyle::SetOptFit(mode) method. This mode has four digits.
  mode = pcev  (default = 0111)
    v = 1;  print name/values of parameters
    e = 1;  print errors (if e=1, v must be 1)
    c = 1;  print Chisquare/Number of degress of freedom
    p = 1;  print Probability

  For example: gStyle->SetOptFit(1011);
  prints the fit probability, parameter names/values, and errors.
  You can change the position of the statistics box with these lines
  (where g is a pointer to the TGraph):

  Root > TPaveStats *st = (TPaveStats*)g->GetListOfFunctions()->FindObject("stats")
  Root > st->SetX1NDC(newx1); //new x start position
  Root > st->SetX2NDC(newx2); //new x end position

Option_t* GetGraphDrawOption(const TGraph *gr) const
 Return the draw option for the TGraph gr in this TMultiGraph
 The return option is the one specified when calling TMultiGraph::Add(gr,option).

void InitGaus(Double_t xmin, Double_t xmax)
*-*-*-*-*-*Compute Initial values of parameters for a gaussian*-*-*-*-*-*-*
*-*        ===================================================

void InitExpo(Double_t xmin, Double_t xmax)
*-*-*-*-*-*Compute Initial values of parameters for an exponential*-*-*-*-*
*-*        =======================================================

void InitPolynom(Double_t xmin, Double_t xmax)
*-*-*-*-*-*Compute Initial values of parameters for a polynom*-*-*-*-*-*-*
*-*        ===================================================

void LeastSquareFit(Int_t m, Double_t *a, Double_t xmin, Double_t xmax)
*-*-*-*-*-*-*-*Least squares lpolynomial fitting without weights*-*-*-*-*-*-*
*-*            =================================================

  m     number of parameters
  a     array of parameters
  first 1st point number to fit (default =0)
  last  last point number to fit (default=fNpoints-1)

   based on CERNLIB routine LSQ: Translated to C++ by Rene Brun



void LeastSquareLinearFit(Int_t ndata, Double_t &a0, Double_t &a1, Int_t &ifail, Double_t xmin, Double_t xmax)
*-*-*-*-*-*-*-*-*-*Least square linear fit without weights*-*-*-*-*-*-*-*-*
*-*                =======================================

  Fit a straight line (a0 + a1*x) to the data in this graph.
  ndata:  number of points to fit
  first:  first point number to fit
  last:   last point to fit O(ndata should be last-first
  ifail:  return parameter indicating the status of the fit (ifail=0, fit is OK)

   extracted from CERNLIB LLSQ: Translated to C++ by Rene Brun


TH1F* GetHistogram() const
    Returns a pointer to the histogram used to draw the axis
    Takes into account the two following cases.
       1- option 'A' was specified in TMultiGraph::Draw. Return fHistogram
       2- user had called TPad::DrawFrame. return pointer to hframe histogram

TF1* GetFunction(const char *name) const
*-*-*-*-*Return pointer to function with name*-*-*-*-*-*-*-*-*-*-*-*-*
*-*      ===================================

 Functions such as TGraph::Fit store the fitted function in the list of
 functions of this graph.

TAxis* GetXaxis() const
 Get x axis of the graph.

TAxis* GetYaxis() const
 Get y axis of the graph.

void Paint(Option_t *option)
 paint all the graphs of this multigraph

void Print(Option_t *option) const
 Print the list of graphs

void RecursiveRemove(TObject *obj)
 Recursively remove this object from a list. Typically implemented
 by classes that can contain mulitple references to a same object.

void SavePrimitive(ofstream &out, Option_t *option)
 Save primitive as a C++ statement(s) on output stream out

void SetMaximum(Double_t maximum)

void SetMinimum(Double_t minimum)



Inline Functions


              TList* GetListOfGraphs() const
              TList* GetListOfFunctions() const
             TClass* Class()
             TClass* IsA() const
                void ShowMembers(TMemberInspector& insp, char* parent)
                void Streamer(TBuffer& b)
                void StreamerNVirtual(TBuffer& b)
         TMultiGraph TMultiGraph(const TMultiGraph&)
        TMultiGraph& operator=(const TMultiGraph&)


Author: Rene Brun 12/10/2000
Last update: root/graf:$Name: $:$Id: TMultiGraph.cxx,v 1.25 2005/09/05 07:25:22 brun Exp $
Copyright (C) 1995-2000, Rene Brun and Fons Rademakers. *


ROOT page - Class index - Class Hierarchy - Top of the page

This page has been automatically generated. If you have any comments or suggestions about the page layout send a mail to ROOT support, or contact the developers with any questions or problems regarding ROOT.