ROOT logo
// @(#)root/foam:$Id: TFoam.cxx 29138 2009-06-22 14:31:00Z brun $
// Author: S. Jadach <mailto:Stanislaw.jadach@ifj.edu.pl>, P.Sawicki <mailto:Pawel.Sawicki@ifj.edu.pl>

//______________________________________________________________________________
//
// FOAM  Version 1.02M
// ===================
// Authors:
//   S. Jadach and P.Sawicki
//   Institute of Nuclear Physics, Cracow, Poland
//   Stanislaw. Jadach@ifj.edu.pl, Pawel.Sawicki@ifj.edu.pl
//
// What is FOAM for?
// =================
// * Suppose you want to generate randomly points (vectors) according to
//   an arbitrary probability distribution  in n dimensions,
//   for which you supply your own subprogram. FOAM can do it for you!
//   Even if your distributions has quite strong peaks and is discontinuous!
// * FOAM generates random points with weight one or with variable weight.
// * FOAM is capable to integrate using efficient "adaptive" MC method.
//   (The distribution does not need to be normalized to one.)
// How does it work?
// =================
// FOAM is the simplified version of the multi-dimensional general purpose
// Monte Carlo event generator (integrator) FOAM.
// It creates hyper-rectangular "foam of cells", which is more dense around its peaks.
// See the following 2-dim. example of the map of 1000 cells for doubly peaked distribution:
//BEGIN_HTML <!--
/* -->
<img src="gif/foam_MapCamel1000.gif">
<!--*/
// -->END_HTML
// FOAM is now fully integrated with the ROOT package.
// The important bonus of the ROOT use is persistency of the FOAM objects!
//
// For more sophisticated problems full version of FOAM may be more appropriate:
//BEGIN_HTML <!--
/* -->
  See <A HREF="http://jadach.home.cern.ch/jadach/Foam/Index.html"> full version of FOAM</A>
<!--*/
// -->END_HTML
// Simple example of the use of FOAM:
// ==================================
// Int_t kanwa(){
//   gSystem->Load("libFoam");
//   TH2D  *hst_xy = new TH2D("hst_xy" ,  "x-y plot", 50,0,1.0, 50,0,1.0);
//   Double_t *MCvect =new Double_t[2]; // 2-dim vector generated in the MC run
//   TRandom3  *PseRan   = new TRandom3();  // Create random number generator
//   PseRan->SetSeed(4357);                // Set seed
//   TFoam   *FoamX    = new TFoam("FoamX");   // Create Simulator
//   FoamX->SetkDim(2);          // No. of dimensions, obligatory!
//   FoamX->SetnCells(500);      // No. of cells, can be omitted, default=2000
//   FoamX->SetRhoInt(Camel2);   // Set 2-dim distribution, included below
//   FoamX->SetPseRan(PseRan);   // Set random number generator
//   FoamX->Initialize();        // Initialize simulator, takes a few seconds...
//   // From now on FoamX is ready to generate events according to Camel2(x,y)
//   for(Long_t loop=0; loop<100000; loop++){
//     FoamX->MakeEvent();          // generate MC event
//     FoamX->GetMCvect( MCvect);   // get generated vector (x,y)
//     Double_t x=MCvect[0];
//     Double_t y=MCvect[1];
//     if(loop<10) cout<<"(x,y) =  ( "<< x <<", "<< y <<" )"<<endl;
//     hst_xy->Fill(x,y);           // fill scattergram
//   }// loop
//   Double_t mcResult, mcError;
//   FoamX->GetIntegMC( mcResult, mcError);  // get MC integral, should be one
//   cout << " mcResult= " << mcResult << " +- " << mcError <<endl;
//   // now hst_xy will be plotted visualizing generated distribution
//   TCanvas *cKanwa = new TCanvas("cKanwa","Canvas for plotting",600,600);
//   cKanwa->cd();
//   hst_xy->Draw("lego2");
// }//kanwa
// Double_t sqr(Double_t x){return x*x;};
// Double_t Camel2(Int_t nDim, Double_t *Xarg){
// // 2-dimensional distribution for FOAM, normalized to one (within 1e-5)
//   Double_t x=Xarg[0];
//   Double_t y=Xarg[1];
//   Double_t GamSq= sqr(0.100e0);
//   Double_t Dist=exp(-(sqr(x-1./3) +sqr(y-1./3))/GamSq)/GamSq/TMath::Pi();
//   Dist        +=exp(-(sqr(x-2./3) +sqr(y-2./3))/GamSq)/GamSq/TMath::Pi();
//   return 0.5*Dist;
// }// Camel2
// Two-dim. histogram of the MC points generated with the above program looks as follows:
//BEGIN_HTML <!--
/* -->
<img src="gif/foam_cKanwa.gif">
<!--*/
// -->END_HTML
// Canonical nine steering parameters of FOAM
// ===========================================
//------------------------------------------------------------------------------
//  Name     | default  | Description
//------------------------------------------------------------------------------
//  kDim     | 0        | Dimension of the integration space. Must be redefined!
//  nCells   | 1000     | No of allocated number of cells,
//  nSampl   | 200      | No. of MC events in the cell MC exploration
//  nBin     | 8        | No. of bins in edge-histogram in cell exploration
//  OptRej   | 1        | OptRej = 0, weighted; OptRej=1, wt=1 MC events
//  OptDrive | 2        | Maximum weight reduction, =1 for variance reduction
//  EvPerBin | 25       | Maximum number of the effective wt=1 events/bin,
//           |          | EvPerBin=0 deactivates this option
//  Chat     | 1        | =0,1,2 is the ``chat level'' in the standard output
//  MaxWtRej | 1.1      | Maximum weight used to get w=1 MC events
//------------------------------------------------------------------------------
// The above can be redefined before calling 'Initialize()' method,
// for instance FoamObject->SetkDim(15) sets dimension of the distribution to 15.
// Only kDim HAS TO BE redefined, the other parameters may be left at their defaults.
// nCell may be increased up to about million cells for wildly peaked distributions.
// Increasing nSampl sometimes helps, but it may cost CPU time.
// MaxWtRej may need to be increased for wild a distribution, while using OptRej=0.
//
// --------------------------------------------------------------------
// Past versions of FOAM: August 2003, v.1.00; September 2003 v.1.01
// Adopted starting from FOAM-2.06 by P. Sawicki
// --------------------------------------------------------------------
// Users of FOAM are kindly requested to cite the following work:
// S. Jadach, Computer Physics Communications 152 (2003) 55.
//
//______________________________________________________________________________

#include "TFoam.h"
#include "TFoamIntegrand.h"
#include "TFoamMaxwt.h"
#include "TFoamVect.h"
#include "TFoamCell.h"
#include "Riostream.h"
#include "TH1.h"
#include "TRefArray.h"
#include "TMethodCall.h"
#include "TRandom.h"
#include "TMath.h"
#include "TInterpreter.h"

ClassImp(TFoam);

//FFFFFF  BoX-FORMATs for nice and flexible outputs
#define BXOPE cout<<\
"FFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFF"<<endl<<\
"F                                                                              F"<<endl
#define BXTXT(text) cout<<\
"F                   "<<setw(40)<<         text           <<"                   F"<<endl
#define BX1I(name,numb,text) cout<<\
"F "<<setw(10)<<name<<" = "<<setw(10)<<numb<<" = "          <<setw(50)<<text<<" F"<<endl
#define BX1F(name,numb,text)     cout<<"F "<<setw(10)<<name<<\
          " = "<<setw(15)<<setprecision(8)<<numb<<"   =    "<<setw(40)<<text<<" F"<<endl
#define BX2F(name,numb,err,text) cout<<"F "<<setw(10)<<name<<\
" = "<<setw(15)<<setprecision(8)<<numb<<" +- "<<setw(15)<<setprecision(8)<<err<<\
                                                      "  = "<<setw(25)<<text<<" F"<<endl
#define BXCLO cout<<\
"F                                                                              F"<<endl<<\
"FFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFF"<<endl
  //FFFFFF  BoX-FORMATs ends here

static const Double_t gHigh= 1.0e150;
static const Double_t gVlow=-1.0e150;

#define SW2 setprecision(7) << setw(12)

//________________________________________________________________________________________________
TFoam::TFoam()
{
  // Default constructor for streamer, user should not use it.
   fDim      = 0;
   fNoAct    = 0;
   fNCells   = 0;
   fRNmax    = 0;
   fMaskDiv  = 0;
   fInhiDiv  = 0;
   fXdivPRD  = 0;
   fCells    = 0;
   fAlpha    = 0;
   fCellsAct = 0;
   fPrimAcu  = 0;
   fHistEdg  = 0;
   fHistWt   = 0;
   fHistDbg  = 0;
   fMCMonit  = 0;
   fRho      = 0;  // Integrand function
   fMCvect   = 0;
   fRvec     = 0;
   fPseRan   = 0;  // generator of pseudorandom numbers

}
//_________________________________________________________________________________________________
TFoam::TFoam(const Char_t* Name)
{
// User constructor, to be employed by the user

   if(strlen(Name)  >129) {
      Error("TFoam","Name too long %s \n",Name);
   }
   fName=Name;                                            // Class name
   fDate="  Release date:  2005.04.10";                   // Release date
   fVersion= "1.02M";                                     // Release version
   fMaskDiv  = 0;             // Dynamic Mask for  cell division, h-cubic
   fInhiDiv  = 0;             // Flag allowing to inhibit cell division in certain projection/edge
   fXdivPRD  = 0;             // Lists of division values encoded in one vector per direction
   fCells    = 0;
   fAlpha    = 0;
   fCellsAct = 0;
   fPrimAcu  = 0;
   fHistEdg  = 0;
   fHistWt   = 0;
   fHistDbg  = 0;
   fDim     = 0;                // dimension of hyp-cubical space
   fNCells   = 1000;             // Maximum number of Cells,    is usually re-set
   fNSampl   = 200;              // No of sampling when dividing cell
   fOptPRD   = 0;                // General Option switch for PRedefined Division, for quick check
   fOptDrive = 2;                // type of Drive =1,2 for TrueVol,Sigma,WtMax
   fChat     = 1;                // Chat=0,1,2 chat level in output, Chat=1 normal level
   fOptRej   = 1;                // OptRej=0, wted events; OptRej=1, wt=1 events
   //------------------------------------------------------
   fNBin     = 8;                // binning of edge-histogram in cell exploration
   fEvPerBin =25;                // maximum no. of EFFECTIVE event per bin, =0 option is inactive
   //------------------------------------------------------
   fNCalls = 0;                  // No of function calls
   fNEffev = 0;                  // Total no of eff. wt=1 events in build=up
   fLastCe =-1;                  // Index of the last cell
   fNoAct  = 0;                  // No of active cells (used in MC generation)
   fWtMin = gHigh;               // Minimal weight
   fWtMax = gVlow;               // Maximal weight
   fMaxWtRej =1.10;              // Maximum weight in rejection for getting wt=1 events
   fPseRan   = 0;                // Initialize private copy of random number generator
   fMCMonit  = 0;                // MC efficiency monitoring
   fRho = 0;                     // pointer to abstract class providing function to integrate
   fMethodCall=0;                // ROOT's pointer to global distribution function
}

//_______________________________________________________________________________________________
TFoam::~TFoam()
{
// Default destructor
//  cout<<" DESTRUCTOR entered "<<endl;
   Int_t i;

   if(fCells!= 0) {
      for(i=0; i<fNCells; i++) delete fCells[i]; // TFoamCell*[]
      delete [] fCells;
   }
   if (fCellsAct) delete fCellsAct ; // WVE FIX LEAK
   delete [] fRvec;    //double[]
   delete [] fAlpha;   //double[]
   delete [] fMCvect;  //double[]
   delete [] fPrimAcu; //double[]
   delete [] fMaskDiv; //int[]
   delete [] fInhiDiv; //int[]
 
   if( fXdivPRD!= 0) {
      for(i=0; i<fDim; i++) delete fXdivPRD[i]; // TFoamVect*[]
      delete [] fXdivPRD;
   }
   delete fMCMonit;
   delete fHistWt;

   // delete histogram arrays
   if (fHistEdg) { 
      fHistEdg->Delete(); 
      delete fHistEdg; 
   }
   if (fHistDbg) { 
      fHistDbg->Delete(); 
      delete fHistDbg; 
   }
}


//_____________________________________________________________________________________________
TFoam::TFoam(const TFoam &From): TObject(From)
{
// Copy Constructor  NOT IMPLEMENTED (NEVER USED)
   Error("TFoam", "COPY CONSTRUCTOR NOT IMPLEMENTED \n");
}

//_____________________________________________________________________________________________
void TFoam::Initialize(TRandom *PseRan, TFoamIntegrand *fun )
{
// Basic initialization of FOAM invoked by the user. Mandatory!
// ============================================================
// This method starts the process of the cell build-up.
// User must invoke Initialize with two arguments or Initialize without arguments.
// This is done BEFORE generating first MC event and AFTER allocating FOAM object
// and reseting (optionally) its internal parameters/switches.
// The overall operational scheme of the FOAM is the following:
//BEGIN_HTML <!--
/* -->
<img src="gif/foam_schema2.gif">
<!--*/
// -->END_HTML
//
// This method invokes several other methods:
// ==========================================
// InitCells initializes memory storage for cells and begins exploration process
// from the root cell. The empty cells are allocated/filled using  CellFill.
// The procedure Grow which loops over cells, picks up the cell with the biggest
// ``driver integral'', see Comp. Phys. Commun. 152 152 (2003) 55 for explanations,
// with the help of PeekMax procedure. The chosen cell is split using Divide.
// Subsequently, the procedure Explore called by the Divide
// (and by InitCells for the root cell) does the most important
// job in the FOAM object build-up: it performs a small MC run for each
// newly allocated daughter cell.
// Explore calculates how profitable the future split of the cell will be
// and defines the optimal cell division geometry with the help of Carver or Varedu
// procedures, for maximum weight or variance optimization respectively.
// All essential results of the exploration are written into
// the explored cell object. At the very end of the foam build-up,
// Finally, MakeActiveList is invoked to create a list of pointers to
// all active cells, for the purpose of the quick access during the MC generation.
// The procedure Explore employs MakeAlpha to generate random coordinates
// inside a given cell with the uniform distribution.
// The above sequence of the procedure calls is depicted in the following figure:
//BEGIN_HTML <!--
/* -->
<img src="gif/foam_Initialize_schema.gif">
<!--*/
// -->END_HTML

   SetPseRan(PseRan);
   SetRho(fun);
   Initialize();
}

