TF2


class description - source file - inheritance tree

class TF2 : public TF1


    public:
TF2 TF2() TF2 TF2(const char* name, const char* formula, Double_t xmin = 0, Double_t xmax = 1, Double_t ymin = 0, Double_t ymax = 1) TF2 TF2(const char* name, void* fcn, Double_t xmin = 0, Double_t xmax = 1, Double_t ymin = 0, Double_t ymax = 1, Int_t npar = 0) TF2 TF2(const char* name, Double_t (*)(Double_t*, Double_t*) fcn, Double_t xmin = 0, Double_t xmax = 1, Double_t ymin = 0, Double_t ymax = 1, Int_t npar = 0) TF2 TF2(const TF2& f2) virtual void ~TF2() static TClass* Class() virtual void Copy(TObject& f2) const virtual Int_t DistancetoPrimitive(Int_t px, Int_t py) virtual void Draw(Option_t* option) virtual TF1* DrawCopy(Option_t* option) const virtual void DrawDerivative(Option_t* = "al") virtual void DrawF2(const char* formula, Double_t xmin, Double_t xmax, Double_t ymin, Double_t ymax, Option_t* option) virtual void DrawIntegral(Option_t* = "al") virtual void ExecuteEvent(Int_t event, Int_t px, Int_t py) virtual Int_t GetContour(Double_t* levels = 0) virtual Double_t GetContourLevel(Int_t level) const Int_t GetNpy() const virtual char* GetObjectInfo(Int_t px, Int_t py) const virtual Double_t GetRandom() virtual Double_t GetRandom(Double_t xmin, Double_t xmax) virtual void GetRandom2(Double_t& xrandom, Double_t& yrandom) virtual void GetRange(Double_t& xmin, Double_t& xmax) const virtual void GetRange(Double_t& xmin, Double_t& ymin, Double_t& xmax, Double_t& ymax) const virtual void GetRange(Double_t& xmin, Double_t& ymin, Double_t& zmin, Double_t& xmax, Double_t& ymax, Double_t& zmax) const virtual Double_t GetSave(const Double_t* x) virtual Double_t GetYmax() const virtual Double_t GetYmin() const virtual Double_t Integral(Double_t a, Double_t b, const Double_t* params = 0, Double_t epsil = 0.000001) virtual Double_t Integral(Double_t ax, Double_t bx, Double_t ay, Double_t by, Double_t epsil = 0.000001) virtual Double_t Integral(Double_t ax, Double_t bx, Double_t ay, Double_t by, Double_t az, Double_t bz, Double_t epsil = 0.000001) virtual TClass* IsA() const virtual Bool_t IsInside(const Double_t* x) const virtual void Paint(Option_t* option) virtual void Save(Double_t xmin, Double_t xmax, Double_t ymin, Double_t ymax, Double_t zmin, Double_t zmax) virtual void SavePrimitive(ofstream& out, Option_t* option) virtual void SetContour(Int_t nlevels = 20, const Double_t* levels = 0) virtual void SetContourLevel(Int_t level, Double_t value) virtual void SetNpy(Int_t npy = 100) virtual void SetRange(Double_t xmin, Double_t xmax) virtual void SetRange(Double_t xmin, Double_t ymin, Double_t xmax, Double_t ymax) virtual void SetRange(Double_t xmin, Double_t ymin, Double_t zmin, Double_t xmax, Double_t ymax, Double_t zmax) virtual void ShowMembers(TMemberInspector& insp, char* parent) virtual void Streamer(TBuffer& b) void StreamerNVirtual(TBuffer& b)

Data Members


    protected:
Double_t fYmin Lower bound for the range in y Double_t fYmax Upper bound for the range in y Int_t fNpy Number of points along y used for the graphical representation TArrayD fContour Array to display contour levels


See also

TF3

Class Description

 a 2-Dim function with parameters
 TF2 graphics function is via the TH1 drawing functions.

      Example of a function

   TF2 *f2 = new TF2("f2","sin(x)*sin(y)/(x*y)",0,5,0,5);
   f2->Draw();

/*

*/


      See TF1 class for the list of functions formats


TF2(): TF1()
*-*-*-*-*-*-*-*-*-*-*F2 default constructor*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
*-*                  ======================

