// @(#)root/hist:$Id$
// Author: Rene Brun   23/08/95

/*************************************************************************
 * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers.               *
 * All rights reserved.                                                  *
 *                                                                       *
 * For the licensing terms see $ROOTSYS/LICENSE.                         *
 * For the list of contributors see $ROOTSYS/README/CREDITS.             *
 *************************************************************************/

#include "TROOT.h"
#include "TF2.h"
#include "TMath.h"
#include "TRandom.h"
#include "TH2.h"
#include "TVirtualPad.h"
#include "TStyle.h"
#include "Riostream.h"
#include "TColor.h"
#include "TVirtualFitter.h"
#include "TClass.h"

ClassImp(TF2)

//______________________________________________________________________________
//
// 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();
//Begin_Html
/*
<img src="gif/function2.gif">
*/
//End_Html
//
//      See TF1 class for the list of functions formats
//

//______________________________________________________________________________
TF2::TF2(): TF1(),fYmin(0),fYmax(0),fNpy(100)
{
//*-*-*-*-*-*-*-*-*-*-*F2 default constructor*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
//*-*                  ======================

}


//______________________________________________________________________________
TF2::TF2(const char *name,const char *formula, Double_t xmin, Double_t xmax, Double_t ymin, Double_t ymax)
      :TF1(name,formula,xmax,xmin)
{
// F2 constructor using a formula definition
//
//  See TFormula constructor for explanation of the formula syntax.
//
//  if formula has the form "fffffff;xxxx;yyyy", it is assumed that
//  the formula string is "fffffff" and "xxxx" and "yyyy" are the
//  titles for the X and Y axis respectively.

   if (ymin < ymax) {
      fYmin   = ymin;
      fYmax   = ymax;
   } else {
      fYmin = ymax;
      fYmax = ymin;
   }
   fNpx    = 30;
   fNpy    = 30;
   fContour.Set(0);
   if (GetNdim() != 2 && xmin < xmax && ymin < ymax) {
      Error("TF2","function: %s/%s has dimension %d instead of 2",name,formula,GetNdim());
      MakeZombie();
   }
}


//______________________________________________________________________________
TF2::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, Int_t ndim)
      : TF1(name, fcn, xmin, xmax, npar,ndim)
{
//*-*-*-*-*-*-*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.
//*-*
//*-* WARNING! A function created with this constructor cannot be Cloned.
//*-*
//*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*

   fYmin   = ymin;
   fYmax   = ymax;
   fNpx    = 30;
   fNpy    = 30;
   fContour.Set(0);

}

//______________________________________________________________________________
TF2::TF2(const char *name, Double_t (*fcn)(const Double_t *, const Double_t *), Double_t xmin, Double_t xmax, Double_t ymin, Double_t ymax, Int_t npar, Int_t ndim)
      : TF1(name, fcn, xmin, xmax, npar,ndim)
{
//*-*-*-*-*-*-*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.
//*-*
//*-* WARNING! A function created with this constructor cannot be Cloned.
//*-*
//*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*

   fYmin   = ymin;
   fYmax   = ymax;
   fNpx    = 30;
   fNpy    = 30;
   fContour.Set(0);

}

//______________________________________________________________________________
TF2::TF2(const char *name, ROOT::Math::ParamFunctor f, Double_t xmin, Double_t xmax, Double_t ymin, Double_t ymax, Int_t npar, Int_t ndim)
      : TF1(name, f, xmin, xmax, npar,ndim)
{
//*-*-*-*-*-*-*F2 constructor using a ParamFunctor, 
//*-*          a functor class implementing operator() (double *, double *)  
//*-*
//*-*   npar is the number of free parameters used by the function
//*-*
//*-* WARNING! A function created with this constructor cannot be Cloned.
//*-*

   fYmin   = ymin;
   fYmax   = ymax;
   fNpx    = 30;
   fNpy    = 30;
   fContour.Set(0);

}

//______________________________________________________________________________
TF2& TF2::operator=(const TF2 &rhs)
{
   // Operator =

   if (this != &rhs) {
      rhs.Copy(*this);
   }
   return *this;
}

//______________________________________________________________________________
TF2::~TF2()
{
//*-*-*-*-*-*-*-*-*-*-*F2 default destructor*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
//*-*                  =====================

}

//______________________________________________________________________________
TF2::TF2(const TF2 &f2) : TF1()
{
   // Copy constructor.

   ((TF2&)f2).Copy(*this);
}

//______________________________________________________________________________
void TF2::Copy(TObject &obj) const
{
//*-*-*-*-*-*-*-*-*-*-*Copy this F2 to a new F2*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
//*-*                  ========================

   TF1::Copy(obj);
   ((TF2&)obj).fYmin    = fYmin;
   ((TF2&)obj).fYmax    = fYmax;
   ((TF2&)obj).fNpy     = fNpy;
   fContour.Copy(((TF2&)obj).fContour);
}

//______________________________________________________________________________
Int_t TF2::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:
//*-*
//*-*
//*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*

   if (!fHistogram) return 9999;
   Int_t distance = fHistogram->DistancetoPrimitive(px,py);
   if (distance <= 1) return distance;

   Double_t x = gPad->PadtoX(gPad->AbsPixeltoX(px));
   Double_t y = gPad->PadtoY(gPad->AbsPixeltoY(py));
   const char *drawOption = GetDrawOption();
   Double_t uxmin,uxmax;
   Double_t uymin,uymax;
   if (gPad->GetView() || strncmp(drawOption,"cont",4) == 0
                       || strncmp(drawOption,"CONT",4) == 0) {
      uxmin=gPad->GetUxmin();
      uxmax=gPad->GetUxmax();
      x = fXmin +(fXmax-fXmin)*(x-uxmin)/(uxmax-uxmin);
      uymin=gPad->GetUymin();
      uymax=gPad->GetUymax();
      y = fYmin +(fYmax-fYmin)*(y-uymin)/(uymax-uymin);
   }
   if (x < fXmin || x > fXmax) return distance;
   if (y < fYmin || y > fYmax) return distance;
   return 0;
}

//______________________________________________________________________________
void TF2::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.
//*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*

   TString opt = option;
   opt.ToLower();
   if (gPad && !opt.Contains("same")) gPad->Clear();

   AppendPad(option);
}

//______________________________________________________________________________
TF1 *TF2::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"
//*-*
//*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*

   TF2 *newf2 = new TF2();
   Copy(*newf2);
   newf2->AppendPad(option);
   newf2->SetBit(kCanDelete);
   return newf2;
}

// remove this function
//______________________________________________________________________________
// void TF2::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*-*-*-*-*-*-*-*
// //*-*                ============================================
// //*-*

//    //if (Compile((char*)formula)) return;

//    SetRange(xmin, ymin, xmax, ymax);

//    Draw(option);

// }

//______________________________________________________________________________
void TF2::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
//*-*
//*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*

   TF1::ExecuteEvent(event, px, py);
}

//______________________________________________________________________________
Int_t TF2::GetContour(Double_t *levels)
{
//*-*-*-*-*-*-*-*Return contour values into array levels*-*-*-*-*-*-*-*-*-*
//*-*            =======================================
//*-*
//*-*  The number of contour levels can be returned by getContourLevel
//*-*
//*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
   Int_t nlevels = fContour.fN;
   if (levels) {
      for (Int_t level=0; level<nlevels; level++) levels[level] = GetContourLevel(level);
   }
   return nlevels;
}

//______________________________________________________________________________
Double_t TF2::GetContourLevel(Int_t level) const
{
//*-*-*-*-*-*-*-*Return the number of contour levels*-*-*-*-*-*-*-*-*-*-*-*-*
//*-*            ===================================
   if (level <0 || level >= fContour.fN) return 0;
   if (fContour.fArray[0] != -9999) return fContour.fArray[level];
   if (fHistogram == 0) return 0;
   return fHistogram->GetContourLevel(level);
}