//_______________________________________________________________________________________
void TFoam::Initialize()
{
// Basic initialization of FOAM invoked by the user.
// IMPORTANT: Random number generator and the distribution object has to be
// provided using SetPseRan and SetRho prior to invoking this initializator!

   Bool_t addStatus = TH1::AddDirectoryStatus();
   TH1::AddDirectory(kFALSE);
   Int_t i;

   if(fChat>0){
      BXOPE;
      BXTXT("****************************************");
      BXTXT("******      TFoam::Initialize    ******");
      BXTXT("****************************************");
      BXTXT(fName);
      BX1F("  Version",fVersion,  fDate);
      BX1I("     kDim",fDim,     " Dimension of the hyper-cubical space             ");
      BX1I("   nCells",fNCells,   " Requested number of Cells (half of them active)  ");
      BX1I("   nSampl",fNSampl,   " No of MC events in exploration of a cell         ");
      BX1I("     nBin",fNBin,     " No of bins in histograms, MC exploration of cell ");
      BX1I(" EvPerBin",fEvPerBin, " Maximum No effective_events/bin, MC exploration  ");
      BX1I(" OptDrive",fOptDrive, " Type of Driver   =1,2 for Sigma,WtMax            ");
      BX1I("   OptRej",fOptRej,   " MC rejection on/off for OptRej=0,1               ");
      BX1F(" MaxWtRej",fMaxWtRej, " Maximum wt in rejection for wt=1 evts");
      BXCLO;
   }

   if(fPseRan==0) Error("Initialize", "Random number generator not set \n");
   if(fRho==0 && fMethodCall==0 ) Error("Initialize", "Distribution function not set \n");
   if(fDim==0) Error("Initialize", "Zero dimension not allowed \n");

   /////////////////////////////////////////////////////////////////////////
   //                   ALLOCATE SMALL LISTS                              //
   //  it is done globally, not for each cell, to save on allocation time //
   /////////////////////////////////////////////////////////////////////////
   fRNmax= fDim+1;
   fRvec = new Double_t[fRNmax];   // Vector of random numbers
   if(fRvec==0)  Error("Initialize", "Cannot initialize buffer fRvec \n");

   if(fDim>0){
      fAlpha = new Double_t[fDim];    // sum<1 for internal parametrization of the simplex
      if(fAlpha==0)  Error("Initialize", "Cannot initialize buffer fAlpha \n" );
   }
   fMCvect = new Double_t[fDim]; // vector generated in the MC run
   if(fMCvect==0)  Error("Initialize", "Cannot initialize buffer fMCvect  \n" );

   //====== List of directions inhibited for division
   if(fInhiDiv == 0){
      fInhiDiv = new Int_t[fDim];
      for(i=0; i<fDim; i++) fInhiDiv[i]=0;
   }
   //====== Dynamic mask used in Explore for edge determination
   if(fMaskDiv == 0){
      fMaskDiv = new Int_t[fDim];
      for(i=0; i<fDim; i++) fMaskDiv[i]=1;
   }
   //====== List of predefined division values in all directions (initialized as empty)
   if(fXdivPRD == 0){
      fXdivPRD = new TFoamVect*[fDim];
      for(i=0; i<fDim; i++)  fXdivPRD[i]=0; // Artificially extended beyond fDim
   }
   //====== Initialize list of histograms
   fHistWt  = new TH1D("HistWt","Histogram of MC weight",100,0.0, 1.5*fMaxWtRej); // MC weight
   fHistEdg = new TObjArray(fDim);           // Initialize list of histograms
   TString hname;
   TString htitle;
   for(i=0;i<fDim;i++){
      hname=fName+TString("_HistEdge_");
      hname+=i;
      htitle=TString("Edge Histogram No. ");
      htitle+=i;
      //cout<<"i= "<<i<<"  hname= "<<hname<<"  htitle= "<<htitle<<endl;
      (*fHistEdg)[i] = new TH1D(hname.Data(),htitle.Data(),fNBin,0.0, 1.0); // Initialize histogram for each edge
      ((TH1D*)(*fHistEdg)[i])->Sumw2();
   }
   //======  extra histograms for debug purposes
   fHistDbg = new TObjArray(fDim);         // Initialize list of histograms
   for(i=0;i<fDim;i++){
      hname=fName+TString("_HistDebug_");
      hname+=i;
      htitle=TString("Debug Histogram ");
      htitle+=i;
      (*fHistDbg)[i] = new TH1D(hname.Data(),htitle.Data(),fNBin,0.0, 1.0); // Initialize histogram for each edge
   }

   // ||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||| //
   //                     BUILD-UP of the FOAM                            //
   // ||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||| //
   //
   //        Define and explore root cell(s)
   InitCells();
   //        PrintCells(); cout<<" ===== after InitCells ====="<<endl;
   Grow();
   //        PrintCells(); cout<<" ===== after Grow      ====="<<endl;

   MakeActiveList(); // Final Preperations for the M.C. generation

   // Preperations for the M.C. generation
   fSumWt  = 0.0;               // M.C. generation sum of Wt
   fSumWt2 = 0.0;               // M.C. generation sum of Wt**2
   fSumOve = 0.0;               // M.C. generation sum of overweighted
   fNevGen = 0.0;               // M.C. generation sum of 1d0
   fWtMax  = gVlow;               // M.C. generation maximum wt
   fWtMin  = gHigh;               // M.C. generation minimum wt
   fMCresult=fCells[0]->GetIntg(); // M.C. Value of INTEGRAL,temporary assignment
   fMCresult=fCells[0]->GetIntg(); // M.C. Value of INTEGRAL,temporary assignment
   fMCerror =fCells[0]->GetIntg(); // M.C. Value of ERROR   ,temporary assignment
   fMCMonit = new TFoamMaxwt(5.0,1000);  // monitoring M.C. efficiency
   //
   if(fChat>0){
      Double_t driver = fCells[0]->GetDriv();
      BXOPE;
      BXTXT("***  TFoam::Initialize FINISHED!!!  ***");
      BX1I("    nCalls",fNCalls,  "Total number of function calls         ");
      BX1F("    XPrime",fPrime,   "Primary total integral                 ");
      BX1F("    XDiver",driver,    "Driver  total integral                 ");
      BX1F("  mcResult",fMCresult,"Estimate of the true MC Integral       ");
      BXCLO;
   }
   if(fChat==2) PrintCells();
   TH1::AddDirectory(addStatus);
} // Initialize

//_______________________________________________________________________________________
void TFoam::InitCells()
{
// Internal subprogram used by Initialize.
// It initializes "root part" of the FOAM of the tree of cells.

   Int_t i;

   fLastCe =-1;                             // Index of the last cell
   if(fCells!= 0) {
      for(i=0; i<fNCells; i++) delete fCells[i];
      delete [] fCells;
   }
   //
   fCells = new TFoamCell*[fNCells];
   for(i=0;i<fNCells;i++){
      fCells[i]= new TFoamCell(fDim); // Allocate BIG list of cells
      fCells[i]->SetSerial(i);
   }
   if(fCells==0) Error("InitCells", "Cannot initialize CELLS \n"  );

   /////////////////////////////////////////////////////////////////////////////
   //              Single Root Hypercube                                      //
   /////////////////////////////////////////////////////////////////////////////
   CellFill(1,   0);  //  0-th cell ACTIVE

   // Exploration of the root cell(s)
   for(Long_t iCell=0; iCell<=fLastCe; iCell++){
      Explore( fCells[iCell] );               // Exploration of root cell(s)
   }
}//InitCells

//_______________________________________________________________________________________
Int_t TFoam::CellFill(Int_t Status, TFoamCell *parent)
{
// Internal subprogram used by Initialize.
// It initializes content of the newly allocated active cell.

   TFoamCell *cell;
   if (fLastCe==fNCells){
      Error( "CellFill", "Too many cells\n");
   }
   fLastCe++;   // 0-th cell is the first
   if (Status==1) fNoAct++;

   cell = fCells[fLastCe];

   cell->Fill(Status, parent, 0, 0);

   cell->SetBest( -1);         // pointer for planning division of the cell
   cell->SetXdiv(0.5);         // factor for division
   Double_t xInt2,xDri2;
   if(parent!=0){
      xInt2  = 0.5*parent->GetIntg();
      xDri2  = 0.5*parent->GetDriv();
      cell->SetIntg(xInt2);
      cell->SetDriv(xDri2);
   }else{
      cell->SetIntg(0.0);
      cell->SetDriv(0.0);
   }
   return fLastCe;
}

//______________________________________________________________________________________
void TFoam::Explore(TFoamCell *cell)
{
// Internal subprogram used by Initialize.
// It explores newly defined cell with help of special short MC sampling.
// As a result, estimates of true and drive volume is defined/determined
// Average and dispersion of the weight distribution will is found along
// each edge and the best edge (minimum dispersion, best maximum weight)
// is memorized for future use.
// The optimal division point for eventual future cell division is
// determined/recorded. Recorded are also minimum and maximum weight etc.
// The volume estimate in all (inactive) parent cells is updated.
// Note that links to parents and initial volume = 1/2 parent has to be
// already defined prior to calling this routine.

   Double_t wt, dx, xBest, yBest;
   Double_t intOld, driOld;

   Long_t iev;
   Double_t nevMC;
   Int_t i, j, k;
   Int_t nProj, kBest;
   Double_t ceSum[5], xproj;

   TFoamVect  cellSize(fDim);
   TFoamVect  cellPosi(fDim);

   cell->GetHcub(cellPosi,cellSize);

   TFoamCell  *parent;

   Double_t *xRand = new Double_t[fDim];

   Double_t *volPart=0;

   cell->CalcVolume();
   dx = cell->GetVolume();
   intOld = cell->GetIntg(); //memorize old values,
   driOld = cell->GetDriv(); //will be needed for correcting parent cells


   /////////////////////////////////////////////////////
   //    Special Short MC sampling to probe cell      //
   /////////////////////////////////////////////////////
   ceSum[0]=0;
   ceSum[1]=0;
   ceSum[2]=0;
   ceSum[3]=gHigh;  //wtmin
   ceSum[4]=gVlow;  //wtmax
   //
   for(i=0;i<fDim;i++) ((TH1D *)(*fHistEdg)[i])->Reset(); // Reset histograms
   fHistWt->Reset();
   //
   // ||||||||||||||||||||||||||BEGIN MC LOOP|||||||||||||||||||||||||||||
   Double_t nevEff=0.;
   for(iev=0;iev<fNSampl;iev++){
      MakeAlpha();               // generate uniformly vector inside hypercube

      if(fDim>0){
      for(j=0; j<fDim; j++)
         xRand[j]= cellPosi[j] +fAlpha[j]*(cellSize[j]);
      }

      wt=dx*Eval(xRand);

      nProj = 0;
      if(fDim>0) {
         for(k=0; k<fDim; k++) {
            xproj =fAlpha[k];
            ((TH1D *)(*fHistEdg)[nProj])->Fill(xproj,wt);
            nProj++;
         }
      }
      //
      fNCalls++;
      ceSum[0] += wt;    // sum of weights
      ceSum[1] += wt*wt; // sum of weights squared
      ceSum[2]++;        // sum of 1
      if (ceSum[3]>wt) ceSum[3]=wt;  // minimum weight;
      if (ceSum[4]<wt) ceSum[4]=wt;  // maximum weight
      // test MC loop exit condition
      nevEff = ceSum[0]*ceSum[0]/ceSum[1];
      if( nevEff >= fNBin*fEvPerBin) break;
   }   // ||||||||||||||||||||||||||END MC LOOP|||||||||||||||||||||||||||||
   //------------------------------------------------------------------
   //---  predefine logics of searching for the best division edge ---
   for(k=0; k<fDim;k++){
      fMaskDiv[k] =1;                       // default is all
      if( fInhiDiv[k]==1) fMaskDiv[k] =0; // inhibit some...
   }
   // Note that predefined division below overrule inhibition above
   kBest=-1;
   Double_t rmin,rmax,rdiv;
   if(fOptPRD) {          // quick check
      for(k=0; k<fDim; k++) {
         rmin= cellPosi[k];
         rmax= cellPosi[k] +cellSize[k];
         if( fXdivPRD[k] != 0) {
            Int_t n= (fXdivPRD[k])->GetDim();
            for(j=0; j<n; j++) {
               rdiv=(*fXdivPRD[k])[j];
               // check predefined divisions is available in this cell
               if( (rmin +1e-99 <rdiv) && (rdiv< rmax -1e-99)) {
                  kBest=k;
                  xBest= (rdiv-cellPosi[k])/cellSize[k] ;
                  goto ee05;
               }
            }
         }
      }//k
   }
   ee05:
   //------------------------------------------------------------------
   fNEffev += (Long_t)nevEff;
   nevMC          = ceSum[2];
   Double_t intTrue = ceSum[0]/(nevMC+0.000001);
   Double_t intDriv=0.;
   Double_t intPrim=0.;

   switch(fOptDrive){
   case 1:                       // VARIANCE REDUCTION
      if(kBest == -1) Varedu(ceSum,kBest,xBest,yBest); // determine the best edge,
      //intDriv =sqrt( ceSum[1]/nevMC -intTrue*intTrue ); // Older ansatz, numerically not bad
      intDriv =sqrt(ceSum[1]/nevMC) -intTrue; // Foam build-up, sqrt(<w**2>) -<w>
      intPrim =sqrt(ceSum[1]/nevMC);          // MC gen. sqrt(<w**2>) =sqrt(<w>**2 +sigma**2)
      break;
   case 2:                       // WTMAX  REDUCTION
      if(kBest == -1) Carver(kBest,xBest,yBest);  // determine the best edge
      intDriv =ceSum[4] -intTrue; // Foam build-up, wtmax-<w>
      intPrim =ceSum[4];          // MC generation, wtmax!
      break;
   default:
      Error("Explore", "Wrong fOptDrive = \n" );
   }//switch
   //=================================================================================
   //hist_Neff_distrib.Fill( fLastCe/2.0+0.01, nevEff+0.01);  //
   //hist_kBest_distrib.Fill( kBest+0.50, 1.0 ); //  debug
   //hist_xBest_distrib.Fill( xBest+0.01, 1.0 ); //  debug
   //=================================================================================
   cell->SetBest(kBest);
   cell->SetXdiv(xBest);
   cell->SetIntg(intTrue);
   cell->SetDriv(intDriv);
   cell->SetPrim(intPrim);
   // correct/update integrals in all parent cells to the top of the tree
   Double_t  parIntg, parDriv;
   for(parent = cell->GetPare(); parent!=0; parent = parent->GetPare()){
      parIntg = parent->GetIntg();
      parDriv = parent->GetDriv();
      parent->SetIntg( parIntg   +intTrue -intOld );
      parent->SetDriv( parDriv   +intDriv -driOld );
   }
   delete [] volPart;
   delete [] xRand;
   //cell->Print();
} // TFoam::Explore