TF2(const char *name,const char *formula, Double_t xmin, Double_t xmax, Double_t ymin, Double_t ymax) :TF1(name,formula,xmin,xmax)
*-*-*-*-*-*-*F2 constructor using a formula definition*-*-*-*-*-*-*-*-*-*-*
*-*          =========================================
*-*
*-*  See TFormula constructor for explanation of the formula syntax.
*-*
*-*
*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*

TF2(const char *name, void *fcn, Double_t xmin, Double_t xmax, Double_t ymin, Double_t ymax, Int_t npar) : TF1(name, fcn, xmin, xmax, npar)
*-*-*-*-*-*-*F2 constructor using a pointer to an interpreted function*-*-*
*-*          =========================================================
*-*
*-*   npar is the number of free parameters used by the function
*-*
*-*  Creates a function of type C between xmin and xmax and ymin,ymax.
*-*  The function is defined with npar parameters
*-*  fcn must be a function of type:
*-*     Double_t fcn(Double_t *x, Double_t *params)
*-*
*-*  This constructor is called for functions of type C by CINT.
*-*
*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*

TF2(const char *name, Double_t (*fcn)(Double_t *, Double_t *), Double_t xmin, Double_t xmax, Double_t ymin, Double_t ymax, Int_t npar) : TF1(name, fcn, xmin, xmax, npar)
*-*-*-*-*-*-*F2 constructor using a pointer to a compiled function*-*-*-*-*
*-*          =====================================================
*-*
*-*   npar is the number of free parameters used by the function
*-*
*-*   This constructor creates a function of type C when invoked
*-*   with the normal C++ compiler.
*-*
*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*

~TF2()
*-*-*-*-*-*-*-*-*-*-*F2 default destructor*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
*-*                  =====================

TF2(const TF2 &f2) : TF1(f2)

void Copy(TObject &obj) const
*-*-*-*-*-*-*-*-*-*-*Copy this F2 to a new F2*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
*-*                  ========================

Int_t DistancetoPrimitive(Int_t px, Int_t py)
*-*-*-*-*-*-*-*-*-*-*Compute distance from point px,py to a function*-*-*-*-*
*-*                  ===============================================
*-*  Compute the closest distance of approach from point px,py to this function.
*-*  The distance is computed in pixels units.
*-*
*-*  Algorithm:
*-*
*-*
*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*

void Draw(Option_t *option)
*-*-*-*-*-*-*-*-*-*-*Draw this function with its current attributes*-*-*-*-*
*-*                  ==============================================
*-* NB. You must use DrawCopy if you want to draw several times the same
*-*     function in the current canvas.
*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*

TF1* DrawCopy(Option_t *option) const
*-*-*-*-*-*-*-*Draw a copy of this function with its current attributes*-*-*
*-*            ========================================================
*-*
*-*  This function MUST be used instead of Draw when you want to draw
*-*  the same function with different parameters settings in the same canvas.
*-*
*-* Possible option values are:
*-*   "SAME"  superimpose on top of existing picture
*-*   "L"     connect all computed points with a straight line
*-*   "C"     connect all computed points with a smooth curve.
*-*
*-* Note that the default value is "F". Therefore to draw on top
*-* of an existing picture, specify option "SL"
*-*
*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*

void DrawF2(const char *formula, Double_t xmin, Double_t ymin, Double_t xmax, Double_t ymax, Option_t *option)
*-*-*-*-*-*-*-*-*-*Draw formula between xmin,ymin and xmax,ymax*-*-*-*-*-*-*-*
*-*                ============================================
*-*

void ExecuteEvent(Int_t event, Int_t px, Int_t py)
*-*-*-*-*-*-*-*-*-*-*Execute action corresponding to one event*-*-*-*
*-*                  =========================================
*-*  This member function is called when a F2 is clicked with the locator
*-*
*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*

Int_t GetContour(Double_t *levels)
*-*-*-*-*-*-*-*Return contour values into array levels*-*-*-*-*-*-*-*-*-*
*-*            =======================================
*-*
*-*  The number of contour levels can be returned by getContourLevel
*-*
*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*

Double_t GetContourLevel(Int_t level) const
*-*-*-*-*-*-*-*Return the number of contour levels*-*-*-*-*-*-*-*-*-*-*-*-*
*-*            ===================================

char* GetObjectInfo(Int_t px, Int_t py) const
   Redefines TObject::GetObjectInfo.
   Displays the function value
   corresponding to cursor position px,py