//______________________________________________________________________________
Double_t TF2::FindMinMax(Double_t *x, Bool_t findmax) const
{
// return minimum/maximum value of the function
// To find the minimum on a range, first set this range via the SetRange function
// If a vector x of coordinate is passed it will be used as starting point for the minimum. 
// In addition on exit x will contain the coordinate values at the minimuma
// If x is NULL or x is inifinity or NaN, first, a grid search is performed to find the initial estimate of the 
// minimum location. The range of the function is divided into fNpx and fNpy
// sub-ranges. If the function is "good" (or "bad"), these values can be changed
// by SetNpx and SetNpy functions
// Then, a minimization is used with starting values found by the grid search
// The minimizer algorithm used (by default Minuit) can be changed by callinga
//  ROOT::Math::Minimizer::SetDefaultMinimizerType("..")
// Other option for the minimizer can be set using the static method of the MinimizerOptions class

   //First do a grid search with step size fNpx and fNpy

   Double_t xx[2]; 
   Double_t rsign = (findmax) ? -1. : 1.;
   TF2 & function = const_cast<TF2&>(*this); // needed since EvalPar is not const
   Double_t xxmin = 0, yymin = 0, zzmin = 0; 
   if (x == NULL || ( (x!= NULL) && ( !TMath::Finite(x[0]) || !TMath::Finite(x[1]) ) ) ){ 
      Double_t dx = (fXmax - fXmin)/fNpx;
      Double_t dy = (fYmax - fYmin)/fNpy;
      xxmin = fXmin;
      yymin = fYmin;
      zzmin = rsign * TMath::Infinity();
      for (Int_t i=0; i<fNpx; i++){
         xx[0]=fXmin + (i+0.5)*dx;
         for (Int_t j=0; j<fNpy; j++){
            xx[1]=fYmin+(j+0.5)*dy;
            Double_t zz = function(xx);
            if (rsign*zz < rsign*zzmin) {xxmin = xx[0], yymin = xx[1]; zzmin = zz;}
         }
      }
      
      xxmin = TMath::Min(fXmax, xxmin);
      yymin = TMath::Min(fYmax, yymin);
   }
   else {
      xxmin = x[0]; 
      yymin = x[1];
      zzmin = function(xx);
   }
   xx[0] = xxmin; 
   xx[1] = yymin; 
      
   double fmin = GetMinMaxNDim(xx,findmax);
   if (rsign*fmin < rsign*zzmin) { 
      if (x) {x[0] = xx[0]; x[1] = xx[1]; }
      return fmin;
   }
   // here if minimization failed
   if (x) { x[0] = xxmin; x[1] = yymin; }
   return zzmin;
}

//______________________________________________________________________________
Double_t TF2::GetMinimumXY(Double_t &x, Double_t &y) const
{
// Compute the X and Y values corresponding to the minimum value of the function
// Return the minimum value of the function
// To find the minimum on a range, first set this range via the SetRange function
// Method:
//   First, a grid search is performed to find the initial estimate of the 
//   minimum location. The range of the function is divided into fNpx and fNpy
//   sub-ranges. If the function is "good" (or "bad"), these values can be changed
//   by SetNpx and SetNpy functions
//   Then, a minimization is used with starting values found by the grid search
//   The minimizer algorithm used (by default Minuit) can be changed by callinga
//   ROOT::Math::Minimizer::SetDefaultMinimizerType("..")
//   Other option for the minimizer can be set using the static method of the MinimizerOptions class
//
//   Note that this method will always do first a grid search in contrast to GetMinimum

   double xx[2] = { 0,0 };
   xx[0] = TMath::QuietNaN();  // to force to do grid search in TF2::FindMinMax
   double fmin = FindMinMax(xx, false);
   x = xx[0]; y = xx[1];
   return fmin;
}

//______________________________________________________________________________
Double_t TF2::GetMaximumXY(Double_t &x, Double_t &y) const
{
// Compute the X and Y values corresponding to the maximum value of the function
// Return the maximum value of the function
//  See TF2::GetMinimumXY

   double xx[2] = { 0,0 };
   xx[0] = TMath::QuietNaN();  // to force to do grid search in TF2::FindMinMax
   double fmax = FindMinMax(xx, true);
   x = xx[0]; y = xx[1];
   return fmax;
}


//______________________________________________________________________________
Double_t TF2::GetMinimum(Double_t *x) const
{
// return minimum/maximum value of the function
// To find the minimum on a range, first set this range via the SetRange function
// If a vector x of coordinate is passed it will be used as starting point for the minimum. 
// In addition on exit x will contain the coordinate values at the minimuma
// If x is NULL or x is inifinity or NaN, first, a grid search is performed to find the initial estimate of the 
// minimum location. The range of the function is divided into fNpx and fNpy
// sub-ranges. If the function is "good" (or "bad"), these values can be changed
// by SetNpx and SetNpy functions
// Then, a minimization is used with starting values found by the grid search
// The minimizer algorithm used (by default Minuit) can be changed by callinga
//  ROOT::Math::Minimizer::SetDefaultMinimizerType("..")
// Other option for the minimizer can be set using the static method of the MinimizerOptions class

   return FindMinMax(x, false); 
}

//______________________________________________________________________________
Double_t TF2::GetMaximum(Double_t *x) const
{
// return maximum value of the function
// See TF2::GetMinimum

   return FindMinMax(x, true); 
}


//______________________________________________________________________________
char *TF2::GetObjectInfo(Int_t px, Int_t py) const
{
//   Redefines TObject::GetObjectInfo.
//   Displays the function value
//   corresponding to cursor position px,py
//
   const char *snull = "";
   if (!gPad) return (char*)snull;
   static char info[64];
   Double_t x = gPad->PadtoX(gPad->AbsPixeltoX(px));
   Double_t y = gPad->PadtoY(gPad->AbsPixeltoY(py));
   const char *drawOption = GetDrawOption();
   Double_t uxmin,uxmax;
   Double_t uymin,uymax;
   if (gPad->GetView() || strncmp(drawOption,"cont",4) == 0
                       || strncmp(drawOption,"CONT",4) == 0) {
      uxmin=gPad->GetUxmin();
      uxmax=gPad->GetUxmax();
      x = fXmin +(fXmax-fXmin)*(x-uxmin)/(uxmax-uxmin);
      uymin=gPad->GetUymin();
      uymax=gPad->GetUymax();
      y = fYmin +(fYmax-fYmin)*(y-uymin)/(uymax-uymin);
   }
   snprintf(info,64,"(x=%g, y=%g, f=%.18g)",x,y,((TF2*)this)->Eval(x,y));
   return info;
}

//______________________________________________________________________________
Double_t TF2::GetRandom()
{
//*-*-*-*-*-*Return a random number following this function shape*-*-*-*-*-*-*
//*-*        ====================================================
//*-*
   Error("GetRandom","cannot be called for TF2/3, use GetRandom2/3 instead");
   return 0;  // not yet implemented
}

//______________________________________________________________________________
Double_t TF2::GetRandom(Double_t, Double_t)
{
//*-*-*-*-*-*Return a random number following this function shape*-*-*-*-*-*-*
//*-*        ====================================================
//*-*
   Error("GetRandom","cannot be called for TF2/3, use GetRandom2/3 instead");
   return 0;  // not yet implemented
}

//______________________________________________________________________________
void TF2::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
//*-*
//*-*
//*-*  IMPORTANT NOTE
//*-*  The integral of the function is computed at fNpx * fNpy points. 
//*-*  If the function has sharp peaks, you should increase the number of 
//*-*  points (SetNpx, SetNpy) such that the peak is correctly tabulated 
//*-*  at several points.

   //  Check if integral array must be build
   Int_t i,j,cell;
   Double_t dx   = (fXmax-fXmin)/fNpx;
   Double_t dy   = (fYmax-fYmin)/fNpy;
   Int_t ncells = fNpx*fNpy;
   if (fIntegral.empty()) {
      fIntegral.resize(ncells+1);
      fIntegral[0] = 0;
      Double_t integ;
      Int_t intNegative = 0;
      cell = 0;
      for (j=0;j<fNpy;j++) {
         for (i=0;i<fNpx;i++) {
            integ = Integral(fXmin+i*dx,fXmin+i*dx+dx,fYmin+j*dy,fYmin+j*dy+dy);
            if (integ < 0) {intNegative++; integ = -integ;}
            fIntegral[cell+1] = fIntegral[cell] + integ;
            cell++;
         }
      }
      if (intNegative > 0) {
         Warning("GetRandom2","function:%s has %d negative values: abs assumed",GetName(),intNegative);
      }
      if (fIntegral[ncells] == 0) {
         Error("GetRandom2","Integral of function is zero");
         return;
      }
      for (i=1;i<=ncells;i++) {  // normalize integral to 1
         fIntegral[i] /= fIntegral[ncells];
      }
   }