//______________________________________________________________________________________
void TFoam::Varedu(Double_t ceSum[5], Int_t &kBest, Double_t &xBest, Double_t &yBest)
{
// Internal subrogram used by Initialize.
// In determines the best edge candidate and the position of the cell division plane
// in case of the variance reduction for future cell division,
// using results of the MC exploration run stored in fHistEdg

   Double_t nent   = ceSum[2];
   Double_t swAll  = ceSum[0];
   Double_t sswAll = ceSum[1];
   Double_t ssw    = sqrt(sswAll)/sqrt(nent);
   //
   Double_t swIn,swOut,sswIn,sswOut,xLo,xUp;
   kBest =-1;
   xBest =0.5;
   yBest =1.0;
   Double_t maxGain=0.0;
   // Now go over all projections kProj
   for(Int_t kProj=0; kProj<fDim; kProj++) {
      if( fMaskDiv[kProj]) {
         // initialize search over bins
         Double_t sigmIn =0.0; Double_t sigmOut =0.0;
         Double_t sswtBest = gHigh;
         Double_t gain =0.0;
         Double_t xMin=0.0; Double_t xMax=0.0;
         // Double loop over all pairs jLo<jUp
         for(Int_t jLo=1; jLo<=fNBin; jLo++) {
            Double_t aswIn=0;  Double_t asswIn=0;
            for(Int_t jUp=jLo; jUp<=fNBin;jUp++) {
               aswIn  +=     ((TH1D *)(*fHistEdg)[kProj])->GetBinContent(jUp);
               asswIn += Sqr(((TH1D *)(*fHistEdg)[kProj])->GetBinError(  jUp));
               xLo=(jLo-1.0)/fNBin;
               xUp=(jUp*1.0)/fNBin;
               swIn  =        aswIn/nent;
               swOut = (swAll-aswIn)/nent;
               sswIn = sqrt(asswIn)       /sqrt(nent*(xUp-xLo))     *(xUp-xLo);
               sswOut= sqrt(sswAll-asswIn)/sqrt(nent*(1.0-xUp+xLo)) *(1.0-xUp+xLo);
               if( (sswIn+sswOut) < sswtBest) {
                  sswtBest = sswIn+sswOut;
                  gain     = ssw-sswtBest;
                  sigmIn   = sswIn -swIn;  // Debug
                  sigmOut  = sswOut-swOut; // Debug
                  xMin    = xLo;
                  xMax    = xUp;
               }
            }//jUp
         }//jLo
         Int_t iLo = (Int_t) (fNBin*xMin);
         Int_t iUp = (Int_t) (fNBin*xMax);
         //----------DEBUG printout
         //cout<<"@@@@@  xMin xMax = "<<xMin   <<" "<<xMax<<"  iLo= "<<iLo<<"  iUp= "<<iUp;
         //cout<<"  sswtBest/ssw= "<<sswtBest/ssw<<"  Gain/ssw= "<< Gain/ssw<<endl;
         //----------DEBUG auxilary Plot
         for(Int_t iBin=1;iBin<=fNBin;iBin++) {
            if( ((iBin-0.5)/fNBin > xMin) && ((iBin-0.5)/fNBin < xMax) ){
               ((TH1D *)(*fHistDbg)[kProj])->SetBinContent(iBin,sigmIn/(xMax-xMin));
            } else {
               ((TH1D *)(*fHistDbg)[kProj])->SetBinContent(iBin,sigmOut/(1-xMax+xMin));
            }
         }
         if(gain>=maxGain) {
            maxGain=gain;
            kBest=kProj; // <--- !!!!! The best edge
            xBest=xMin;
            yBest=xMax;
            if(iLo == 0)     xBest=yBest; // The best division point
            if(iUp == fNBin) yBest=xBest; // this is not really used
         }
      }
   } //kProj
   //----------DEBUG printout
   //cout<<"@@@@@@@>>>>> kBest= "<<kBest<<"  maxGain/ssw= "<< maxGain/ssw<<endl;
   if( (kBest >= fDim) || (kBest<0) ) Error("Varedu", "Something wrong with kBest \n" );
}          //TFoam::Varedu

//________________________________________________________________________________________
void TFoam::Carver(Int_t &kBest, Double_t &xBest, Double_t &yBest)
{
// Internal subrogram used by Initialize.
// Determines the best edge-candidate and the position of the division plane
// for the future cell division, in the case of the optimization of the maximum weight.
// It exploits results of the cell MC exploration run stored in fHistEdg.

   Int_t    kProj,iBin;
   Double_t carve,carvTot,carvMax,carvOne,binMax,binTot,primTot,primMax;
   Int_t    jLow,jUp,iLow,iUp;
   Double_t theBin;
   Int_t    jDivi; // TEST
   Int_t j;

   Double_t *bins  = new Double_t[fNBin];      // bins of histogram for single  PROJECTION
   if(bins==0)    Error("Carver", "Cannot initialize buffer Bins \n" );

   kBest =-1;
   xBest =0.5;
   yBest =1.0;
   carvMax = gVlow;
   primMax = gVlow;
   for(kProj=0; kProj<fDim; kProj++)
      if( fMaskDiv[kProj] ){
      //if( kProj==1 ){
      //cout<<"==================== Carver histogram: kProj ="<<kProj<<"==================="<<endl;
      //((TH1D *)(*fHistEdg)[kProj])->Print("all");
      binMax = gVlow;
      for(iBin=0; iBin<fNBin;iBin++){
         bins[iBin]= ((TH1D *)(*fHistEdg)[kProj])->GetBinContent(iBin+1);
         binMax = TMath::Max( binMax, bins[iBin]);       // Maximum content/bin
      }
      if(binMax < 0 ) {       //case of empty cell
         delete [] bins;
         return;
      }
      carvTot = 0.0;
      binTot  = 0.0;
      for(iBin=0;iBin<fNBin;iBin++){
         carvTot = carvTot + (binMax-bins[iBin]);     // Total Carve (more stable)
         binTot  +=bins[iBin];
      }
      primTot = binMax*fNBin;
      //cout <<"Carver:  CarvTot "<<CarvTot<< "    primTot "<<primTot<<endl;
      jLow =0;
      jUp  =fNBin-1;
      carvOne = gVlow;
      Double_t yLevel = gVlow;
      for(iBin=0; iBin<fNBin;iBin++) {
         theBin = bins[iBin];
         //-----  walk to the left and find first bin > theBin
         iLow = iBin;
         for(j=iBin; j>-1; j-- ) {
            if(theBin< bins[j]) break;
            iLow = j;
         }
         //iLow = iBin;
         //if(iLow>0)     while( (theBin >= bins[iLow-1])&&(iLow >0) ){iLow--;} // horror!!!
         //------ walk to the right and find first bin > theBin
         iUp  = iBin;
         for(j=iBin; j<fNBin; j++) {
            if(theBin< bins[j]) break;
            iUp = j;
         }
         //iUp  = iBin;
         //if(iUp<fNBin-1) while( (theBin >= bins[iUp+1])&&( iUp<fNBin-1 ) ){iUp++;} // horror!!!
         //
         carve = (iUp-iLow+1)*(binMax-theBin);
         if( carve > carvOne) {
            carvOne = carve;
            jLow = iLow;
            jUp  = iUp;
            yLevel = theBin;
         }
      }//iBin
      if( carvTot > carvMax) {
         carvMax   = carvTot;
         primMax   = primTot;
         //cout <<"Carver:   primMax "<<primMax<<endl;
         kBest = kProj;    // Best edge
         xBest = ((Double_t)(jLow))/fNBin;
         yBest = ((Double_t)(jUp+1))/fNBin;
         if(jLow == 0 )       xBest = yBest;
         if(jUp  == fNBin-1) yBest = xBest;
         // division ratio in units of 1/fNBin, testing
         jDivi = jLow;
         if(jLow == 0 )     jDivi=jUp+1;
      }
      //======  extra histograms for debug purposes
      //cout<<"kProj= "<<kProj<<" jLow= "<<jLow<<" jUp= "<<jUp<<endl;
      for(iBin=0;    iBin<fNBin;  iBin++)
         ((TH1D *)(*fHistDbg)[kProj])->SetBinContent(iBin+1,binMax);
      for(iBin=jLow; iBin<jUp+1;   iBin++)
         ((TH1D *)(*fHistDbg)[kProj])->SetBinContent(iBin+1,yLevel);
   }//kProj
   if( (kBest >= fDim) || (kBest<0) ) Error("Carver", "Something wrong with kBest \n" );
   delete [] bins;
}          //TFoam::Carver

//______________________________________________________________________________________________
void TFoam::MakeAlpha()
{
// Internal subrogram used by Initialize.
// Provides random vector Alpha  0< Alpha(i) < 1
   Int_t k;
   if(fDim<1) return;

   // simply generate and load kDim uniform random numbers
   fPseRan->RndmArray(fDim,fRvec);   // kDim random numbers needed
   for(k=0; k<fDim; k++) fAlpha[k] = fRvec[k];
} //MakeAlpha


//_____________________________________________________________________________________________
void TFoam::Grow()
{
// Internal subrogram used by Initialize.
// It grow new cells by the binary division process.

   Long_t iCell;
   TFoamCell* newCell;

   while ( (fLastCe+2) < fNCells ) {  // this condition also checked inside Divide
      iCell   = PeekMax();            // peek up cell with maximum driver integral
      if( (iCell<0) || (iCell>fLastCe) ) Error("Grow", "Wrong iCell \n");
      newCell = fCells[iCell];

      if(fLastCe !=0) {
         Int_t kEcho=10;
         if(fLastCe>=10000) kEcho=100;
         if( (fLastCe%kEcho)==0) {
	   if (fChat>0) {
	     if(fDim<10)
	       cout<<fDim<<flush;
             else
	       cout<<"."<<flush;
	     if( (fLastCe%(100*kEcho))==0)  cout<<"|"<<fLastCe<<endl<<flush;
	   }
	 }
      }
      if( Divide( newCell )==0) break;  // and divide it into two
   }
   if (fChat>0) {
     cout<<endl<<flush;
   }
   CheckAll(0);   // set arg=1 for more info
}// Grow

//_____________________________________________________________________________________________
Long_t  TFoam::PeekMax()
{
// Internal subprogram used by Initialize.
// It finds cell with maximal driver integral for the purpose of the division.

   Long_t  i;
   Long_t iCell = -1;
   Double_t  drivMax, driv;

   drivMax = gVlow;
   for(i=0; i<=fLastCe; i++) {//without root
      if( fCells[i]->GetStat() == 1 ) {
         driv =  TMath::Abs( fCells[i]->GetDriv());
         //cout<<"PeekMax: Driv = "<<driv<<endl;
         if(driv > drivMax) {
            drivMax = driv;
            iCell = i;
         }
      }
   }
   //  cout << 'TFoam_PeekMax: iCell=' << iCell << endl;
   if (iCell == -1)
      cout << "STOP in TFoam::PeekMax: not found iCell=" <<  iCell << endl;
   return(iCell);
}                 // TFoam_PeekMax

//_____________________________________________________________________________________________
Int_t TFoam::Divide(TFoamCell *cell)
{
// Internal subrogram used by Initialize.
// It divides cell iCell into two daughter cells.
// The iCell is retained and tagged as inactive, daughter cells are appended
// at the end of the buffer.
// New vertex is added to list of vertices.
// List of active cells is updated, iCell removed, two daughters added
// and their properties set with help of MC sampling (TFoam_Explore)
// Returns Code RC=-1 of buffer limit is reached,  fLastCe=fnBuf.

   Double_t xdiv;
   Int_t   kBest;

   if(fLastCe+1 >= fNCells) Error("Divide", "Buffer limit is reached, fLastCe=fnBuf \n");

   cell->SetStat(0); // reset cell as inactive
   fNoAct--;

   xdiv  = cell->GetXdiv();
   kBest = cell->GetBest();
   if( kBest<0 || kBest>=fDim ) Error("Divide", "Wrong kBest \n");

   //////////////////////////////////////////////////////////////////
   //           define two daughter cells (active)                 //
   //////////////////////////////////////////////////////////////////

   Int_t d1 = CellFill(1,   cell);
   Int_t d2 = CellFill(1,   cell);
   cell->SetDau0((fCells[d1]));
   cell->SetDau1((fCells[d2]));
   Explore( (fCells[d1]) );
   Explore( (fCells[d2]) );
   return 1;
} // TFoam_Divide