Double_t GetRandom()
*-*-*-*-*-*Return a random number following this function shape*-*-*-*-*-*-*
*-*        ====================================================
*-*

Double_t GetRandom(Double_t, Double_t)
*-*-*-*-*-*Return a random number following this function shape*-*-*-*-*-*-*
*-*        ====================================================
*-*

void GetRandom2(Double_t &xrandom, Double_t &yrandom)
*-*-*-*-*-*Return 2 random numbers following this function shape*-*-*-*-*-*
*-*        =====================================================
*-*
*-*   The distribution contained in this TF2 function is integrated
*-*   over the cell contents.
*-*   It is normalized to 1.
*-*   Getting the two random numbers implies:
*-*     - Generating a random number between 0 and 1 (say r1)
*-*     - Look in which cell in the normalized integral r1 corresponds to
*-*     - make a linear interpolation in the returned cell
*-*

void GetRange(Double_t &xmin, Double_t &ymin, Double_t &xmax, Double_t &ymax) const
*-*-*-*-*-*-*-*-*-*-*Return range of a 2-D function*-*-*-*-*-*-*-*-*-*-*-*-*
*-*                  ==============================

void GetRange(Double_t &xmin, Double_t &ymin, Double_t &zmin, Double_t &xmax, Double_t &ymax, Double_t &zmax) const
*-*-*-*-*-*-*-*-*-*-*Return range of function*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
*-*                  ========================

Double_t GetSave(const Double_t *xx)
 Get value corresponding to X in array of fSave values

Double_t Integral(Double_t ax, Double_t bx, Double_t ay, Double_t by, Double_t epsilon)
 Return Integral of a 2d function in range [ax,bx],[ay,by]


Bool_t IsInside(const Double_t *x) const
 Return kTRUE is the point is inside the function range

void Paint(Option_t *option)
*-*-*-*-*-*-*-*-*Paint this 2-D function with its current attributes*-*-*-*-*
*-*              ===================================================

void Save(Double_t xmin, Double_t xmax, Double_t ymin, Double_t ymax, Double_t, Double_t)
 Save values of function in array fSave

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

void SetContour(Int_t nlevels, const Double_t *levels)
*-*-*-*-*-*-*-*Set the number and values of contour levels*-*-*-*-*-*-*-*-*
*-*            ===========================================

  By default the number of contour levels is set to 20.

  if argument levels = 0 or missing, equidistant contours are computed


void SetContourLevel(Int_t level, Double_t value)
*-*-*-*-*-*-*-*-*-*-*Set value for one contour level*-*-*-*-*-*-*-*-*-*-*-*
*-*                  ===============================

void SetNpy(Int_t npy)
*-*-*-*-*-*-*-*Set the number of points used to draw the function*-*-*-*-*-*
*-*            ==================================================

void SetRange(Double_t xmin, Double_t ymin, Double_t xmax, Double_t ymax)
*-*-*-*-*-*Initialize the upper and lower bounds to draw the function*-*-*-*
*-*        ==========================================================

void Streamer(TBuffer &R__b)
 Stream an object of class TF2.



Inline Functions


               void DrawDerivative(Option_t* = "al")
               void DrawIntegral(Option_t* = "al")
              Int_t GetNpy() const
               void GetRange(Double_t& xmin, Double_t& ymin, Double_t& zmin, Double_t& xmax, Double_t& ymax, Double_t& zmax) const
           Double_t GetYmin() const
           Double_t GetYmax() const
           Double_t Integral(Double_t ax, Double_t bx, Double_t ay, Double_t by, Double_t epsil = 0.000001)
           Double_t Integral(Double_t ax, Double_t bx, Double_t ay, Double_t by, Double_t az, Double_t bz, Double_t epsil = 0.000001)
               void SetRange(Double_t xmin, Double_t ymin, Double_t xmax, Double_t ymax)
               void SetRange(Double_t xmin, Double_t ymin, Double_t zmin, Double_t xmax, Double_t ymax, Double_t zmax)
            TClass* Class()
            TClass* IsA() const
               void ShowMembers(TMemberInspector& insp, char* parent)
               void StreamerNVirtual(TBuffer& b)


Author: Rene Brun 23/08/95
Last update: root/hist:$Name: $:$Id: TF2.cxx,v 1.20 2003/04/20 20:03:04 brun Exp $
Copyright (C) 1995-2000, Rene Brun and Fons Rademakers. *


ROOT page - Class index - 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.