// return random numbers
   Double_t r,ddx,ddy,dxint;
   r     = gRandom->Rndm();
   cell  = TMath::BinarySearch(ncells,fIntegral.data(),r);
   dxint = fIntegral[cell+1] - fIntegral[cell];
   if (dxint > 0) ddx = dx*(r - fIntegral[cell])/dxint;
   else           ddx = 0;
   ddy = dy*gRandom->Rndm();
   j   = cell/fNpx;
   i   = cell%fNpx;
   xrandom = fXmin +dx*i +ddx;
   yrandom = fYmin +dy*j +ddy;
}

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

   xmin = fXmin;
   xmax = fXmax;
   ymin = fYmin;
   ymax = fYmax;
}

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

   xmin = fXmin;
   xmax = fXmax;
   ymin = fYmin;
   ymax = fYmax;
   zmin = 0;
   zmax = 0;
}


//______________________________________________________________________________
Double_t TF2::GetSave(const Double_t *xx)
{
    // Get value corresponding to X in array of fSave values

   //if (fNsave <= 0) return 0;
   if (fSave.empty()) return 0;
   Int_t fNsave = fSave.size(); 
   Int_t np = fNsave - 6;
   Double_t xmin = Double_t(fSave[np+0]);
   Double_t xmax = Double_t(fSave[np+1]);
   Double_t ymin = Double_t(fSave[np+2]);
   Double_t ymax = Double_t(fSave[np+3]);
   Int_t npx     = Int_t(fSave[np+4]);
   Int_t npy     = Int_t(fSave[np+5]);
   Double_t x    = Double_t(xx[0]);
   Double_t dx   = (xmax-xmin)/npx;
   if (x < xmin || x > xmax) return 0;
   if (dx <= 0) return 0;
   Double_t y    = Double_t(xx[1]);
   Double_t dy   = (ymax-ymin)/npy;
   if (y < ymin || y > ymax) return 0;
   if (dy <= 0) return 0;

   //we make a bilinear interpolation using the 4 points surrounding x,y
   Int_t ibin    = Int_t((x-xmin)/dx);
   Int_t jbin    = Int_t((y-ymin)/dy);
   Double_t xlow = xmin + ibin*dx;
   Double_t ylow = ymin + jbin*dy;
   Double_t t    = (x-xlow)/dx;
   Double_t u    = (y-ylow)/dy;
   Int_t k1      = jbin*(npx+1) + ibin;
   Int_t k2      = jbin*(npx+1) + ibin +1;
   Int_t k3      = (jbin+1)*(npx+1) + ibin +1;
   Int_t k4      = (jbin+1)*(npx+1) + ibin;
   Double_t z    = (1-t)*(1-u)*fSave[k1] +t*(1-u)*fSave[k2] +t*u*fSave[k3] + (1-t)*u*fSave[k4];
   return z;
}

//______________________________________________________________________________
Double_t TF2::Integral(Double_t ax, Double_t bx, Double_t ay, Double_t by, Double_t epsrel)
{
// Return Integral of a 2d function in range [ax,bx],[ay,by]
// with desired relative accuracy (default value of eps is 1.e-9)
//
   Double_t a[2], b[2];
   a[0] = ax;
   b[0] = bx;
   a[1] = ay;
   b[1] = by;
   Double_t relerr  = 0;
   Int_t n = 2;
   Int_t maxpts = 20*fNpx*fNpy;
   Int_t nfnevl,ifail;
   Double_t result = IntegralMultiple(n,a,b,maxpts,epsrel,epsrel,relerr,nfnevl,ifail);
   if (ifail > 0) {
      Warning("Integral","failed code=%d, maxpts=%d, epsrel=%g, nfnevl=%d, relerr=%g ",ifail,maxpts,epsrel,nfnevl,relerr);
   }
   return result;
}

//______________________________________________________________________________
Bool_t TF2::IsInside(const Double_t *x) const
{
// Return kTRUE is the point is inside the function range

   if (x[0] < fXmin || x[0] > fXmax) return kFALSE;
   if (x[1] < fYmin || x[1] > fYmax) return kFALSE;
   return kTRUE;
}

//______________________________________________________________________________
TH1* TF2::CreateHistogram()
{
   // Create a histogram from function.
   // always created it, even if it is already existing

   Int_t i,j,bin;
   Double_t dx, dy;
   Double_t xv[2];


   Double_t *parameters = GetParameters();
   TH2F* h = new TH2F("Func",(char*)GetTitle(),fNpx,fXmin,fXmax,fNpy,fYmin,fYmax);
   h->SetDirectory(0);

   InitArgs(xv,parameters);
   dx = (fXmax - fXmin)/Double_t(fNpx);
   dy = (fYmax - fYmin)/Double_t(fNpy);
   for (i=1;i<=fNpx;i++) {
      xv[0] = fXmin + (Double_t(i) - 0.5)*dx;
      for (j=1;j<=fNpy;j++) {
         xv[1] = fYmin + (Double_t(j) - 0.5)*dy;
         bin   = j*(fNpx + 2) + i;
         h->SetBinContent(bin,EvalPar(xv,parameters));
      }
   }
   h->Fill(fXmin-1,fYmin-1,0);  //This call to force fNentries non zero

   Double_t *levels = fContour.GetArray();
   if (levels && levels[0] == -9999) levels = 0;
   h->SetMinimum(fMinimum);
   h->SetMaximum(fMaximum);
   h->SetContour(fContour.fN, levels);
   h->SetLineColor(GetLineColor());
   h->SetLineStyle(GetLineStyle());
   h->SetLineWidth(GetLineWidth());
   h->SetFillColor(GetFillColor());
   h->SetFillStyle(GetFillStyle());
   h->SetMarkerColor(GetMarkerColor());
   h->SetMarkerStyle(GetMarkerStyle());
   h->SetMarkerSize(GetMarkerSize());
   h->SetStats(0);

   return h;
}

//______________________________________________________________________________
void TF2::Paint(Option_t *option)
{
//*-*-*-*-*-*-*-*-*Paint this 2-D function with its current attributes*-*-*-*-*
//*-*              ===================================================

   Int_t i,j,bin;
   Double_t dx, dy;
   Double_t xv[2];
   Double_t *parameters = GetParameters();
   TString opt = option;
   opt.ToLower();

//*-*-  Create a temporary histogram and fill each channel with the function value
   if (!fHistogram) {
      fHistogram = new TH2F("Func",(char*)GetTitle(),fNpx,fXmin,fXmax,fNpy,fYmin,fYmax);
      if (!fHistogram) return;
      fHistogram->SetDirectory(0);
   }
   InitArgs(xv,parameters);
   dx = (fXmax - fXmin)/Double_t(fNpx);
   dy = (fYmax - fYmin)/Double_t(fNpy);
   for (i=1;i<=fNpx;i++) {
      xv[0] = fXmin + (Double_t(i) - 0.5)*dx;
      for (j=1;j<=fNpy;j++) {
         xv[1] = fYmin + (Double_t(j) - 0.5)*dy;
         bin   = j*(fNpx + 2) + i;
         fHistogram->SetBinContent(bin,EvalPar(xv,parameters));
      }
   }
   ((TH2F*)fHistogram)->Fill(fXmin-1,fYmin-1,0);  //This call to force fNentries non zero

//*-*- Copy Function attributes to histogram attributes
   Double_t *levels = fContour.GetArray();
   if (levels && levels[0] == -9999) levels = 0;
   fHistogram->SetMinimum(fMinimum);
   fHistogram->SetMaximum(fMaximum);
   fHistogram->SetContour(fContour.fN, levels);
   fHistogram->SetLineColor(GetLineColor());
   fHistogram->SetLineStyle(GetLineStyle());
   fHistogram->SetLineWidth(GetLineWidth());
   fHistogram->SetFillColor(GetFillColor());
   fHistogram->SetFillStyle(GetFillStyle());
   fHistogram->SetMarkerColor(GetMarkerColor());
   fHistogram->SetMarkerStyle(GetMarkerStyle());
   fHistogram->SetMarkerSize(GetMarkerSize());
   fHistogram->SetStats(0);

//*-*-  Draw the histogram
   if (!gPad) return;
   if (opt.Length() == 0)  fHistogram->Paint("cont3");
   else if (opt == "same") fHistogram->Paint("cont2same");
   else                    fHistogram->Paint(option);
}