//_________________________________________________________________________________________
void TFoam::MakeActiveList()
{
// Internal subrogram used by Initialize.
// It finds out number of active cells fNoAct,
// creates list of active cell fCellsAct and primary cumulative fPrimAcu.
// They are used during the MC generation to choose randomly an active cell.

   Long_t n, iCell;
   Double_t sum;
   // flush previous result
   if(fPrimAcu  != 0) delete [] fPrimAcu;
   if(fCellsAct != 0) delete fCellsAct;

   // Allocate tables of active cells
   fCellsAct = new TRefArray();

   // Count Active cells and find total Primary
   // Fill-in tables of active cells

   fPrime = 0.0; n = 0;
   for(iCell=0; iCell<=fLastCe; iCell++) { 
      if (fCells[iCell]->GetStat()==1) {
         fPrime += fCells[iCell]->GetPrim();
         fCellsAct->Add(fCells[iCell]);
         n++;
      }
   }

   if(fNoAct != n)  Error("MakeActiveList", "Wrong fNoAct               \n"  );
   if(fPrime == 0.) Error("MakeActiveList", "Integrand function is zero  \n"  );

   fPrimAcu  = new  Double_t[fNoAct]; // cumulative primary for MC generation
   if( fCellsAct==0 || fPrimAcu==0 ) Error("MakeActiveList", "Cant allocate fCellsAct or fPrimAcu \n");

   sum =0.0;
   for(iCell=0; iCell<fNoAct; iCell++) {
      sum = sum + ( (TFoamCell *) (fCellsAct->At(iCell)) )->GetPrim()/fPrime;
      fPrimAcu[iCell]=sum;
   }

} //MakeActiveList

//__________________________________________________________________________________________
void TFoam::ResetPseRan(TRandom *PseRan)
{
// User may optionally reset random number generator using this method
// Usually it is done when FOAM object is restored from the disk.
// IMPORTANT: this method deletes existing  random number generator registered in the FOAM object.
// In particular such an object is created by the streamer during the disk-read operation.

   if(fPseRan) {
      Info("ResetPseRan", "Resetting random number generator  \n");
      delete fPseRan;
   }
   SetPseRan(PseRan);
}

//__________________________________________________________________________________________
void TFoam::SetRho(TFoamIntegrand *fun)
{
// User may use this method to set (register) random number generator used by
// the given instance of the FOAM event generator. Note that single r.n. generator
// may serve several FOAM objects.

   if (fun)
      fRho=fun;
   else
      Error("SetRho", "Bad function \n" );
}

//__________________________________________________________________________________________
void TFoam::ResetRho(TFoamIntegrand *fun)
{
// User may optionally reset the distribution using this method
// Usually it is done when FOAM object is restored from the disk.
// IMPORTANT: this method deletes existing  distribution object registered in the FOAM object.
// In particular such an object is created by the streamer diring the disk-read operation.
// This method is used only in very special cases, because the distribution in most cases
// should be "owned" by the FOAM object and should not be replaced by another one after initialization.

   if(fRho) {
      Info("ResetRho", "!!! Resetting distribution function  !!!\n");
      delete fRho;
   }
   SetRho(fun);
}

//__________________________________________________________________________________________
void TFoam::SetRhoInt(void *fun)
{
// User may use this to set pointer to the global function (not descending
// from TFoamIntegrand) serving as a distribution for FOAM.
// It is useful for simple interactive applications.
// Note that persistency for FOAM object will not work in the case of such
// a distribution.

   const char *namefcn = gCint->Getp2f2funcname(fun); //name of integrand function
   if(namefcn) {
      fMethodCall=new TMethodCall();
      fMethodCall->InitWithPrototype(namefcn, "Int_t, Double_t *");
   }
   fRho=0;
}

//__________________________________________________________________________________________
Double_t TFoam::Eval(Double_t *xRand)
{
// Internal subprogram.
// Evaluates distribution to be generated.

   Double_t result;

   if(!fRho) {   //interactive mode
      Long_t paramArr[3];
      paramArr[0]=(Long_t)fDim;
      paramArr[1]=(Long_t)xRand;
      fMethodCall->SetParamPtrs(paramArr);
      fMethodCall->Execute(result);
   } else {       //compiled mode
      result=fRho->Density(fDim,xRand);
   }

   return result;
}

//___________________________________________________________________________________________
void TFoam::GenerCel2(TFoamCell *&pCell)
{
// Internal subprogram.
// Return randomly chosen active cell with probability equal to its
// contribution into total driver integral using interpolation search.

   Long_t  lo, hi, hit;
   Double_t fhit, flo, fhi;
   Double_t random;

   random=fPseRan->Rndm();
   lo  = 0;              hi =fNoAct-1;
   flo = fPrimAcu[lo];  fhi=fPrimAcu[hi];
   while(lo+1<hi) {
      hit = lo + (Int_t)( (hi-lo)*(random-flo)/(fhi-flo)+0.5);
      if (hit<=lo)
         hit = lo+1;
      else if(hit>=hi)
         hit = hi-1;
      fhit=fPrimAcu[hit];
      if (fhit>random) {
         hi = hit;
         fhi = fhit;
      } else {
         lo = hit;
         flo = fhit;
      }
   }
   if (fPrimAcu[lo]>random)
      pCell = (TFoamCell *) fCellsAct->At(lo);
   else
      pCell = (TFoamCell *) fCellsAct->At(hi);
}       // TFoam::GenerCel2


//___________________________________________________________________________________________
void TFoam::MakeEvent(void)
{
// User subprogram.
// It generates randomly point/vector according to user-defined distribution.
// Prior initialization with help of Initialize() is mandatory.
// Generated MC point/vector is available using GetMCvect and the MC weight with GetMCwt.
// MC point is generated with wt=1 or with variable weight, see OptRej switch.

   Int_t      j;
   Double_t   wt,dx,mcwt;
   TFoamCell *rCell;
   //
   //********************** MC LOOP STARS HERE **********************
ee0:
   GenerCel2(rCell);   // choose randomly one cell

   MakeAlpha();

   TFoamVect  cellPosi(fDim); TFoamVect  cellSize(fDim);
   rCell->GetHcub(cellPosi,cellSize);
   for(j=0; j<fDim; j++)
      fMCvect[j]= cellPosi[j] +fAlpha[j]*cellSize[j];
   dx = rCell->GetVolume();      // Cartesian volume of the Cell
   //  weight average normalized to PRIMARY integral over the cell

   wt=dx*Eval(fMCvect);

   mcwt = wt / rCell->GetPrim();  // PRIMARY controls normalization
   fNCalls++;
   fMCwt   =  mcwt;
   // accumulation of statistics for the main MC weight
   fSumWt  += mcwt;           // sum of Wt
   fSumWt2 += mcwt*mcwt;      // sum of Wt**2
   fNevGen++;                 // sum of 1d0
   fWtMax  =  TMath::Max(fWtMax, mcwt);   // maximum wt
   fWtMin  =  TMath::Min(fWtMin, mcwt);   // minimum wt
   fMCMonit->Fill(mcwt);
   fHistWt->Fill(mcwt,1.0);          // histogram
   //*******  Optional rejection ******
   if(fOptRej == 1) {
      Double_t random;
      random=fPseRan->Rndm();
      if( fMaxWtRej*random > fMCwt) goto ee0;  // Wt=1 events, internal rejection
      if( fMCwt<fMaxWtRej ) {
         fMCwt = 1.0;                  // normal Wt=1 event
      } else {
         fMCwt = fMCwt/fMaxWtRej;    // weight for overweighted events! kept for debug
         fSumOve += fMCwt-fMaxWtRej; // contribution of overweighted
      }
   }
   //********************** MC LOOP ENDS HERE **********************
} // MakeEvent

//_________________________________________________________________________________
void TFoam::GetMCvect(Double_t *MCvect)
{
// User may get generated MC point/vector with help of this method

   for ( Int_t k=0 ; k<fDim ; k++) *(MCvect +k) = fMCvect[k];
}//GetMCvect

//___________________________________________________________________________________
Double_t TFoam::GetMCwt(void)
{
// User may get weight MC weight using this method

   return(fMCwt);
}
//___________________________________________________________________________________
void TFoam::GetMCwt(Double_t &mcwt)
{
// User may get weight MC weight using this method

   mcwt=fMCwt;
}

//___________________________________________________________________________________
Double_t TFoam::MCgenerate(Double_t *MCvect)
{
// User subprogram which generates MC event and returns MC weight

   MakeEvent();
   GetMCvect(MCvect);
   return(fMCwt);
}//MCgenerate

//___________________________________________________________________________________
void TFoam::GetIntegMC(Double_t &mcResult, Double_t &mcError)
{
// User subprogram.
// It provides the value of the integral calculated from the averages of the MC run
// May be called after (or during) the MC run.

   Double_t mCerelat;
   mcResult = 0.0;
   mCerelat = 1.0;
   if (fNevGen>0) {
      mcResult = fPrime*fSumWt/fNevGen;
      mCerelat = sqrt( fSumWt2/(fSumWt*fSumWt) - 1/fNevGen);
   }
   mcError = mcResult *mCerelat;
}//GetIntegMC

//____________________________________________________________________________________
void  TFoam::GetIntNorm(Double_t& IntNorm, Double_t& Errel )
{
// User subprogram.
// It returns NORMALIZATION integral to be combined with the average weights
// and content of the histograms in order to get proper absolute normalization
// of the integrand and distributions.
// It can be called after initialization, before or during the MC run.

   if(fOptRej == 1) {    // Wt=1 events, internal rejection
      Double_t intMC,errMC;
      GetIntegMC(intMC,errMC);
      IntNorm = intMC;
      Errel   = errMC;
   } else {                // Wted events, NO internal rejection
      IntNorm = fPrime;
      Errel   = 0;
   }
}//GetIntNorm

//______________________________________________________________________________________
void  TFoam::GetWtParams(Double_t eps, Double_t &aveWt, Double_t &wtMax, Double_t &sigma)
{
// May be called optionally after the MC run.
// Returns various parameters of the MC weight for efficiency evaluation

   Double_t mCeff, wtLim;
   fMCMonit->GetMCeff(eps, mCeff, wtLim);
   wtMax = wtLim;
   aveWt = fSumWt/fNevGen;
   sigma = sqrt( fSumWt2/fNevGen -aveWt*aveWt );
}//GetmCeff

//_______________________________________________________________________________________
void TFoam::Finalize(Double_t& IntNorm, Double_t& Errel)
{
// May be called optionally by the user after the MC run.
// It provides normalization and also prints some information/statistics on the MC run.

   GetIntNorm(IntNorm,Errel);
   Double_t mcResult,mcError;
   GetIntegMC(mcResult,mcError);
   Double_t mCerelat= mcError/mcResult;
   //
   if(fChat>0) {
      Double_t eps = 0.0005;
      Double_t mCeff, mcEf2, wtMax, aveWt, sigma;
      GetWtParams(eps, aveWt, wtMax, sigma);
      mCeff=0;
      if(wtMax>0.0) mCeff=aveWt/wtMax;
      mcEf2 = sigma/aveWt;
      Double_t driver = fCells[0]->GetDriv();
      //
      BXOPE;
      BXTXT("****************************************");
      BXTXT("******     TFoam::Finalize       ******");
      BXTXT("****************************************");
      BX1I("    NevGen",fNevGen, "Number of generated events in the MC generation   ");
      BX1I("    nCalls",fNCalls, "Total number of function calls                    ");
      BXTXT("----------------------------------------");
      BX1F("     AveWt",aveWt,    "Average MC weight                      ");
      BX1F("     WtMin",fWtMin,  "Minimum MC weight (absolute)           ");
      BX1F("     WtMax",fWtMax,  "Maximum MC weight (absolute)           ");
      BXTXT("----------------------------------------");
      BX1F("    XPrime",fPrime,  "Primary total integral, R_prime        ");
      BX1F("    XDiver",driver,   "Driver  total integral, R_loss         ");
      BXTXT("----------------------------------------");
      BX2F("    IntMC", mcResult,  mcError,      "Result of the MC Integral");
      BX1F(" mCerelat", mCerelat,  "Relative error of the MC integral      ");
      BX1F(" <w>/WtMax",mCeff,     "MC efficiency, acceptance rate");
      BX1F(" Sigma/<w>",mcEf2,     "MC efficiency, variance/ave_wt");
      BX1F("     WtMax",wtMax,     "WtMax(esp= 0.0005)            ");
      BX1F("     Sigma",sigma,     "variance of MC weight         ");
      if(fOptRej==1) {
         Double_t avOve=fSumOve/fSumWt;
         BX1F("<OveW>/<W>",avOve,     "Contrib. of events wt>MaxWtRej");
      }
      BXCLO;
   }
}  // Finalize

//_____________________________________________________________________________________
void  TFoam::SetInhiDiv(Int_t iDim, Int_t InhiDiv)
{
// This can be called before Initialize, after setting kDim
// It defines which variables are excluded in the process of the cell division.
// For example 'FoamX->SetInhiDiv(1, 1);' inhibits division of y-variable.
// The resulting map of cells in 2-dim. case will look as follows:
//BEGIN_HTML <!--
/* -->
<img src="gif/foam_Map2.gif">
<!--*/
// -->END_HTML

   if(fDim==0) Error("TFoam","SetInhiDiv: fDim=0 \n");
   if(fInhiDiv == 0) {
      fInhiDiv = new Int_t[ fDim ];
      for(Int_t i=0; i<fDim; i++) fInhiDiv[i]=0;
   }
   //
   if( ( 0<=iDim) && (iDim<fDim)) {
      fInhiDiv[iDim] = InhiDiv;
   } else
      Error("SetInhiDiv:","Wrong iDim \n");
}//SetInhiDiv

//______________________________________________________________________________________
void  TFoam::SetXdivPRD(Int_t iDim, Int_t len, Double_t xDiv[])
{
// This should be called before Initialize, after setting  kDim
// It predefines values of the cell division for certain variable iDim.
// For example setting 3 predefined division lines using:
//     xDiv[0]=0.30; xDiv[1]=0.40; xDiv[2]=0.65;
//     FoamX->SetXdivPRD(0,3,xDiv);
// results in the following 2-dim. pattern of the cells:
//BEGIN_HTML <!--
/* -->
<img src="gif/foam_Map3.gif">
<!--*/
// -->END_HTML

   Int_t i;

   if(fDim<=0)  Error("SetXdivPRD", "fDim=0 \n");
   if(   len<1 )  Error("SetXdivPRD", "len<1 \n");
   // allocate list of pointers, if it was not done before
   if(fXdivPRD == 0) {
      fXdivPRD = new TFoamVect*[fDim];
      for(i=0; i<fDim; i++)  fXdivPRD[i]=0;
   }
  // set division list for direction iDim in H-cubic space!!!
   if( ( 0<=iDim) && (iDim<fDim)) {
      fOptPRD =1;      // !!!!
      if( fXdivPRD[iDim] != 0)
         Error("SetXdivPRD", "Second allocation of XdivPRD not allowed \n");
      fXdivPRD[iDim] = new TFoamVect(len); // allocate list of division points
      for(i=0; i<len; i++) {
         (*fXdivPRD[iDim])[i]=xDiv[i]; // set list of division points
      }
   } else {
      Error("SetXdivPRD", "Wrong iDim  \n");
   }
   // Priting predefined division points
   cout<<" SetXdivPRD, idim= "<<iDim<<"  len= "<<len<<"   "<<endl;
   for(i=0; i<len; i++) {
      cout<< (*fXdivPRD[iDim])[i] <<"  ";
   }
   cout<<endl;
   for(i=0; i<len; i++)  cout<< xDiv[i] <<"   ";
   cout<<endl;
   //
}//SetXdivPRD

//_______________________________________________________________________________________
void TFoam::CheckAll(Int_t level)
{
//  User utility, miscellaneous and debug.
//  Checks all pointers in the tree of cells. This is useful autodiagnostic.
//  level=0, no printout, failures causes STOP
//  level=1, printout, failures lead to WARNINGS only

   Int_t errors, warnings;
   TFoamCell *cell;
   Long_t iCell;

   errors = 0; warnings = 0;
   if (level==1) cout << "///////////////////////////// FOAM_Checks /////////////////////////////////" << endl;
   for(iCell=1; iCell<=fLastCe; iCell++) {
      cell = fCells[iCell];
      //  checking general rules
      if( ((cell->GetDau0()==0) && (cell->GetDau1()!=0) ) ||
         ((cell->GetDau1()==0) && (cell->GetDau0()!=0) ) ) {
         errors++;
         if (level==1) Error("CheckAll","ERROR: Cell's no %d has only one daughter \n",iCell);
      }
      if( (cell->GetDau0()==0) && (cell->GetDau1()==0) && (cell->GetStat()==0) ) {
         errors++;
         if (level==1) Error("CheckAll","ERROR: Cell's no %d  has no daughter and is inactive \n",iCell);
      }
      if( (cell->GetDau0()!=0) && (cell->GetDau1()!=0) && (cell->GetStat()==1) ) {
         errors++;
         if (level==1) Error("CheckAll","ERROR: Cell's no %d has two daughters and is active \n",iCell);
      }

      // checking parents
      if( (cell->GetPare())!=fCells[0] ) { // not child of the root
         if ( (cell != cell->GetPare()->GetDau0()) && (cell != cell->GetPare()->GetDau1()) ) {
            errors++;
            if (level==1) Error("CheckAll","ERROR: Cell's no %d parent not pointing to this cell\n ",iCell);
         }
      }

      // checking daughters
      if(cell->GetDau0()!=0) {
         if(cell != (cell->GetDau0())->GetPare()) {
            errors++;
            if (level==1)  Error("CheckAll","ERROR: Cell's no %d daughter 0 not pointing to this cell \n",iCell);
         }
      }
      if(cell->GetDau1()!=0) {
         if(cell != (cell->GetDau1())->GetPare()) {
            errors++;
            if (level==1) Error("CheckAll","ERROR: Cell's no %d daughter 1 not pointing to this cell \n",iCell);
         }
      }
   }// loop after cells;

   // Check for empty cells
   for(iCell=0; iCell<=fLastCe; iCell++) {
      cell = fCells[iCell];
      if( (cell->GetStat()==1) && (cell->GetDriv()==0) ) {
         warnings++;
         if(level==1) Warning("CheckAll", "Warning: Cell no. %d is active but empty \n", iCell);
      }
   }
   // summary
   if(level==1){
      Info("CheckAll","Check has found %d errors and %d warnings \n",errors, warnings);
   }
   if(errors>0){
      Info("CheckAll","Check - found total %d  errors \n",errors);
   }
} // Check

//________________________________________________________________________________________
void TFoam::PrintCells(void)
{
// Prints geometry of ALL cells of the FOAM

   Long_t iCell;

   for(iCell=0; iCell<=fLastCe; iCell++) {
      cout<<"Cell["<<iCell<<"]={ ";
      //cout<<"  "<< fCells[iCell]<<"  ";  // extra DEBUG
      cout<<endl;
      fCells[iCell]->Print("1");
      cout<<"}"<<endl;
   }
}

//_________________________________________________________________________________________
void TFoam::RootPlot2dim(Char_t *filename)
{
// Debugging tool which plots 2-dimensional cells as rectangles
// in C++ format readable for root

   ofstream outfile(filename, ios::out);
   Double_t   x1,y1,x2,y2,x,y;
   Long_t    iCell;
   Double_t offs =0.1;
   Double_t lpag   =1-2*offs;
   outfile<<"{" << endl;
   outfile<<"cMap = new TCanvas(\"Map1\",\" Cell Map \",600,600);"<<endl;
   //
   outfile<<"TBox*a=new TBox();"<<endl;
   outfile<<"a->SetFillStyle(0);"<<endl;  // big frame
   outfile<<"a->SetLineWidth(4);"<<endl;
   outfile<<"a->SetLineColor(2);"<<endl;
   outfile<<"a->DrawBox("<<offs<<","<<offs<<","<<(offs+lpag)<<","<<(offs+lpag)<<");"<<endl;
   //
   outfile<<"TText*t=new TText();"<<endl;  // text for numbering
   outfile<<"t->SetTextColor(4);"<<endl;
   if(fLastCe<51)
      outfile<<"t->SetTextSize(0.025);"<<endl;  // text for numbering
   else if(fLastCe<251)
      outfile<<"t->SetTextSize(0.015);"<<endl;
   else
      outfile<<"t->SetTextSize(0.008);"<<endl;
   //
   outfile<<"TBox*b=new TBox();"<<endl;  // single cell
   outfile <<"b->SetFillStyle(0);"<<endl;
   //
   if(fDim==2 && fLastCe<=2000) {
      TFoamVect  cellPosi(fDim); TFoamVect  cellSize(fDim);
      outfile << "// =========== Rectangular cells  ==========="<< endl;
      for(iCell=1; iCell<=fLastCe; iCell++) {
         if( fCells[iCell]->GetStat() == 1) {
            fCells[iCell]->GetHcub(cellPosi,cellSize);
            x1 = offs+lpag*(        cellPosi[0]); y1 = offs+lpag*(        cellPosi[1]);
            x2 = offs+lpag*(cellPosi[0]+cellSize[0]); y2 = offs+lpag*(cellPosi[1]+cellSize[1]);
            //     cell rectangle
            if(fLastCe<=2000)
            outfile<<"b->DrawBox("<<x1<<","<<y1<<","<<x2<<","<<y2<<");"<<endl;
            //     cell number
            if(fLastCe<=250) {
               x = offs+lpag*(cellPosi[0]+0.5*cellSize[0]); y = offs+lpag*(cellPosi[1]+0.5*cellSize[1]);
               outfile<<"t->DrawText("<<x<<","<<y<<","<<"\""<<iCell<<"\""<<");"<<endl;
            }
         }
      }
      outfile<<"// ============== End Rectangles ==========="<< endl;
   }//kDim=2
   //
   //
   outfile << "}" << endl;
   outfile.close();
}

void TFoam::LinkCells()
{
// Void function for backward compatibility

   Info("LinkCells", "VOID function for backward compatibility \n");
   return;
}