//______________________________________________________________________________
void TF2::Save(Double_t xmin, Double_t xmax, Double_t ymin, Double_t ymax, Double_t, Double_t)
{
    // Save values of function in array fSave
   if (!fSave.empty()) fSave.clear(); 
   //if (fSave != 0) {delete [] fSave; fSave = 0;}
   Int_t nsave = (fNpx+1)*(fNpy+1);
   Int_t fNsave = nsave+6;
   if (fNsave <= 6) {fNsave=0; return;}
   //fSave  = new Double_t[fNsave];
   fSave.resize(fNsave);
   Int_t i,j,k=0;
   Double_t dx = (xmax-xmin)/fNpx;
   Double_t dy = (ymax-ymin)/fNpy;
   if (dx <= 0) {
      dx = (fXmax-fXmin)/fNpx;
      xmin = fXmin +0.5*dx;
      xmax = fXmax -0.5*dx;
   }
   if (dy <= 0) {
      dy = (fYmax-fYmin)/fNpy;
      ymin = fYmin +0.5*dy;
      ymax = fYmax -0.5*dy;
   }
   Double_t xv[2];
   Double_t *parameters = GetParameters();
   InitArgs(xv,parameters);
   for (j=0;j<=fNpy;j++) {
      xv[1]    = ymin + dy*j;
      for (i=0;i<=fNpx;i++) {
         xv[0]    = xmin + dx*i;
         fSave[k] = EvalPar(xv,parameters);
         k++;
      }
   }
   fSave[nsave+0] = xmin;
   fSave[nsave+1] = xmax;
   fSave[nsave+2] = ymin;
   fSave[nsave+3] = ymax;
   fSave[nsave+4] = fNpx;
   fSave[nsave+5] = fNpy;
}

//______________________________________________________________________________
void TF2::SavePrimitive(std::ostream &out, Option_t *option /*= ""*/)
{
   // Save primitive as a C++ statement(s) on output stream out

   char quote = '"';
   out<<"   "<<std::endl;
   if (gROOT->ClassSaved(TF2::Class())) {
      out<<"   ";
   } else {
      out<<"   TF2 *";
   }
   if (!fMethodCall) {
      out<<GetName()<<" = new TF2("<<quote<<GetName()<<quote<<","<<quote<<GetTitle()<<quote<<","<<fXmin<<","<<fXmax<<","<<fYmin<<","<<fYmax<<");"<<std::endl;
   } else {
      out<<GetName()<<" = new TF2("<<quote<<GetName()<<quote<<","<<GetTitle()<<","<<fXmin<<","<<fXmax<<","<<fYmin<<","<<fYmax<<","<<GetNpar()<<");"<<std::endl;
   }

   if (GetFillColor() != 0) {
      if (GetFillColor() > 228) {
         TColor::SaveColor(out, GetFillColor());
         out<<"   "<<GetName()<<"->SetFillColor(ci);" << std::endl;
      } else 
         out<<"   "<<GetName()<<"->SetFillColor("<<GetFillColor()<<");"<<std::endl;
   }
   if (GetFillStyle() != 1001) {
      out<<"   "<<GetName()<<"->SetFillStyle("<<GetFillStyle()<<");"<<std::endl;
   }
   if (GetMarkerColor() != 1) {
      if (GetMarkerColor() > 228) {
         TColor::SaveColor(out, GetMarkerColor());
         out<<"   "<<GetName()<<"->SetMarkerColor(ci);" << std::endl;
      } else 
         out<<"   "<<GetName()<<"->SetMarkerColor("<<GetMarkerColor()<<");"<<std::endl;
   }
   if (GetMarkerStyle() != 1) {
      out<<"   "<<GetName()<<"->SetMarkerStyle("<<GetMarkerStyle()<<");"<<std::endl;
   }
   if (GetMarkerSize() != 1) {
      out<<"   "<<GetName()<<"->SetMarkerSize("<<GetMarkerSize()<<");"<<std::endl;
   }
   if (GetLineColor() != 1) {
      if (GetLineColor() > 228) {
         TColor::SaveColor(out, GetLineColor());
         out<<"   "<<GetName()<<"->SetLineColor(ci);" << std::endl;
      } else 
         out<<"   "<<GetName()<<"->SetLineColor("<<GetLineColor()<<");"<<std::endl;
   }
   if (GetLineWidth() != 4) {
      out<<"   "<<GetName()<<"->SetLineWidth("<<GetLineWidth()<<");"<<std::endl;
   }
   if (GetLineStyle() != 1) {
      out<<"   "<<GetName()<<"->SetLineStyle("<<GetLineStyle()<<");"<<std::endl;
   }
   if (GetNpx() != 100) {
      out<<"   "<<GetName()<<"->SetNpx("<<GetNpx()<<");"<<std::endl;
   }
   if (GetChisquare() != 0) {
      out<<"   "<<GetName()<<"->SetChisquare("<<GetChisquare()<<");"<<std::endl;
   }
   Double_t parmin, parmax;
   for (Int_t i=0;i<GetNpar();i++) {
      out<<"   "<<GetName()<<"->SetParameter("<<i<<","<<GetParameter(i)<<");"<<std::endl;
      out<<"   "<<GetName()<<"->SetParError("<<i<<","<<GetParError(i)<<");"<<std::endl;
      GetParLimits(i,parmin,parmax);
      out<<"   "<<GetName()<<"->SetParLimits("<<i<<","<<parmin<<","<<parmax<<");"<<std::endl;
   }
   out<<"   "<<GetName()<<"->Draw("
      <<quote<<option<<quote<<");"<<std::endl;
}



//______________________________________________________________________________
void TF2::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
   //

   Int_t level;
   if (nlevels <=0 ) {
      fContour.Set(0);
      return;
   }
   fContour.Set(nlevels);

   //*-*-  Contour levels are specified
   if (levels) {
      for (level=0; level<nlevels; level++) fContour.fArray[level] = levels[level];
   } else {
      fContour.fArray[0] = -9999; // means not defined at this point
   }
}


//______________________________________________________________________________
void TF2::SetContourLevel(Int_t level, Double_t value)
{
   //*-*-*-*-*-*-*-*-*-*-*Set value for one contour level*-*-*-*-*-*-*-*-*-*-*-*
   //*-*                  ===============================
   if (level <0 || level >= fContour.fN) return;
   fContour.fArray[level] = value;
}

//______________________________________________________________________________
void TF2::SetNpy(Int_t npy)
{
   // Set the number of points used to draw the function
   //
   // The default number of points along x is 30 for 2-d/3-d functions.
   // You can increase this value to get a better resolution when drawing
   // pictures with sharp peaks or to get a better result when using TF2::GetRandom2   
   // the minimum number of points is 4, the maximum is 10000 for 2-d/3-d functions

   if (npy < 4) {
      Warning("SetNpy","Number of points must be >=4 && <= 10000, fNpy set to 4");
      fNpy = 4;
   } else if(npy > 10000) {
      Warning("SetNpy","Number of points must be >=4 && <= 10000, fNpy set to 10000");
      fNpy = 10000;
   } else {
      fNpy = npy;
   }
   Update();
}

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

   fXmin = xmin;
   fXmax = xmax;
   fYmin = ymin;
   fYmax = ymax;
   Update();
}

//______________________________________________________________________________
void TF2::Streamer(TBuffer &R__b)
{
   // Stream an object of class TF2.

   if (R__b.IsReading()) {
      UInt_t R__s, R__c;
      Version_t R__v = R__b.ReadVersion(&R__s, &R__c);
      if (R__v > 3) {
         R__b.ReadClassBuffer(TF2::Class(), this, R__v, R__s, R__c);
         return;
      }
      //====process old versions before automatic schema evolution
      Int_t nlevels;
      TF1::Streamer(R__b);
      if (R__v < 3) {
         Float_t ymin,ymax;
         R__b >> ymin; fYmin = ymin;
         R__b >> ymax; fYmax = ymax;
      } else {
         R__b >> fYmin;
         R__b >> fYmax;
      }
      R__b >> fNpy;
      R__b >> nlevels;
      if (R__v < 3) {
         Float_t *contour = 0;
         Int_t n = R__b.ReadArray(contour);
         fContour.Set(n);
         for (Int_t i=0;i<n;i++) fContour.fArray[i] = contour[i];
         delete [] contour;
      } else {
         fContour.Streamer(R__b);
      }
      R__b.CheckByteCount(R__s, R__c, TF2::IsA());
      //====end of old versions

   } else {
      Int_t saved = 0;
      if (fType > 0 && fSave.empty()) { saved = 1; Save(fXmin,fXmax,fYmin,fYmax,0,0);}

      R__b.WriteClassBuffer(TF2::Class(),this);

      if (saved) {fSave.clear(); }
   }
}

//______________________________________________________________________________
Double_t TF2::Moment2(Double_t nx, Double_t ax, Double_t bx, Double_t ny, Double_t ay, Double_t by, Double_t epsilon)
{
// Return x^nx * y^ny moment of a 2d function in range [ax,bx],[ay,by]
//   Author: Gene Van Buren <gene@bnl.gov>

   Double_t norm = Integral(ax,bx,ay,by,epsilon);
   if (norm == 0) {
      Error("Moment2", "Integral zero over range");
      return 0;
   }

   TF2 fnc("TF2_ExpValHelper",Form("%s*pow(x,%f)*pow(y,%f)",GetName(),nx,ny));
   return fnc.Integral(ax,bx,ay,by,epsilon)/norm;
}