////////////////////////////////////////////////////////////////////////////////
//       End of Class TFoam                                                   //
////////////////////////////////////////////////////////////////////////////////

 TFoam.cxx:1
 TFoam.cxx:2
 TFoam.cxx:3
 TFoam.cxx:4
 TFoam.cxx:5
 TFoam.cxx:6
 TFoam.cxx:7
 TFoam.cxx:8
 TFoam.cxx:9
 TFoam.cxx:10
 TFoam.cxx:11
 TFoam.cxx:12
 TFoam.cxx:13
 TFoam.cxx:14
 TFoam.cxx:15
 TFoam.cxx:16
 TFoam.cxx:17
 TFoam.cxx:18
 TFoam.cxx:19
 TFoam.cxx:20
 TFoam.cxx:21
 TFoam.cxx:22
 TFoam.cxx:23
 TFoam.cxx:24
 TFoam.cxx:25
 TFoam.cxx:26
 TFoam.cxx:27
 TFoam.cxx:28
 TFoam.cxx:29
 TFoam.cxx:30
 TFoam.cxx:31
 TFoam.cxx:32
 TFoam.cxx:33
 TFoam.cxx:34
 TFoam.cxx:35
 TFoam.cxx:36
 TFoam.cxx:37
 TFoam.cxx:38
 TFoam.cxx:39
 TFoam.cxx:40
 TFoam.cxx:41
 TFoam.cxx:42
 TFoam.cxx:43
 TFoam.cxx:44
 TFoam.cxx:45
 TFoam.cxx:46
 TFoam.cxx:47
 TFoam.cxx:48
 TFoam.cxx:49
 TFoam.cxx:50
 TFoam.cxx:51
 TFoam.cxx:52
 TFoam.cxx:53
 TFoam.cxx:54
 TFoam.cxx:55
 TFoam.cxx:56
 TFoam.cxx:57
 TFoam.cxx:58
 TFoam.cxx:59
 TFoam.cxx:60
 TFoam.cxx:61
 TFoam.cxx:62
 TFoam.cxx:63
 TFoam.cxx:64
 TFoam.cxx:65
 TFoam.cxx:66
 TFoam.cxx:67
 TFoam.cxx:68
 TFoam.cxx:69
 TFoam.cxx:70
 TFoam.cxx:71
 TFoam.cxx:72
 TFoam.cxx:73
 TFoam.cxx:74
 TFoam.cxx:75
 TFoam.cxx:76
 TFoam.cxx:77
 TFoam.cxx:78
 TFoam.cxx:79
 TFoam.cxx:80
 TFoam.cxx:81
 TFoam.cxx:82
 TFoam.cxx:83
 TFoam.cxx:84
 TFoam.cxx:85
 TFoam.cxx:86
 TFoam.cxx:87
 TFoam.cxx:88
 TFoam.cxx:89
 TFoam.cxx:90
 TFoam.cxx:91
 TFoam.cxx:92
 TFoam.cxx:93
 TFoam.cxx:94
 TFoam.cxx:95
 TFoam.cxx:96
 TFoam.cxx:97
 TFoam.cxx:98
 TFoam.cxx:99
 TFoam.cxx:100
 TFoam.cxx:101
 TFoam.cxx:102
 TFoam.cxx:103
 TFoam.cxx:104
 TFoam.cxx:105
 TFoam.cxx:106
 TFoam.cxx:107
 TFoam.cxx:108
 TFoam.cxx:109
 TFoam.cxx:110
 TFoam.cxx:111
 TFoam.cxx:112
 TFoam.cxx:113
 TFoam.cxx:114
 TFoam.cxx:115
 TFoam.cxx:116
 TFoam.cxx:117
 TFoam.cxx:118
 TFoam.cxx:119
 TFoam.cxx:120
 TFoam.cxx:121
 TFoam.cxx:122
 TFoam.cxx:123
 TFoam.cxx:124
 TFoam.cxx:125
 TFoam.cxx:126
 TFoam.cxx:127
 TFoam.cxx:128
 TFoam.cxx:129
 TFoam.cxx:130
 TFoam.cxx:131
 TFoam.cxx:132
 TFoam.cxx:133
 TFoam.cxx:134
 TFoam.cxx:135
 TFoam.cxx:136
 TFoam.cxx:137
 TFoam.cxx:138
 TFoam.cxx:139
 TFoam.cxx:140
 TFoam.cxx:141
 TFoam.cxx:142
 TFoam.cxx:143
 TFoam.cxx:144
 TFoam.cxx:145
 TFoam.cxx:146
 TFoam.cxx:147
 TFoam.cxx:148
 TFoam.cxx:149
 TFoam.cxx:150
 TFoam.cxx:151
 TFoam.cxx:152
 TFoam.cxx:153
 TFoam.cxx:154
 TFoam.cxx:155
 TFoam.cxx:156
 TFoam.cxx:157
 TFoam.cxx:158
 TFoam.cxx:159
 TFoam.cxx:160
 TFoam.cxx:161
 TFoam.cxx:162
 TFoam.cxx:163
 TFoam.cxx:164
 TFoam.cxx:165
 TFoam.cxx:166
 TFoam.cxx:167
 TFoam.cxx:168
 TFoam.cxx:169
 TFoam.cxx:170
 TFoam.cxx:171
 TFoam.cxx:172
 TFoam.cxx:173
 TFoam.cxx:174
 TFoam.cxx:175
 TFoam.cxx:176
 TFoam.cxx:177
 TFoam.cxx:178
 TFoam.cxx:179
 TFoam.cxx:180
 TFoam.cxx:181
 TFoam.cxx:182
 TFoam.cxx:183
 TFoam.cxx:184
 TFoam.cxx:185
 TFoam.cxx:186
 TFoam.cxx:187
 TFoam.cxx:188
 TFoam.cxx:189
 TFoam.cxx:190
 TFoam.cxx:191
 TFoam.cxx:192
 TFoam.cxx:193
 TFoam.cxx:194
 TFoam.cxx:195
 TFoam.cxx:196
 TFoam.cxx:197
 TFoam.cxx:198
 TFoam.cxx:199
 TFoam.cxx:200
 TFoam.cxx:201
 TFoam.cxx:202
 TFoam.cxx:203
 TFoam.cxx:204
 TFoam.cxx:205
 TFoam.cxx:206
 TFoam.cxx:207
 TFoam.cxx:208
 TFoam.cxx:209
 TFoam.cxx:210
 TFoam.cxx:211
 TFoam.cxx:212
 TFoam.cxx:213
 TFoam.cxx:214
 TFoam.cxx:215
 TFoam.cxx:216
 TFoam.cxx:217
 TFoam.cxx:218
 TFoam.cxx:219
 TFoam.cxx:220
 TFoam.cxx:221
 TFoam.cxx:222
 TFoam.cxx:223
 TFoam.cxx:224
 TFoam.cxx:225
 TFoam.cxx:226
 TFoam.cxx:227
 TFoam.cxx:228
 TFoam.cxx:229
 TFoam.cxx:230
 TFoam.cxx:231
 TFoam.cxx:232
 TFoam.cxx:233
 TFoam.cxx:234
 TFoam.cxx:235
 TFoam.cxx:236
 TFoam.cxx:237
 TFoam.cxx:238
 TFoam.cxx:239
 TFoam.cxx:240
 TFoam.cxx:241
 TFoam.cxx:242
 TFoam.cxx:243
 TFoam.cxx:244
 TFoam.cxx:245
 TFoam.cxx:246
 TFoam.cxx:247
 TFoam.cxx:248
 TFoam.cxx:249
 TFoam.cxx:250
 TFoam.cxx:251
 TFoam.cxx:252
 TFoam.cxx:253
 TFoam.cxx:254
 TFoam.cxx:255
 TFoam.cxx:256
 TFoam.cxx:257
 TFoam.cxx:258
 TFoam.cxx:259
 TFoam.cxx:260
 TFoam.cxx:261
 TFoam.cxx:262
 TFoam.cxx:263
 TFoam.cxx:264
 TFoam.cxx:265
 TFoam.cxx:266
 TFoam.cxx:267
 TFoam.cxx:268
 TFoam.cxx:269
 TFoam.cxx:270
 TFoam.cxx:271
 TFoam.cxx:272
 TFoam.cxx:273
 TFoam.cxx:274
 TFoam.cxx:275
 TFoam.cxx:276
 TFoam.cxx:277
 TFoam.cxx:278
 TFoam.cxx:279
 TFoam.cxx:280
 TFoam.cxx:281
 TFoam.cxx:282
 TFoam.cxx:283
 TFoam.cxx:284
 TFoam.cxx:285
 TFoam.cxx:286
 TFoam.cxx:287
 TFoam.cxx:288
 TFoam.cxx:289
 TFoam.cxx:290
 TFoam.cxx:291
 TFoam.cxx:292
 TFoam.cxx:293
 TFoam.cxx:294
 TFoam.cxx:295
 TFoam.cxx:296
 TFoam.cxx:297
 TFoam.cxx:298
 TFoam.cxx:299
 TFoam.cxx:300
 TFoam.cxx:301
 TFoam.cxx:302
 TFoam.cxx:303
 TFoam.cxx:304
 TFoam.cxx:305
 TFoam.cxx:306
 TFoam.cxx:307
 TFoam.cxx:308
 TFoam.cxx:309
 TFoam.cxx:310
 TFoam.cxx:311
 TFoam.cxx:312
 TFoam.cxx:313
 TFoam.cxx:314
 TFoam.cxx:315
 TFoam.cxx:316
 TFoam.cxx:317
 TFoam.cxx:318
 TFoam.cxx:319
 TFoam.cxx:320
 TFoam.cxx:321
 TFoam.cxx:322
 TFoam.cxx:323
 TFoam.cxx:324
 TFoam.cxx:325
 TFoam.cxx:326
 TFoam.cxx:327
 TFoam.cxx:328
 TFoam.cxx:329
 TFoam.cxx:330
 TFoam.cxx:331
 TFoam.cxx:332
 TFoam.cxx:333
 TFoam.cxx:334
 TFoam.cxx:335
 TFoam.cxx:336
 TFoam.cxx:337
 TFoam.cxx:338
 TFoam.cxx:339
 TFoam.cxx:340
 TFoam.cxx:341
 TFoam.cxx:342
 TFoam.cxx:343
 TFoam.cxx:344
 TFoam.cxx:345
 TFoam.cxx:346
 TFoam.cxx:347
 TFoam.cxx:348
 TFoam.cxx:349
 TFoam.cxx:350
 TFoam.cxx:351
 TFoam.cxx:352
 TFoam.cxx:353
 TFoam.cxx:354
 TFoam.cxx:355
 TFoam.cxx:356
 TFoam.cxx:357
 TFoam.cxx:358
 TFoam.cxx:359
 TFoam.cxx:360
 TFoam.cxx:361
 TFoam.cxx:362
 TFoam.cxx:363
 TFoam.cxx:364
 TFoam.cxx:365
 TFoam.cxx:366
 TFoam.cxx:367
 TFoam.cxx:368
 TFoam.cxx:369
 TFoam.cxx:370
 TFoam.cxx:371
 TFoam.cxx:372
 TFoam.cxx:373
 TFoam.cxx:374
 TFoam.cxx:375
 TFoam.cxx:376
 TFoam.cxx:377
 TFoam.cxx:378
 TFoam.cxx:379
 TFoam.cxx:380
 TFoam.cxx:381
 TFoam.cxx:382
 TFoam.cxx:383
 TFoam.cxx:384
 TFoam.cxx:385
 TFoam.cxx:386
 TFoam.cxx:387
 TFoam.cxx:388
 TFoam.cxx:389
 TFoam.cxx:390
 TFoam.cxx:391
 TFoam.cxx:392
 TFoam.cxx:393
 TFoam.cxx:394
 TFoam.cxx:395
 TFoam.cxx:396
 TFoam.cxx:397
 TFoam.cxx:398
 TFoam.cxx:399
 TFoam.cxx:400
 TFoam.cxx:401
 TFoam.cxx:402
 TFoam.cxx:403
 TFoam.cxx:404
 TFoam.cxx:405
 TFoam.cxx:406
 TFoam.cxx:407
 TFoam.cxx:408
 TFoam.cxx:409
 TFoam.cxx:410
 TFoam.cxx:411
 TFoam.cxx:412
 TFoam.cxx:413
 TFoam.cxx:414
 TFoam.cxx:415
 TFoam.cxx:416
 TFoam.cxx:417
 TFoam.cxx:418
 TFoam.cxx:419
 TFoam.cxx:420
 TFoam.cxx:421
 TFoam.cxx:422
 TFoam.cxx:423
 TFoam.cxx:424
 TFoam.cxx:425
 TFoam.cxx:426
 TFoam.cxx:427
 TFoam.cxx:428
 TFoam.cxx:429
 TFoam.cxx:430
 TFoam.cxx:431
 TFoam.cxx:432
 TFoam.cxx:433
 TFoam.cxx:434
 TFoam.cxx:435
 TFoam.cxx:436
 TFoam.cxx:437
 TFoam.cxx:438
 TFoam.cxx:439
 TFoam.cxx:440
 TFoam.cxx:441
 TFoam.cxx:442
 TFoam.cxx:443
 TFoam.cxx:444
 TFoam.cxx:445
 TFoam.cxx:446
 TFoam.cxx:447
 TFoam.cxx:448
 TFoam.cxx:449
 TFoam.cxx:450
 TFoam.cxx:451
 TFoam.cxx:452
 TFoam.cxx:453
 TFoam.cxx:454
 TFoam.cxx:455
 TFoam.cxx:456
 TFoam.cxx:457
 TFoam.cxx:458
 TFoam.cxx:459
 TFoam.cxx:460
 TFoam.cxx:461
 TFoam.cxx:462
 TFoam.cxx:463
 TFoam.cxx:464
 TFoam.cxx:465
 TFoam.cxx:466
 TFoam.cxx:467
 TFoam.cxx:468
 TFoam.cxx:469
 TFoam.cxx:470
 TFoam.cxx:471
 TFoam.cxx:472
 TFoam.cxx:473
 TFoam.cxx:474
 TFoam.cxx:475
 TFoam.cxx:476
 TFoam.cxx:477
 TFoam.cxx:478
 TFoam.cxx:479
 TFoam.cxx:480
 TFoam.cxx:481
 TFoam.cxx:482
 TFoam.cxx:483
 TFoam.cxx:484
 TFoam.cxx:485
 TFoam.cxx:486
 TFoam.cxx:487
 TFoam.cxx:488
 TFoam.cxx:489
 TFoam.cxx:490
 TFoam.cxx:491
 TFoam.cxx:492
 TFoam.cxx:493
 TFoam.cxx:494
 TFoam.cxx:495
 TFoam.cxx:496
 TFoam.cxx:497
 TFoam.cxx:498
 TFoam.cxx:499
 TFoam.cxx:500
 TFoam.cxx:501
 TFoam.cxx:502
 TFoam.cxx:503
 TFoam.cxx:504
 TFoam.cxx:505
 TFoam.cxx:506
 TFoam.cxx:507
 TFoam.cxx:508
 TFoam.cxx:509
 TFoam.cxx:510
 TFoam.cxx:511
 TFoam.cxx:512
 TFoam.cxx:513
 TFoam.cxx:514
 TFoam.cxx:515
 TFoam.cxx:516
 TFoam.cxx:517
 TFoam.cxx:518
 TFoam.cxx:519
 TFoam.cxx:520
 TFoam.cxx:521
 TFoam.cxx:522
 TFoam.cxx:523
 TFoam.cxx:524
 TFoam.cxx:525
 TFoam.cxx:526
 TFoam.cxx:527
 TFoam.cxx:528
 TFoam.cxx:529
 TFoam.cxx:530
 TFoam.cxx:531
 TFoam.cxx:532
 TFoam.cxx:533
 TFoam.cxx:534
 TFoam.cxx:535
 TFoam.cxx:536
 TFoam.cxx:537
 TFoam.cxx:538
 TFoam.cxx:539
 TFoam.cxx:540
 TFoam.cxx:541
 TFoam.cxx:542
 TFoam.cxx:543
 TFoam.cxx:544
 TFoam.cxx:545
 TFoam.cxx:546
 TFoam.cxx:547
 TFoam.cxx:548
 TFoam.cxx:549
 TFoam.cxx:550
 TFoam.cxx:551
 TFoam.cxx:552
 TFoam.cxx:553
 TFoam.cxx:554
 TFoam.cxx:555
 TFoam.cxx:556
 TFoam.cxx:557
 TFoam.cxx:558
 TFoam.cxx:559
 TFoam.cxx:560
 TFoam.cxx:561
 TFoam.cxx:562
 TFoam.cxx:563
 TFoam.cxx:564
 TFoam.cxx:565
 TFoam.cxx:566
 TFoam.cxx:567
 TFoam.cxx:568
 TFoam.cxx:569
 TFoam.cxx:570
 TFoam.cxx:571
 TFoam.cxx:572
 TFoam.cxx:573
 TFoam.cxx:574
 TFoam.cxx:575
 TFoam.cxx:576
 TFoam.cxx:577
 TFoam.cxx:578
 TFoam.cxx:579
 TFoam.cxx:580
 TFoam.cxx:581
 TFoam.cxx:582
 TFoam.cxx:583
 TFoam.cxx:584
 TFoam.cxx:585
 TFoam.cxx:586
 TFoam.cxx:587
 TFoam.cxx:588
 TFoam.cxx:589
 TFoam.cxx:590
 TFoam.cxx:591
 TFoam.cxx:592
 TFoam.cxx:593
 TFoam.cxx:594
 TFoam.cxx:595
 TFoam.cxx:596
 TFoam.cxx:597
 TFoam.cxx:598
 TFoam.cxx:599
 TFoam.cxx:600
 TFoam.cxx:601
 TFoam.cxx:602
 TFoam.cxx:603
 TFoam.cxx:604
 TFoam.cxx:605
 TFoam.cxx:606
 TFoam.cxx:607
 TFoam.cxx:608
 TFoam.cxx:609
 TFoam.cxx:610
 TFoam.cxx:611
 TFoam.cxx:612
 TFoam.cxx:613
 TFoam.cxx:614
 TFoam.cxx:615
 TFoam.cxx:616
 TFoam.cxx:617
 TFoam.cxx:618
 TFoam.cxx:619
 TFoam.cxx:620
 TFoam.cxx:621
 TFoam.cxx:622
 TFoam.cxx:623
 TFoam.cxx:624
 TFoam.cxx:625
 TFoam.cxx:626
 TFoam.cxx:627
 TFoam.cxx:628
 TFoam.cxx:629
 TFoam.cxx:630
 TFoam.cxx:631
 TFoam.cxx:632
 TFoam.cxx:633
 TFoam.cxx:634
 TFoam.cxx:635
 TFoam.cxx:636
 TFoam.cxx:637
 TFoam.cxx:638
 TFoam.cxx:639
 TFoam.cxx:640
 TFoam.cxx:641
 TFoam.cxx:642
 TFoam.cxx:643
 TFoam.cxx:644
 TFoam.cxx:645
 TFoam.cxx:646
 TFoam.cxx:647
 TFoam.cxx:648
 TFoam.cxx:649
 TFoam.cxx:650
 TFoam.cxx:651
 TFoam.cxx:652
 TFoam.cxx:653
 TFoam.cxx:654
 TFoam.cxx:655
 TFoam.cxx:656
 TFoam.cxx:657
 TFoam.cxx:658
 TFoam.cxx:659
 TFoam.cxx:660
 TFoam.cxx:661
 TFoam.cxx:662
 TFoam.cxx:663
 TFoam.cxx:664
 TFoam.cxx:665
 TFoam.cxx:666
 TFoam.cxx:667
 TFoam.cxx:668
 TFoam.cxx:669
 TFoam.cxx:670
 TFoam.cxx:671
 TFoam.cxx:672
 TFoam.cxx:673
 TFoam.cxx:674
 TFoam.cxx:675
 TFoam.cxx:676
 TFoam.cxx:677
 TFoam.cxx:678
 TFoam.cxx:679
 TFoam.cxx:680
 TFoam.cxx:681
 TFoam.cxx:682
 TFoam.cxx:683
 TFoam.cxx:684
 TFoam.cxx:685
 TFoam.cxx:686
 TFoam.cxx:687
 TFoam.cxx:688
 TFoam.cxx:689
 TFoam.cxx:690
 TFoam.cxx:691
 TFoam.cxx:692
 TFoam.cxx:693
 TFoam.cxx:694
 TFoam.cxx:695
 TFoam.cxx:696
 TFoam.cxx:697
 TFoam.cxx:698
 TFoam.cxx:699
 TFoam.cxx:700
 TFoam.cxx:701
 TFoam.cxx:702
 TFoam.cxx:703
 TFoam.cxx:704
 TFoam.cxx:705
 TFoam.cxx:706
 TFoam.cxx:707
 TFoam.cxx:708
 TFoam.cxx:709
 TFoam.cxx:710
 TFoam.cxx:711
 TFoam.cxx:712
 TFoam.cxx:713
 TFoam.cxx:714
 TFoam.cxx:715
 TFoam.cxx:716
 TFoam.cxx:717
 TFoam.cxx:718
 TFoam.cxx:719
 TFoam.cxx:720
 TFoam.cxx:721
 TFoam.cxx:722
 TFoam.cxx:723
 TFoam.cxx:724
 TFoam.cxx:725
 TFoam.cxx:726
 TFoam.cxx:727
 TFoam.cxx:728
 TFoam.cxx:729
 TFoam.cxx:730
 TFoam.cxx:731
 TFoam.cxx:732
 TFoam.cxx:733
 TFoam.cxx:734
 TFoam.cxx:735
 TFoam.cxx:736
 TFoam.cxx:737
 TFoam.cxx:738
 TFoam.cxx:739
 TFoam.cxx:740
 TFoam.cxx:741
 TFoam.cxx:742
 TFoam.cxx:743
 TFoam.cxx:744
 TFoam.cxx:745
 TFoam.cxx:746
 TFoam.cxx:747
 TFoam.cxx:748
 TFoam.cxx:749
 TFoam.cxx:750
 TFoam.cxx:751
 TFoam.cxx:752
 TFoam.cxx:753
 TFoam.cxx:754
 TFoam.cxx:755
 TFoam.cxx:756
 TFoam.cxx:757
 TFoam.cxx:758
 TFoam.cxx:759
 TFoam.cxx:760
 TFoam.cxx:761
 TFoam.cxx:762
 TFoam.cxx:763
 TFoam.cxx:764
 TFoam.cxx:765
 TFoam.cxx:766
 TFoam.cxx:767
 TFoam.cxx:768
 TFoam.cxx:769
 TFoam.cxx:770
 TFoam.cxx:771
 TFoam.cxx:772
 TFoam.cxx:773
 TFoam.cxx:774
 TFoam.cxx:775
 TFoam.cxx:776
 TFoam.cxx:777
 TFoam.cxx:778
 TFoam.cxx:779
 TFoam.cxx:780
 TFoam.cxx:781
 TFoam.cxx:782
 TFoam.cxx:783
 TFoam.cxx:784
 TFoam.cxx:785
 TFoam.cxx:786
 TFoam.cxx:787
 TFoam.cxx:788
 TFoam.cxx:789
 TFoam.cxx:790
 TFoam.cxx:791
 TFoam.cxx:792
 TFoam.cxx:793
 TFoam.cxx:794
 TFoam.cxx:795
 TFoam.cxx:796
 TFoam.cxx:797
 TFoam.cxx:798
 TFoam.cxx:799
 TFoam.cxx:800
 TFoam.cxx:801
 TFoam.cxx:802
 TFoam.cxx:803
 TFoam.cxx:804
 TFoam.cxx:805
 TFoam.cxx:806
 TFoam.cxx:807
 TFoam.cxx:808
 TFoam.cxx:809
 TFoam.cxx:810
 TFoam.cxx:811
 TFoam.cxx:812
 TFoam.cxx:813
 TFoam.cxx:814
 TFoam.cxx:815
 TFoam.cxx:816
 TFoam.cxx:817
 TFoam.cxx:818
 TFoam.cxx:819
 TFoam.cxx:820
 TFoam.cxx:821
 TFoam.cxx:822
 TFoam.cxx:823
 TFoam.cxx:824
 TFoam.cxx:825
 TFoam.cxx:826
 TFoam.cxx:827
 TFoam.cxx:828
 TFoam.cxx:829
 TFoam.cxx:830
 TFoam.cxx:831
 TFoam.cxx:832
 TFoam.cxx:833
 TFoam.cxx:834
 TFoam.cxx:835
 TFoam.cxx:836
 TFoam.cxx:837
 TFoam.cxx:838
 TFoam.cxx:839
 TFoam.cxx:840
 TFoam.cxx:841
 TFoam.cxx:842
 TFoam.cxx:843
 TFoam.cxx:844
 TFoam.cxx:845
 TFoam.cxx:846
 TFoam.cxx:847
 TFoam.cxx:848
 TFoam.cxx:849
 TFoam.cxx:850
 TFoam.cxx:851
 TFoam.cxx:852
 TFoam.cxx:853
 TFoam.cxx:854
 TFoam.cxx:855
 TFoam.cxx:856
 TFoam.cxx:857
 TFoam.cxx:858
 TFoam.cxx:859
 TFoam.cxx:860
 TFoam.cxx:861
 TFoam.cxx:862
 TFoam.cxx:863
 TFoam.cxx:864
 TFoam.cxx:865
 TFoam.cxx:866
 TFoam.cxx:867
 TFoam.cxx:868
 TFoam.cxx:869
 TFoam.cxx:870
 TFoam.cxx:871
 TFoam.cxx:872
 TFoam.cxx:873
 TFoam.cxx:874
 TFoam.cxx:875
 TFoam.cxx:876
 TFoam.cxx:877
 TFoam.cxx:878
 TFoam.cxx:879
 TFoam.cxx:880
 TFoam.cxx:881
 TFoam.cxx:882
 TFoam.cxx:883
 TFoam.cxx:884
 TFoam.cxx:885
 TFoam.cxx:886
 TFoam.cxx:887
 TFoam.cxx:888
 TFoam.cxx:889
 TFoam.cxx:890
 TFoam.cxx:891
 TFoam.cxx:892
 TFoam.cxx:893
 TFoam.cxx:894
 TFoam.cxx:895
 TFoam.cxx:896
 TFoam.cxx:897
 TFoam.cxx:898
 TFoam.cxx:899
 TFoam.cxx:900
 TFoam.cxx:901
 TFoam.cxx:902
 TFoam.cxx:903
 TFoam.cxx:904
 TFoam.cxx:905
 TFoam.cxx:906
 TFoam.cxx:907
 TFoam.cxx:908
 TFoam.cxx:909
 TFoam.cxx:910
 TFoam.cxx:911
 TFoam.cxx:912
 TFoam.cxx:913
 TFoam.cxx:914
 TFoam.cxx:915
 TFoam.cxx:916
 TFoam.cxx:917
 TFoam.cxx:918
 TFoam.cxx:919
 TFoam.cxx:920
 TFoam.cxx:921
 TFoam.cxx:922
 TFoam.cxx:923
 TFoam.cxx:924
 TFoam.cxx:925
 TFoam.cxx:926
 TFoam.cxx:927
 TFoam.cxx:928
 TFoam.cxx:929
 TFoam.cxx:930
 TFoam.cxx:931
 TFoam.cxx:932
 TFoam.cxx:933
 TFoam.cxx:934
 TFoam.cxx:935
 TFoam.cxx:936
 TFoam.cxx:937
 TFoam.cxx:938
 TFoam.cxx:939
 TFoam.cxx:940
 TFoam.cxx:941
 TFoam.cxx:942
 TFoam.cxx:943
 TFoam.cxx:944
 TFoam.cxx:945
 TFoam.cxx:946
 TFoam.cxx:947
 TFoam.cxx:948
 TFoam.cxx:949
 TFoam.cxx:950
 TFoam.cxx:951
 TFoam.cxx:952
 TFoam.cxx:953
 TFoam.cxx:954
 TFoam.cxx:955
 TFoam.cxx:956
 TFoam.cxx:957
 TFoam.cxx:958
 TFoam.cxx:959
 TFoam.cxx:960
 TFoam.cxx:961
 TFoam.cxx:962
 TFoam.cxx:963
 TFoam.cxx:964
 TFoam.cxx:965
 TFoam.cxx:966
 TFoam.cxx:967
 TFoam.cxx:968
 TFoam.cxx:969
 TFoam.cxx:970
 TFoam.cxx:971
 TFoam.cxx:972
 TFoam.cxx:973
 TFoam.cxx:974
 TFoam.cxx:975
 TFoam.cxx:976
 TFoam.cxx:977
 TFoam.cxx:978
 TFoam.cxx:979
 TFoam.cxx:980
 TFoam.cxx:981
 TFoam.cxx:982
 TFoam.cxx:983
 TFoam.cxx:984
 TFoam.cxx:985
 TFoam.cxx:986
 TFoam.cxx:987
 TFoam.cxx:988
 TFoam.cxx:989
 TFoam.cxx:990
 TFoam.cxx:991
 TFoam.cxx:992
 TFoam.cxx:993
 TFoam.cxx:994
 TFoam.cxx:995
 TFoam.cxx:996
 TFoam.cxx:997
 TFoam.cxx:998
 TFoam.cxx:999
 TFoam.cxx:1000
 TFoam.cxx:1001
 TFoam.cxx:1002
 TFoam.cxx:1003
 TFoam.cxx:1004
 TFoam.cxx:1005
 TFoam.cxx:1006
 TFoam.cxx:1007
 TFoam.cxx:1008
 TFoam.cxx:1009
 TFoam.cxx:1010
 TFoam.cxx:1011
 TFoam.cxx:1012
 TFoam.cxx:1013
 TFoam.cxx:1014
 TFoam.cxx:1015
 TFoam.cxx:1016
 TFoam.cxx:1017
 TFoam.cxx:1018
 TFoam.cxx:1019
 TFoam.cxx:1020
 TFoam.cxx:1021
 TFoam.cxx:1022
 TFoam.cxx:1023
 TFoam.cxx:1024
 TFoam.cxx:1025
 TFoam.cxx:1026
 TFoam.cxx:1027
 TFoam.cxx:1028
 TFoam.cxx:1029
 TFoam.cxx:1030
 TFoam.cxx:1031
 TFoam.cxx:1032
 TFoam.cxx:1033
 TFoam.cxx:1034
 TFoam.cxx:1035
 TFoam.cxx:1036
 TFoam.cxx:1037
 TFoam.cxx:1038
 TFoam.cxx:1039
 TFoam.cxx:1040
 TFoam.cxx:1041
 TFoam.cxx:1042
 TFoam.cxx:1043
 TFoam.cxx:1044
 TFoam.cxx:1045
 TFoam.cxx:1046
 TFoam.cxx:1047
 TFoam.cxx:1048
 TFoam.cxx:1049
 TFoam.cxx:1050
 TFoam.cxx:1051
 TFoam.cxx:1052
 TFoam.cxx:1053
 TFoam.cxx:1054
 TFoam.cxx:1055
 TFoam.cxx:1056
 TFoam.cxx:1057
 TFoam.cxx:1058
 TFoam.cxx:1059
 TFoam.cxx:1060
 TFoam.cxx:1061
 TFoam.cxx:1062
 TFoam.cxx:1063
 TFoam.cxx:1064
 TFoam.cxx:1065
 TFoam.cxx:1066
 TFoam.cxx:1067
 TFoam.cxx:1068
 TFoam.cxx:1069
 TFoam.cxx:1070
 TFoam.cxx:1071
 TFoam.cxx:1072
 TFoam.cxx:1073
 TFoam.cxx:1074
 TFoam.cxx:1075
 TFoam.cxx:1076
 TFoam.cxx:1077
 TFoam.cxx:1078
 TFoam.cxx:1079
 TFoam.cxx:1080
 TFoam.cxx:1081
 TFoam.cxx:1082
 TFoam.cxx:1083
 TFoam.cxx:1084
 TFoam.cxx:1085
 TFoam.cxx:1086
 TFoam.cxx:1087
 TFoam.cxx:1088
 TFoam.cxx:1089
 TFoam.cxx:1090
 TFoam.cxx:1091
 TFoam.cxx:1092
 TFoam.cxx:1093
 TFoam.cxx:1094
 TFoam.cxx:1095
 TFoam.cxx:1096
 TFoam.cxx:1097
 TFoam.cxx:1098
 TFoam.cxx:1099
 TFoam.cxx:1100
 TFoam.cxx:1101
 TFoam.cxx:1102
 TFoam.cxx:1103
 TFoam.cxx:1104
 TFoam.cxx:1105
 TFoam.cxx:1106
 TFoam.cxx:1107
 TFoam.cxx:1108
 TFoam.cxx:1109
 TFoam.cxx:1110
 TFoam.cxx:1111
 TFoam.cxx:1112
 TFoam.cxx:1113
 TFoam.cxx:1114
 TFoam.cxx:1115
 TFoam.cxx:1116
 TFoam.cxx:1117
 TFoam.cxx:1118
 TFoam.cxx:1119
 TFoam.cxx:1120
 TFoam.cxx:1121
 TFoam.cxx:1122
 TFoam.cxx:1123
 TFoam.cxx:1124
 TFoam.cxx:1125
 TFoam.cxx:1126
 TFoam.cxx:1127
 TFoam.cxx:1128
 TFoam.cxx:1129
 TFoam.cxx:1130
 TFoam.cxx:1131
 TFoam.cxx:1132
 TFoam.cxx:1133
 TFoam.cxx:1134
 TFoam.cxx:1135
 TFoam.cxx:1136
 TFoam.cxx:1137
 TFoam.cxx:1138
 TFoam.cxx:1139
 TFoam.cxx:1140
 TFoam.cxx:1141
 TFoam.cxx:1142
 TFoam.cxx:1143
 TFoam.cxx:1144
 TFoam.cxx:1145
 TFoam.cxx:1146
 TFoam.cxx:1147
 TFoam.cxx:1148
 TFoam.cxx:1149
 TFoam.cxx:1150
 TFoam.cxx:1151
 TFoam.cxx:1152
 TFoam.cxx:1153
 TFoam.cxx:1154
 TFoam.cxx:1155
 TFoam.cxx:1156
 TFoam.cxx:1157
 TFoam.cxx:1158
 TFoam.cxx:1159
 TFoam.cxx:1160
 TFoam.cxx:1161
 TFoam.cxx:1162
 TFoam.cxx:1163
 TFoam.cxx:1164
 TFoam.cxx:1165
 TFoam.cxx:1166
 TFoam.cxx:1167
 TFoam.cxx:1168
 TFoam.cxx:1169
 TFoam.cxx:1170
 TFoam.cxx:1171
 TFoam.cxx:1172
 TFoam.cxx:1173
 TFoam.cxx:1174
 TFoam.cxx:1175
 TFoam.cxx:1176
 TFoam.cxx:1177
 TFoam.cxx:1178
 TFoam.cxx:1179
 TFoam.cxx:1180
 TFoam.cxx:1181
 TFoam.cxx:1182
 TFoam.cxx:1183
 TFoam.cxx:1184
 TFoam.cxx:1185
 TFoam.cxx:1186
 TFoam.cxx:1187
 TFoam.cxx:1188
 TFoam.cxx:1189
 TFoam.cxx:1190
 TFoam.cxx:1191
 TFoam.cxx:1192
 TFoam.cxx:1193
 TFoam.cxx:1194
 TFoam.cxx:1195
 TFoam.cxx:1196
 TFoam.cxx:1197
 TFoam.cxx:1198
 TFoam.cxx:1199
 TFoam.cxx:1200
 TFoam.cxx:1201
 TFoam.cxx:1202
 TFoam.cxx:1203
 TFoam.cxx:1204
 TFoam.cxx:1205
 TFoam.cxx:1206
 TFoam.cxx:1207
 TFoam.cxx:1208
 TFoam.cxx:1209
 TFoam.cxx:1210
 TFoam.cxx:1211
 TFoam.cxx:1212
 TFoam.cxx:1213
 TFoam.cxx:1214
 TFoam.cxx:1215
 TFoam.cxx:1216
 TFoam.cxx:1217
 TFoam.cxx:1218
 TFoam.cxx:1219
 TFoam.cxx:1220
 TFoam.cxx:1221
 TFoam.cxx:1222
 TFoam.cxx:1223
 TFoam.cxx:1224
 TFoam.cxx:1225
 TFoam.cxx:1226
 TFoam.cxx:1227
 TFoam.cxx:1228
 TFoam.cxx:1229
 TFoam.cxx:1230
 TFoam.cxx:1231
 TFoam.cxx:1232
 TFoam.cxx:1233
 TFoam.cxx:1234
 TFoam.cxx:1235
 TFoam.cxx:1236
 TFoam.cxx:1237
 TFoam.cxx:1238
 TFoam.cxx:1239
 TFoam.cxx:1240
 TFoam.cxx:1241
 TFoam.cxx:1242
 TFoam.cxx:1243
 TFoam.cxx:1244
 TFoam.cxx:1245
 TFoam.cxx:1246
 TFoam.cxx:1247
 TFoam.cxx:1248
 TFoam.cxx:1249
 TFoam.cxx:1250
 TFoam.cxx:1251
 TFoam.cxx:1252
 TFoam.cxx:1253
 TFoam.cxx:1254
 TFoam.cxx:1255
 TFoam.cxx:1256
 TFoam.cxx:1257
 TFoam.cxx:1258
 TFoam.cxx:1259
 TFoam.cxx:1260
 TFoam.cxx:1261
 TFoam.cxx:1262
 TFoam.cxx:1263
 TFoam.cxx:1264
 TFoam.cxx:1265
 TFoam.cxx:1266
 TFoam.cxx:1267
 TFoam.cxx:1268
 TFoam.cxx:1269
 TFoam.cxx:1270
 TFoam.cxx:1271
 TFoam.cxx:1272
 TFoam.cxx:1273
 TFoam.cxx:1274
 TFoam.cxx:1275
 TFoam.cxx:1276
 TFoam.cxx:1277
 TFoam.cxx:1278
 TFoam.cxx:1279
 TFoam.cxx:1280
 TFoam.cxx:1281
 TFoam.cxx:1282
 TFoam.cxx:1283
 TFoam.cxx:1284
 TFoam.cxx:1285
 TFoam.cxx:1286
 TFoam.cxx:1287
 TFoam.cxx:1288
 TFoam.cxx:1289
 TFoam.cxx:1290
 TFoam.cxx:1291
 TFoam.cxx:1292
 TFoam.cxx:1293
 TFoam.cxx:1294
 TFoam.cxx:1295
 TFoam.cxx:1296
 TFoam.cxx:1297
 TFoam.cxx:1298
 TFoam.cxx:1299
 TFoam.cxx:1300
 TFoam.cxx:1301
 TFoam.cxx:1302
 TFoam.cxx:1303
 TFoam.cxx:1304
 TFoam.cxx:1305
 TFoam.cxx:1306
 TFoam.cxx:1307
 TFoam.cxx:1308
 TFoam.cxx:1309
 TFoam.cxx:1310
 TFoam.cxx:1311
 TFoam.cxx:1312
 TFoam.cxx:1313
 TFoam.cxx:1314
 TFoam.cxx:1315
 TFoam.cxx:1316
 TFoam.cxx:1317
 TFoam.cxx:1318
 TFoam.cxx:1319
 TFoam.cxx:1320
 TFoam.cxx:1321
 TFoam.cxx:1322
 TFoam.cxx:1323
 TFoam.cxx:1324
 TFoam.cxx:1325
 TFoam.cxx:1326
 TFoam.cxx:1327
 TFoam.cxx:1328
 TFoam.cxx:1329
 TFoam.cxx:1330
 TFoam.cxx:1331
 TFoam.cxx:1332
 TFoam.cxx:1333
 TFoam.cxx:1334
 TFoam.cxx:1335
 TFoam.cxx:1336
 TFoam.cxx:1337
 TFoam.cxx:1338
 TFoam.cxx:1339
 TFoam.cxx:1340
 TFoam.cxx:1341
 TFoam.cxx:1342
 TFoam.cxx:1343
 TFoam.cxx:1344
 TFoam.cxx:1345
 TFoam.cxx:1346
 TFoam.cxx:1347
 TFoam.cxx:1348
 TFoam.cxx:1349
 TFoam.cxx:1350
 TFoam.cxx:1351
 TFoam.cxx:1352
 TFoam.cxx:1353
 TFoam.cxx:1354
 TFoam.cxx:1355
 TFoam.cxx:1356
 TFoam.cxx:1357
 TFoam.cxx:1358
 TFoam.cxx:1359
 TFoam.cxx:1360
 TFoam.cxx:1361
 TFoam.cxx:1362
 TFoam.cxx:1363
 TFoam.cxx:1364
 TFoam.cxx:1365
 TFoam.cxx:1366
 TFoam.cxx:1367
 TFoam.cxx:1368
 TFoam.cxx:1369
 TFoam.cxx:1370
 TFoam.cxx:1371
 TFoam.cxx:1372
 TFoam.cxx:1373
 TFoam.cxx:1374
 TFoam.cxx:1375
 TFoam.cxx:1376
 TFoam.cxx:1377
 TFoam.cxx:1378
 TFoam.cxx:1379
 TFoam.cxx:1380
 TFoam.cxx:1381
 TFoam.cxx:1382
 TFoam.cxx:1383
 TFoam.cxx:1384
 TFoam.cxx:1385
 TFoam.cxx:1386
 TFoam.cxx:1387
 TFoam.cxx:1388
 TFoam.cxx:1389
 TFoam.cxx:1390
 TFoam.cxx:1391
 TFoam.cxx:1392
 TFoam.cxx:1393
 TFoam.cxx:1394
 TFoam.cxx:1395
 TFoam.cxx:1396
 TFoam.cxx:1397
 TFoam.cxx:1398
 TFoam.cxx:1399
 TFoam.cxx:1400
 TFoam.cxx:1401
 TFoam.cxx:1402
 TFoam.cxx:1403
 TFoam.cxx:1404
 TFoam.cxx:1405
 TFoam.cxx:1406
 TFoam.cxx:1407
 TFoam.cxx:1408
 TFoam.cxx:1409
 TFoam.cxx:1410
 TFoam.cxx:1411
 TFoam.cxx:1412
 TFoam.cxx:1413
 TFoam.cxx:1414
 TFoam.cxx:1415
 TFoam.cxx:1416
 TFoam.cxx:1417
 TFoam.cxx:1418
 TFoam.cxx:1419
 TFoam.cxx:1420
 TFoam.cxx:1421
 TFoam.cxx:1422
 TFoam.cxx:1423
 TFoam.cxx:1424
 TFoam.cxx:1425
 TFoam.cxx:1426
 TFoam.cxx:1427
 TFoam.cxx:1428
 TFoam.cxx:1429
 TFoam.cxx:1430
 TFoam.cxx:1431
 TFoam.cxx:1432
 TFoam.cxx:1433
 TFoam.cxx:1434
 TFoam.cxx:1435
 TFoam.cxx:1436
 TFoam.cxx:1437
 TFoam.cxx:1438
 TFoam.cxx:1439
 TFoam.cxx:1440
 TFoam.cxx:1441
 TFoam.cxx:1442
 TFoam.cxx:1443
 TFoam.cxx:1444
 TFoam.cxx:1445
 TFoam.cxx:1446
 TFoam.cxx:1447
 TFoam.cxx:1448
 TFoam.cxx:1449
 TFoam.cxx:1450
 TFoam.cxx:1451
 TFoam.cxx:1452
 TFoam.cxx:1453
 TFoam.cxx:1454
 TFoam.cxx:1455
 TFoam.cxx:1456
 TFoam.cxx:1457
 TFoam.cxx:1458
 TFoam.cxx:1459
 TFoam.cxx:1460
 TFoam.cxx:1461
 TFoam.cxx:1462
 TFoam.cxx:1463
 TFoam.cxx:1464
 TFoam.cxx:1465
 TFoam.cxx:1466
 TFoam.cxx:1467
 TFoam.cxx:1468
 TFoam.cxx:1469
 TFoam.cxx:1470
 TFoam.cxx:1471
 TFoam.cxx:1472
 TFoam.cxx:1473
 TFoam.cxx:1474
 TFoam.cxx:1475
 TFoam.cxx:1476
 TFoam.cxx:1477
 TFoam.cxx:1478
 TFoam.cxx:1479
 TFoam.cxx:1480
 TFoam.cxx:1481
 TFoam.cxx:1482
 TFoam.cxx:1483
 TFoam.cxx:1484
 TFoam.cxx:1485
 TFoam.cxx:1486
 TFoam.cxx:1487
 TFoam.cxx:1488
 TFoam.cxx:1489
 TFoam.cxx:1490
 TFoam.cxx:1491
 TFoam.cxx:1492
 TFoam.cxx:1493
 TFoam.cxx:1494
 TFoam.cxx:1495
 TFoam.cxx:1496
 TFoam.cxx:1497
 TFoam.cxx:1498
 TFoam.cxx:1499
 TFoam.cxx:1500
 TFoam.cxx:1501
 TFoam.cxx:1502
 TFoam.cxx:1503
 TFoam.cxx:1504
 TFoam.cxx:1505
 TFoam.cxx:1506
 TFoam.cxx:1507
 TFoam.cxx:1508
 TFoam.cxx:1509
 TFoam.cxx:1510
 TFoam.cxx:1511
 TFoam.cxx:1512
 TFoam.cxx:1513
 TFoam.cxx:1514
 TFoam.cxx:1515
 TFoam.cxx:1516
 TFoam.cxx:1517
 TFoam.cxx:1518
 TFoam.cxx:1519
 TFoam.cxx:1520
 TFoam.cxx:1521
 TFoam.cxx:1522
 TFoam.cxx:1523
 TFoam.cxx:1524
 TFoam.cxx:1525
 TFoam.cxx:1526
 TFoam.cxx:1527
 TFoam.cxx:1528
 TFoam.cxx:1529
 TFoam.cxx:1530
 TFoam.cxx:1531
 TFoam.cxx:1532
 TFoam.cxx:1533