//______________________________________________________________________________
Double_t TF2::CentralMoment2(Double_t nx, Double_t ax, Double_t bx, Double_t ny, Double_t ay, Double_t by, Double_t epsilon)
{
// Return x^nx * y^ny central moment of a 2d function in range [ax,bx],[ay,by]
//   Author: Gene Van Buren <gene@bnl.gov>

   Double_t norm = Integral(ax,bx,ay,by,epsilon);
   if (norm == 0) {
      Error("CentralMoment2", "Integral zero over range");
      return 0;
   }

   Double_t xbar = 0;
   Double_t ybar = 0;
   if (nx!=0) {
      TF2 fncx("TF2_ExpValHelperx",Form("%s*x",GetName()));
      xbar = fncx.Integral(ax,bx,ay,by,epsilon)/norm;
   }
   if (ny!=0) {
      TF2 fncx("TF2_ExpValHelpery",Form("%s*y",GetName()));
      ybar = fncx.Integral(ax,bx,ay,by,epsilon)/norm;
   }
   TF2 fnc("TF2_ExpValHelper",Form("%s*pow(x-%f,%f)*pow(y-%f,%f)",GetName(),xbar,nx,ybar,ny));
   return fnc.Integral(ax,bx,ay,by,epsilon)/norm;
}

 TF2.cxx:1
 TF2.cxx:2
 TF2.cxx:3
 TF2.cxx:4
 TF2.cxx:5
 TF2.cxx:6
 TF2.cxx:7
 TF2.cxx:8
 TF2.cxx:9
 TF2.cxx:10
 TF2.cxx:11
 TF2.cxx:12
 TF2.cxx:13
 TF2.cxx:14
 TF2.cxx:15
 TF2.cxx:16
 TF2.cxx:17
 TF2.cxx:18
 TF2.cxx:19
 TF2.cxx:20
 TF2.cxx:21
 TF2.cxx:22
 TF2.cxx:23
 TF2.cxx:24
 TF2.cxx:25
 TF2.cxx:26
 TF2.cxx:27
 TF2.cxx:28
 TF2.cxx:29
 TF2.cxx:30
 TF2.cxx:31
 TF2.cxx:32
 TF2.cxx:33
 TF2.cxx:34
 TF2.cxx:35
 TF2.cxx:36
 TF2.cxx:37
 TF2.cxx:38
 TF2.cxx:39
 TF2.cxx:40
 TF2.cxx:41
 TF2.cxx:42
 TF2.cxx:43
 TF2.cxx:44
 TF2.cxx:45
 TF2.cxx:46
 TF2.cxx:47
 TF2.cxx:48
 TF2.cxx:49
 TF2.cxx:50
 TF2.cxx:51
 TF2.cxx:52
 TF2.cxx:53
 TF2.cxx:54
 TF2.cxx:55
 TF2.cxx:56
 TF2.cxx:57
 TF2.cxx:58
 TF2.cxx:59
 TF2.cxx:60
 TF2.cxx:61
 TF2.cxx:62
 TF2.cxx:63
 TF2.cxx:64
 TF2.cxx:65
 TF2.cxx:66
 TF2.cxx:67
 TF2.cxx:68
 TF2.cxx:69
 TF2.cxx:70
 TF2.cxx:71
 TF2.cxx:72
 TF2.cxx:73
 TF2.cxx:74
 TF2.cxx:75
 TF2.cxx:76
 TF2.cxx:77
 TF2.cxx:78
 TF2.cxx:79
 TF2.cxx:80
 TF2.cxx:81
 TF2.cxx:82
 TF2.cxx:83
 TF2.cxx:84
 TF2.cxx:85
 TF2.cxx:86
 TF2.cxx:87
 TF2.cxx:88
 TF2.cxx:89
 TF2.cxx:90
 TF2.cxx:91
 TF2.cxx:92
 TF2.cxx:93
 TF2.cxx:94
 TF2.cxx:95
 TF2.cxx:96
 TF2.cxx:97
 TF2.cxx:98
 TF2.cxx:99
 TF2.cxx:100
 TF2.cxx:101
 TF2.cxx:102
 TF2.cxx:103
 TF2.cxx:104
 TF2.cxx:105
 TF2.cxx:106
 TF2.cxx:107
 TF2.cxx:108
 TF2.cxx:109
 TF2.cxx:110
 TF2.cxx:111
 TF2.cxx:112
 TF2.cxx:113
 TF2.cxx:114
 TF2.cxx:115
 TF2.cxx:116
 TF2.cxx:117
 TF2.cxx:118
 TF2.cxx:119
 TF2.cxx:120
 TF2.cxx:121
 TF2.cxx:122
 TF2.cxx:123
 TF2.cxx:124
 TF2.cxx:125
 TF2.cxx:126
 TF2.cxx:127
 TF2.cxx:128
 TF2.cxx:129
 TF2.cxx:130
 TF2.cxx:131
 TF2.cxx:132
 TF2.cxx:133
 TF2.cxx:134
 TF2.cxx:135
 TF2.cxx:136
 TF2.cxx:137
 TF2.cxx:138
 TF2.cxx:139
 TF2.cxx:140
 TF2.cxx:141
 TF2.cxx:142
 TF2.cxx:143
 TF2.cxx:144
 TF2.cxx:145
 TF2.cxx:146
 TF2.cxx:147
 TF2.cxx:148
 TF2.cxx:149
 TF2.cxx:150
 TF2.cxx:151
 TF2.cxx:152
 TF2.cxx:153
 TF2.cxx:154
 TF2.cxx:155
 TF2.cxx:156
 TF2.cxx:157
 TF2.cxx:158
 TF2.cxx:159
 TF2.cxx:160
 TF2.cxx:161
 TF2.cxx:162
 TF2.cxx:163
 TF2.cxx:164
 TF2.cxx:165
 TF2.cxx:166
 TF2.cxx:167
 TF2.cxx:168
 TF2.cxx:169
 TF2.cxx:170
 TF2.cxx:171
 TF2.cxx:172
 TF2.cxx:173
 TF2.cxx:174
 TF2.cxx:175
 TF2.cxx:176
 TF2.cxx:177
 TF2.cxx:178
 TF2.cxx:179
 TF2.cxx:180
 TF2.cxx:181
 TF2.cxx:182
 TF2.cxx:183
 TF2.cxx:184
 TF2.cxx:185
 TF2.cxx:186
 TF2.cxx:187
 TF2.cxx:188
 TF2.cxx:189
 TF2.cxx:190
 TF2.cxx:191
 TF2.cxx:192
 TF2.cxx:193
 TF2.cxx:194
 TF2.cxx:195
 TF2.cxx:196
 TF2.cxx:197
 TF2.cxx:198
 TF2.cxx:199
 TF2.cxx:200
 TF2.cxx:201
 TF2.cxx:202
 TF2.cxx:203
 TF2.cxx:204
 TF2.cxx:205
 TF2.cxx:206
 TF2.cxx:207
 TF2.cxx:208
 TF2.cxx:209
 TF2.cxx:210
 TF2.cxx:211
 TF2.cxx:212
 TF2.cxx:213
 TF2.cxx:214
 TF2.cxx:215
 TF2.cxx:216
 TF2.cxx:217
 TF2.cxx:218
 TF2.cxx:219
 TF2.cxx:220
 TF2.cxx:221
 TF2.cxx:222
 TF2.cxx:223
 TF2.cxx:224
 TF2.cxx:225
 TF2.cxx:226
 TF2.cxx:227
 TF2.cxx:228
 TF2.cxx:229
 TF2.cxx:230
 TF2.cxx:231
 TF2.cxx:232
 TF2.cxx:233
 TF2.cxx:234
 TF2.cxx:235
 TF2.cxx:236
 TF2.cxx:237
 TF2.cxx:238
 TF2.cxx:239
 TF2.cxx:240
 TF2.cxx:241
 TF2.cxx:242
 TF2.cxx:243
 TF2.cxx:244
 TF2.cxx:245
 TF2.cxx:246
 TF2.cxx:247
 TF2.cxx:248
 TF2.cxx:249
 TF2.cxx:250
 TF2.cxx:251
 TF2.cxx:252
 TF2.cxx:253
 TF2.cxx:254
 TF2.cxx:255
 TF2.cxx:256
 TF2.cxx:257
 TF2.cxx:258
 TF2.cxx:259
 TF2.cxx:260
 TF2.cxx:261
 TF2.cxx:262
 TF2.cxx:263
 TF2.cxx:264
 TF2.cxx:265
 TF2.cxx:266
 TF2.cxx:267
 TF2.cxx:268
 TF2.cxx:269
 TF2.cxx:270
 TF2.cxx:271
 TF2.cxx:272
 TF2.cxx:273
 TF2.cxx:274
 TF2.cxx:275
 TF2.cxx:276
 TF2.cxx:277
 TF2.cxx:278
 TF2.cxx:279
 TF2.cxx:280
 TF2.cxx:281
 TF2.cxx:282
 TF2.cxx:283
 TF2.cxx:284
 TF2.cxx:285
 TF2.cxx:286
 TF2.cxx:287
 TF2.cxx:288
 TF2.cxx:289
 TF2.cxx:290
 TF2.cxx:291
 TF2.cxx:292
 TF2.cxx:293
 TF2.cxx:294
 TF2.cxx:295
 TF2.cxx:296
 TF2.cxx:297
 TF2.cxx:298
 TF2.cxx:299
 TF2.cxx:300
 TF2.cxx:301
 TF2.cxx:302
 TF2.cxx:303
 TF2.cxx:304
 TF2.cxx:305
 TF2.cxx:306
 TF2.cxx:307
 TF2.cxx:308
 TF2.cxx:309
 TF2.cxx:310
 TF2.cxx:311
 TF2.cxx:312
 TF2.cxx:313
 TF2.cxx:314
 TF2.cxx:315
 TF2.cxx:316
 TF2.cxx:317
 TF2.cxx:318
 TF2.cxx:319
 TF2.cxx:320
 TF2.cxx:321
 TF2.cxx:322
 TF2.cxx:323
 TF2.cxx:324
 TF2.cxx:325
 TF2.cxx:326
 TF2.cxx:327
 TF2.cxx:328
 TF2.cxx:329
 TF2.cxx:330
 TF2.cxx:331
 TF2.cxx:332
 TF2.cxx:333
 TF2.cxx:334
 TF2.cxx:335
 TF2.cxx:336
 TF2.cxx:337
 TF2.cxx:338
 TF2.cxx:339
 TF2.cxx:340
 TF2.cxx:341
 TF2.cxx:342
 TF2.cxx:343
 TF2.cxx:344
 TF2.cxx:345
 TF2.cxx:346
 TF2.cxx:347
 TF2.cxx:348
 TF2.cxx:349
 TF2.cxx:350
 TF2.cxx:351
 TF2.cxx:352
 TF2.cxx:353
 TF2.cxx:354
 TF2.cxx:355
 TF2.cxx:356
 TF2.cxx:357
 TF2.cxx:358
 TF2.cxx:359
 TF2.cxx:360
 TF2.cxx:361
 TF2.cxx:362
 TF2.cxx:363
 TF2.cxx:364
 TF2.cxx:365
 TF2.cxx:366
 TF2.cxx:367
 TF2.cxx:368
 TF2.cxx:369
 TF2.cxx:370
 TF2.cxx:371
 TF2.cxx:372
 TF2.cxx:373
 TF2.cxx:374
 TF2.cxx:375
 TF2.cxx:376
 TF2.cxx:377
 TF2.cxx:378
 TF2.cxx:379
 TF2.cxx:380
 TF2.cxx:381
 TF2.cxx:382
 TF2.cxx:383
 TF2.cxx:384
 TF2.cxx:385
 TF2.cxx:386
 TF2.cxx:387
 TF2.cxx:388
 TF2.cxx:389
 TF2.cxx:390
 TF2.cxx:391
 TF2.cxx:392
 TF2.cxx:393
 TF2.cxx:394
 TF2.cxx:395
 TF2.cxx:396
 TF2.cxx:397
 TF2.cxx:398
 TF2.cxx:399
 TF2.cxx:400
 TF2.cxx:401
 TF2.cxx:402
 TF2.cxx:403
 TF2.cxx:404
 TF2.cxx:405
 TF2.cxx:406
 TF2.cxx:407
 TF2.cxx:408
 TF2.cxx:409
 TF2.cxx:410
 TF2.cxx:411
 TF2.cxx:412
 TF2.cxx:413
 TF2.cxx:414
 TF2.cxx:415
 TF2.cxx:416
 TF2.cxx:417
 TF2.cxx:418
 TF2.cxx:419
 TF2.cxx:420
 TF2.cxx:421
 TF2.cxx:422
 TF2.cxx:423
 TF2.cxx:424
 TF2.cxx:425
 TF2.cxx:426
 TF2.cxx:427
 TF2.cxx:428
 TF2.cxx:429
 TF2.cxx:430
 TF2.cxx:431
 TF2.cxx:432
 TF2.cxx:433
 TF2.cxx:434
 TF2.cxx:435
 TF2.cxx:436
 TF2.cxx:437
 TF2.cxx:438
 TF2.cxx:439
 TF2.cxx:440
 TF2.cxx:441
 TF2.cxx:442
 TF2.cxx:443
 TF2.cxx:444
 TF2.cxx:445
 TF2.cxx:446
 TF2.cxx:447
 TF2.cxx:448
 TF2.cxx:449
 TF2.cxx:450
 TF2.cxx:451
 TF2.cxx:452
 TF2.cxx:453
 TF2.cxx:454
 TF2.cxx:455
 TF2.cxx:456
 TF2.cxx:457
 TF2.cxx:458
 TF2.cxx:459
 TF2.cxx:460
 TF2.cxx:461
 TF2.cxx:462
 TF2.cxx:463
 TF2.cxx:464
 TF2.cxx:465
 TF2.cxx:466
 TF2.cxx:467
 TF2.cxx:468
 TF2.cxx:469
 TF2.cxx:470
 TF2.cxx:471
 TF2.cxx:472
 TF2.cxx:473
 TF2.cxx:474
 TF2.cxx:475
 TF2.cxx:476
 TF2.cxx:477
 TF2.cxx:478
 TF2.cxx:479
 TF2.cxx:480
 TF2.cxx:481
 TF2.cxx:482
 TF2.cxx:483
 TF2.cxx:484
 TF2.cxx:485
 TF2.cxx:486
 TF2.cxx:487
 TF2.cxx:488
 TF2.cxx:489
 TF2.cxx:490
 TF2.cxx:491
 TF2.cxx:492
 TF2.cxx:493
 TF2.cxx:494
 TF2.cxx:495
 TF2.cxx:496
 TF2.cxx:497
 TF2.cxx:498
 TF2.cxx:499
 TF2.cxx:500
 TF2.cxx:501
 TF2.cxx:502
 TF2.cxx:503
 TF2.cxx:504
 TF2.cxx:505
 TF2.cxx:506
 TF2.cxx:507
 TF2.cxx:508
 TF2.cxx:509
 TF2.cxx:510
 TF2.cxx:511
 TF2.cxx:512
 TF2.cxx:513
 TF2.cxx:514
 TF2.cxx:515
 TF2.cxx:516
 TF2.cxx:517
 TF2.cxx:518
 TF2.cxx:519
 TF2.cxx:520
 TF2.cxx:521
 TF2.cxx:522
 TF2.cxx:523
 TF2.cxx:524
 TF2.cxx:525
 TF2.cxx:526
 TF2.cxx:527
 TF2.cxx:528
 TF2.cxx:529
 TF2.cxx:530
 TF2.cxx:531
 TF2.cxx:532
 TF2.cxx:533
 TF2.cxx:534
 TF2.cxx:535
 TF2.cxx:536
 TF2.cxx:537
 TF2.cxx:538
 TF2.cxx:539
 TF2.cxx:540
 TF2.cxx:541
 TF2.cxx:542
 TF2.cxx:543
 TF2.cxx:544
 TF2.cxx:545
 TF2.cxx:546
 TF2.cxx:547
 TF2.cxx:548
 TF2.cxx:549
 TF2.cxx:550
 TF2.cxx:551
 TF2.cxx:552
 TF2.cxx:553
 TF2.cxx:554
 TF2.cxx:555
 TF2.cxx:556
 TF2.cxx:557
 TF2.cxx:558
 TF2.cxx:559
 TF2.cxx:560
 TF2.cxx:561
 TF2.cxx:562
 TF2.cxx:563
 TF2.cxx:564
 TF2.cxx:565
 TF2.cxx:566
 TF2.cxx:567
 TF2.cxx:568
 TF2.cxx:569
 TF2.cxx:570
 TF2.cxx:571
 TF2.cxx:572
 TF2.cxx:573
 TF2.cxx:574
 TF2.cxx:575
 TF2.cxx:576
 TF2.cxx:577
 TF2.cxx:578
 TF2.cxx:579
 TF2.cxx:580
 TF2.cxx:581
 TF2.cxx:582
 TF2.cxx:583
 TF2.cxx:584
 TF2.cxx:585
 TF2.cxx:586
 TF2.cxx:587
 TF2.cxx:588
 TF2.cxx:589
 TF2.cxx:590
 TF2.cxx:591
 TF2.cxx:592
 TF2.cxx:593
 TF2.cxx:594
 TF2.cxx:595
 TF2.cxx:596
 TF2.cxx:597
 TF2.cxx:598
 TF2.cxx:599
 TF2.cxx:600
 TF2.cxx:601
 TF2.cxx:602
 TF2.cxx:603
 TF2.cxx:604
 TF2.cxx:605
 TF2.cxx:606
 TF2.cxx:607
 TF2.cxx:608
 TF2.cxx:609
 TF2.cxx:610
 TF2.cxx:611
 TF2.cxx:612
 TF2.cxx:613
 TF2.cxx:614
 TF2.cxx:615
 TF2.cxx:616
 TF2.cxx:617
 TF2.cxx:618
 TF2.cxx:619
 TF2.cxx:620
 TF2.cxx:621
 TF2.cxx:622
 TF2.cxx:623
 TF2.cxx:624
 TF2.cxx:625
 TF2.cxx:626
 TF2.cxx:627
 TF2.cxx:628
 TF2.cxx:629
 TF2.cxx:630
 TF2.cxx:631
 TF2.cxx:632
 TF2.cxx:633
 TF2.cxx:634
 TF2.cxx:635
 TF2.cxx:636
 TF2.cxx:637
 TF2.cxx:638
 TF2.cxx:639
 TF2.cxx:640
 TF2.cxx:641
 TF2.cxx:642
 TF2.cxx:643
 TF2.cxx:644
 TF2.cxx:645
 TF2.cxx:646
 TF2.cxx:647
 TF2.cxx:648
 TF2.cxx:649
 TF2.cxx:650
 TF2.cxx:651
 TF2.cxx:652
 TF2.cxx:653
 TF2.cxx:654
 TF2.cxx:655
 TF2.cxx:656
 TF2.cxx:657
 TF2.cxx:658
 TF2.cxx:659
 TF2.cxx:660
 TF2.cxx:661
 TF2.cxx:662
 TF2.cxx:663
 TF2.cxx:664
 TF2.cxx:665
 TF2.cxx:666
 TF2.cxx:667
 TF2.cxx:668
 TF2.cxx:669
 TF2.cxx:670
 TF2.cxx:671
 TF2.cxx:672
 TF2.cxx:673
 TF2.cxx:674
 TF2.cxx:675
 TF2.cxx:676
 TF2.cxx:677
 TF2.cxx:678
 TF2.cxx:679
 TF2.cxx:680
 TF2.cxx:681
 TF2.cxx:682
 TF2.cxx:683
 TF2.cxx:684
 TF2.cxx:685
 TF2.cxx:686
 TF2.cxx:687
 TF2.cxx:688
 TF2.cxx:689
 TF2.cxx:690
 TF2.cxx:691
 TF2.cxx:692
 TF2.cxx:693
 TF2.cxx:694
 TF2.cxx:695
 TF2.cxx:696
 TF2.cxx:697
 TF2.cxx:698
 TF2.cxx:699
 TF2.cxx:700
 TF2.cxx:701
 TF2.cxx:702
 TF2.cxx:703
 TF2.cxx:704
 TF2.cxx:705
 TF2.cxx:706
 TF2.cxx:707
 TF2.cxx:708
 TF2.cxx:709
 TF2.cxx:710
 TF2.cxx:711
 TF2.cxx:712
 TF2.cxx:713
 TF2.cxx:714
 TF2.cxx:715
 TF2.cxx:716
 TF2.cxx:717
 TF2.cxx:718
 TF2.cxx:719
 TF2.cxx:720
 TF2.cxx:721
 TF2.cxx:722
 TF2.cxx:723
 TF2.cxx:724
 TF2.cxx:725
 TF2.cxx:726
 TF2.cxx:727
 TF2.cxx:728
 TF2.cxx:729
 TF2.cxx:730
 TF2.cxx:731
 TF2.cxx:732
 TF2.cxx:733
 TF2.cxx:734
 TF2.cxx:735
 TF2.cxx:736
 TF2.cxx:737
 TF2.cxx:738
 TF2.cxx:739
 TF2.cxx:740
 TF2.cxx:741
 TF2.cxx:742
 TF2.cxx:743
 TF2.cxx:744
 TF2.cxx:745
 TF2.cxx:746
 TF2.cxx:747
 TF2.cxx:748
 TF2.cxx:749
 TF2.cxx:750
 TF2.cxx:751
 TF2.cxx:752
 TF2.cxx:753
 TF2.cxx:754
 TF2.cxx:755
 TF2.cxx:756
 TF2.cxx:757
 TF2.cxx:758
 TF2.cxx:759
 TF2.cxx:760
 TF2.cxx:761
 TF2.cxx:762
 TF2.cxx:763
 TF2.cxx:764
 TF2.cxx:765
 TF2.cxx:766
 TF2.cxx:767
 TF2.cxx:768
 TF2.cxx:769
 TF2.cxx:770
 TF2.cxx:771
 TF2.cxx:772
 TF2.cxx:773
 TF2.cxx:774
 TF2.cxx:775
 TF2.cxx:776
 TF2.cxx:777
 TF2.cxx:778
 TF2.cxx:779
 TF2.cxx:780
 TF2.cxx:781
 TF2.cxx:782
 TF2.cxx:783
 TF2.cxx:784
 TF2.cxx:785
 TF2.cxx:786
 TF2.cxx:787
 TF2.cxx:788
 TF2.cxx:789
 TF2.cxx:790
 TF2.cxx:791
 TF2.cxx:792
 TF2.cxx:793
 TF2.cxx:794
 TF2.cxx:795
 TF2.cxx:796
 TF2.cxx:797
 TF2.cxx:798
 TF2.cxx:799
 TF2.cxx:800
 TF2.cxx:801
 TF2.cxx:802
 TF2.cxx:803
 TF2.cxx:804
 TF2.cxx:805
 TF2.cxx:806
 TF2.cxx:807
 TF2.cxx:808
 TF2.cxx:809
 TF2.cxx:810
 TF2.cxx:811
 TF2.cxx:812
 TF2.cxx:813
 TF2.cxx:814
 TF2.cxx:815
 TF2.cxx:816
 TF2.cxx:817
 TF2.cxx:818
 TF2.cxx:819
 TF2.cxx:820
 TF2.cxx:821
 TF2.cxx:822
 TF2.cxx:823
 TF2.cxx:824
 TF2.cxx:825
 TF2.cxx:826
 TF2.cxx:827
 TF2.cxx:828
 TF2.cxx:829
 TF2.cxx:830
 TF2.cxx:831
 TF2.cxx:832
 TF2.cxx:833
 TF2.cxx:834
 TF2.cxx:835
 TF2.cxx:836
 TF2.cxx:837
 TF2.cxx:838
 TF2.cxx:839
 TF2.cxx:840
 TF2.cxx:841
 TF2.cxx:842
 TF2.cxx:843
 TF2.cxx:844
 TF2.cxx:845
 TF2.cxx:846
 TF2.cxx:847
 TF2.cxx:848
 TF2.cxx:849
 TF2.cxx:850
 TF2.cxx:851
 TF2.cxx:852
 TF2.cxx:853
 TF2.cxx:854
 TF2.cxx:855
 TF2.cxx:856
 TF2.cxx:857
 TF2.cxx:858
 TF2.cxx:859
 TF2.cxx:860
 TF2.cxx:861
 TF2.cxx:862
 TF2.cxx:863
 TF2.cxx:864
 TF2.cxx:865
 TF2.cxx:866
 TF2.cxx:867
 TF2.cxx:868
 TF2.cxx:869
 TF2.cxx:870
 TF2.cxx:871
 TF2.cxx:872
 TF2.cxx:873
 TF2.cxx:874
 TF2.cxx:875
 TF2.cxx:876
 TF2.cxx:877
 TF2.cxx:878
 TF2.cxx:879
 TF2.cxx:880
 TF2.cxx:881
 TF2.cxx:882
 TF2.cxx:883
 TF2.cxx:884
 TF2.cxx:885
 TF2.cxx:886
 TF2.cxx:887
 TF2.cxx:888
 TF2.cxx:889
 TF2.cxx:890
 TF2.cxx:891
 TF2.cxx:892
 TF2.cxx:893
 TF2.cxx:894
 TF2.cxx:895
 TF2.cxx:896
 TF2.cxx:897
 TF2.cxx:898
 TF2.cxx:899
 TF2.cxx:900
 TF2.cxx:901
 TF2.cxx:902
 TF2.cxx:903
 TF2.cxx:904
 TF2.cxx:905
 TF2.cxx:906
 TF2.cxx:907
 TF2.cxx:908
 TF2.cxx:909
 TF2.cxx:910
 TF2.cxx:911
 TF2.cxx:912
 TF2.cxx:913
 TF2.cxx:914
 TF2.cxx:915
 TF2.cxx:916
 TF2.cxx:917
 TF2.cxx:918
 TF2.cxx:919
 TF2.cxx:920
 TF2.cxx:921
 TF2.cxx:922
 TF2.cxx:923
 TF2.cxx:924
 TF2.cxx:925
 TF2.cxx:926
 TF2.cxx:927
 TF2.cxx:928
 TF2.cxx:929
 TF2.cxx:930
 TF2.cxx:931
 TF2.cxx:932
 TF2.cxx:933
 TF2.cxx:934
 TF2.cxx:935
 TF2.cxx:936
 TF2.cxx:937
 TF2.cxx:938
 TF2.cxx:939
 TF2.cxx:940
 TF2.cxx:941
 TF2.cxx:942
 TF2.cxx:943
 TF2.cxx:944
 TF2.cxx:945
 TF2.cxx:946
 TF2.cxx:947
 TF2.cxx:948
 TF2.cxx:949
 TF2.cxx:950
 TF2.cxx:951
 TF2.cxx:952
 TF2.cxx:953
 TF2.cxx:954
 TF2.cxx:955
 TF2.cxx:956
 TF2.cxx:957
 TF2.cxx:958
 TF2.cxx:959
 TF2.cxx:960
 TF2.cxx:961
 TF2.cxx:962
 TF2.cxx:963
 TF2.cxx:964
 TF2.cxx:965
 TF2.cxx:966
 TF2.cxx:967
 TF2.cxx:968
 TF2.cxx:969
 TF2.cxx:970
 TF2.cxx:971
 TF2.cxx:972
 TF2.cxx:973
 TF2.cxx:974
 TF2.cxx:975
 TF2.cxx:976
 TF2.cxx:977
 TF2.cxx:978
 TF2.cxx:979
 TF2.cxx:980
 TF2.cxx:981
 TF2.cxx:982
 TF2.cxx:983
 TF2.cxx:984
 TF2.cxx:985
 TF2.cxx:986
 TF2.cxx:987
 TF2.cxx:988
 TF2.cxx:989
 TF2.cxx:990
 TF2.cxx:991
 TF2.cxx:992
 TF2.cxx:993
 TF2.cxx:994
 TF2.cxx:995
 TF2.cxx:996
 TF2.cxx:997
 TF2.cxx:998
 TF2.cxx:999
 TF2.cxx:1000
 TF2.cxx:1001
 TF2.cxx:1002
 TF2.cxx:1003
 TF2.cxx:1004
 TF2.cxx:1005
 TF2.cxx:1006
 TF2.cxx:1007
 TF2.cxx:1008
 TF2.cxx:1009
 TF2.cxx:1010
 TF2.cxx:1011
 TF2.cxx:1012
 TF2.cxx:1013
 TF2.cxx:1014
 TF2.cxx:1015
 TF2.cxx:1016
 TF2.cxx:1017
 TF2.cxx:1018
 TF2.cxx:1019
 TF2.cxx:1020
 TF2.cxx:1021
 TF2.cxx:1022
 TF2.cxx:1023
 TF2.cxx:1024
 TF2.cxx:1025
 TF2.cxx:1026
 TF2.cxx:1027
 TF2.cxx:1028
 TF2.cxx:1029
 TF2.cxx:1030
 TF2.cxx:1031
 TF2.cxx:1032
 TF2.cxx:1033
 TF2.cxx:1034
 TF2.cxx:1035
 TF2.cxx:1036
 TF2.cxx:1037
 TF2.cxx:1038