TFoam.cxx

Go to the documentation of this file.
00001 // @(#)root/foam:$Id: TFoam.cxx 36037 2010-10-02 08:06:04Z moneta $
00002 // Author: S. Jadach <mailto:Stanislaw.jadach@ifj.edu.pl>, P.Sawicki <mailto:Pawel.Sawicki@ifj.edu.pl>
00003 
00004 //______________________________________________________________________________
00005 //
00006 // FOAM  Version 1.02M
00007 // ===================
00008 // Authors:
00009 //   S. Jadach and P.Sawicki
00010 //   Institute of Nuclear Physics, Cracow, Poland
00011 //   Stanislaw. Jadach@ifj.edu.pl, Pawel.Sawicki@ifj.edu.pl
00012 //
00013 // What is FOAM for?
00014 // =================
00015 // * Suppose you want to generate randomly points (vectors) according to
00016 //   an arbitrary probability distribution  in n dimensions,
00017 //   for which you supply your own subprogram. FOAM can do it for you!
00018 //   Even if your distributions has quite strong peaks and is discontinuous!
00019 // * FOAM generates random points with weight one or with variable weight.
00020 // * FOAM is capable to integrate using efficient "adaptive" MC method.
00021 //   (The distribution does not need to be normalized to one.)
00022 // How does it work?
00023 // =================
00024 // FOAM is the simplified version of the multi-dimensional general purpose
00025 // Monte Carlo event generator (integrator) FOAM.
00026 // It creates hyper-rectangular "foam of cells", which is more dense around its peaks.
00027 // See the following 2-dim. example of the map of 1000 cells for doubly peaked distribution:
00028 //BEGIN_HTML <!--
00029 /* -->
00030 <img src="gif/foam_MapCamel1000.gif">
00031 <!--*/
00032 // -->END_HTML
00033 // FOAM is now fully integrated with the ROOT package.
00034 // The important bonus of the ROOT use is persistency of the FOAM objects!
00035 //
00036 // For more sophisticated problems full version of FOAM may be more appropriate:
00037 //BEGIN_HTML <!--
00038 /* -->
00039   See <A HREF="http://jadach.home.cern.ch/jadach/Foam/Index.html"> full version of FOAM</A>
00040 <!--*/
00041 // -->END_HTML
00042 // Simple example of the use of FOAM:
00043 // ==================================
00044 // Int_t kanwa(){
00045 //   gSystem->Load("libFoam");
00046 //   TH2D  *hst_xy = new TH2D("hst_xy" ,  "x-y plot", 50,0,1.0, 50,0,1.0);
00047 //   Double_t *MCvect =new Double_t[2]; // 2-dim vector generated in the MC run
00048 //   TRandom3  *PseRan   = new TRandom3();  // Create random number generator
00049 //   PseRan->SetSeed(4357);                // Set seed
00050 //   TFoam   *FoamX    = new TFoam("FoamX");   // Create Simulator
00051 //   FoamX->SetkDim(2);          // No. of dimensions, obligatory!
00052 //   FoamX->SetnCells(500);      // No. of cells, can be omitted, default=2000
00053 //   FoamX->SetRhoInt(Camel2);   // Set 2-dim distribution, included below
00054 //   FoamX->SetPseRan(PseRan);   // Set random number generator
00055 //   FoamX->Initialize();        // Initialize simulator, takes a few seconds...
00056 //   // From now on FoamX is ready to generate events according to Camel2(x,y)
00057 //   for(Long_t loop=0; loop<100000; loop++){
00058 //     FoamX->MakeEvent();          // generate MC event
00059 //     FoamX->GetMCvect( MCvect);   // get generated vector (x,y)
00060 //     Double_t x=MCvect[0];
00061 //     Double_t y=MCvect[1];
00062 //     if(loop<10) cout<<"(x,y) =  ( "<< x <<", "<< y <<" )"<<endl;
00063 //     hst_xy->Fill(x,y);           // fill scattergram
00064 //   }// loop
00065 //   Double_t mcResult, mcError;
00066 //   FoamX->GetIntegMC( mcResult, mcError);  // get MC integral, should be one
00067 //   cout << " mcResult= " << mcResult << " +- " << mcError <<endl;
00068 //   // now hst_xy will be plotted visualizing generated distribution
00069 //   TCanvas *cKanwa = new TCanvas("cKanwa","Canvas for plotting",600,600);
00070 //   cKanwa->cd();
00071 //   hst_xy->Draw("lego2");
00072 // }//kanwa
00073 // Double_t sqr(Double_t x){return x*x;};
00074 // Double_t Camel2(Int_t nDim, Double_t *Xarg){
00075 // // 2-dimensional distribution for FOAM, normalized to one (within 1e-5)
00076 //   Double_t x=Xarg[0];
00077 //   Double_t y=Xarg[1];
00078 //   Double_t GamSq= sqr(0.100e0);
00079 //   Double_t Dist=exp(-(sqr(x-1./3) +sqr(y-1./3))/GamSq)/GamSq/TMath::Pi();
00080 //   Dist        +=exp(-(sqr(x-2./3) +sqr(y-2./3))/GamSq)/GamSq/TMath::Pi();
00081 //   return 0.5*Dist;
00082 // }// Camel2
00083 // Two-dim. histogram of the MC points generated with the above program looks as follows:
00084 //BEGIN_HTML <!--
00085 /* -->
00086 <img src="gif/foam_cKanwa.gif">
00087 <!--*/
00088 // -->END_HTML
00089 // Canonical nine steering parameters of FOAM
00090 // ===========================================
00091 //------------------------------------------------------------------------------
00092 //  Name     | default  | Description
00093 //------------------------------------------------------------------------------
00094 //  kDim     | 0        | Dimension of the integration space. Must be redefined!
00095 //  nCells   | 1000     | No of allocated number of cells,
00096 //  nSampl   | 200      | No. of MC events in the cell MC exploration
00097 //  nBin     | 8        | No. of bins in edge-histogram in cell exploration
00098 //  OptRej   | 1        | OptRej = 0, weighted; OptRej=1, wt=1 MC events
00099 //  OptDrive | 2        | Maximum weight reduction, =1 for variance reduction
00100 //  EvPerBin | 25       | Maximum number of the effective wt=1 events/bin,
00101 //           |          | EvPerBin=0 deactivates this option
00102 //  Chat     | 1        | =0,1,2 is the ``chat level'' in the standard output
00103 //  MaxWtRej | 1.1      | Maximum weight used to get w=1 MC events
00104 //------------------------------------------------------------------------------
00105 // The above can be redefined before calling 'Initialize()' method,
00106 // for instance FoamObject->SetkDim(15) sets dimension of the distribution to 15.
00107 // Only kDim HAS TO BE redefined, the other parameters may be left at their defaults.
00108 // nCell may be increased up to about million cells for wildly peaked distributions.
00109 // Increasing nSampl sometimes helps, but it may cost CPU time.
00110 // MaxWtRej may need to be increased for wild a distribution, while using OptRej=0.
00111 //
00112 // --------------------------------------------------------------------
00113 // Past versions of FOAM: August 2003, v.1.00; September 2003 v.1.01
00114 // Adopted starting from FOAM-2.06 by P. Sawicki
00115 // --------------------------------------------------------------------
00116 // Users of FOAM are kindly requested to cite the following work:
00117 // S. Jadach, Computer Physics Communications 152 (2003) 55.
00118 //
00119 //______________________________________________________________________________
00120 
00121 #include "TFoam.h"
00122 #include "TFoamIntegrand.h"
00123 #include "TFoamMaxwt.h"
00124 #include "TFoamVect.h"
00125 #include "TFoamCell.h"
00126 #include "Riostream.h"
00127 #include "TH1.h"
00128 #include "TRefArray.h"
00129 #include "TMethodCall.h"
00130 #include "TRandom.h"
00131 #include "TMath.h"
00132 #include "TInterpreter.h"
00133 
00134 ClassImp(TFoam);
00135 
00136 //FFFFFF  BoX-FORMATs for nice and flexible outputs
00137 #define BXOPE cout<<\
00138 "FFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFF"<<endl<<\
00139 "F                                                                              F"<<endl
00140 #define BXTXT(text) cout<<\
00141 "F                   "<<setw(40)<<         text           <<"                   F"<<endl
00142 #define BX1I(name,numb,text) cout<<\
00143 "F "<<setw(10)<<name<<" = "<<setw(10)<<numb<<" = "          <<setw(50)<<text<<" F"<<endl
00144 #define BX1F(name,numb,text)     cout<<"F "<<setw(10)<<name<<\
00145           " = "<<setw(15)<<setprecision(8)<<numb<<"   =    "<<setw(40)<<text<<" F"<<endl
00146 #define BX2F(name,numb,err,text) cout<<"F "<<setw(10)<<name<<\
00147 " = "<<setw(15)<<setprecision(8)<<numb<<" +- "<<setw(15)<<setprecision(8)<<err<<\
00148                                                       "  = "<<setw(25)<<text<<" F"<<endl
00149 #define BXCLO cout<<\
00150 "F                                                                              F"<<endl<<\
00151 "FFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFF"<<endl
00152   //FFFFFF  BoX-FORMATs ends here
00153 
00154 static const Double_t gHigh= 1.0e150;
00155 static const Double_t gVlow=-1.0e150;
00156 
00157 #define SW2 setprecision(7) << setw(12)
00158 
00159 //________________________________________________________________________________________________
00160 TFoam::TFoam() : 
00161    fDim(0), fNCells(0), fRNmax(0), 
00162    fOptDrive(0), fChat(0), fOptRej(0), 
00163    fNBin(0), fNSampl(0), fEvPerBin(0), 
00164    fMaskDiv(0), fInhiDiv(0), fOptPRD(0), fXdivPRD(0), 
00165    fNoAct(0), fLastCe(0), fCells(0), 
00166    fMCMonit(0), fMaxWtRej(0), fCellsAct(0), fPrimAcu(0), 
00167    fHistEdg(0), fHistDbg(0), fHistWt(0), 
00168    fMCvect(0), fMCwt(0), fRvec(0), 
00169    fRho(0), fMethodCall(0), fPseRan(0), 
00170    fNCalls(0), fNEffev(0), 
00171    fSumWt(0), fSumWt2(0), 
00172    fSumOve(0), fNevGen(0), 
00173    fWtMax(0), fWtMin(0), 
00174    fPrime(0), fMCresult(0), fMCerror(0), 
00175    fAlpha(0)
00176 {
00177   // Default constructor for streamer, user should not use it.
00178 }
00179 //_________________________________________________________________________________________________
00180 TFoam::TFoam(const Char_t* Name) :
00181    fDim(0), fNCells(0), fRNmax(0), 
00182    fOptDrive(0), fChat(0), fOptRej(0), 
00183    fNBin(0), fNSampl(0), fEvPerBin(0), 
00184    fMaskDiv(0), fInhiDiv(0), fOptPRD(0), fXdivPRD(0), 
00185    fNoAct(0), fLastCe(0), fCells(0), 
00186    fMCMonit(0), fMaxWtRej(0), fCellsAct(0), fPrimAcu(0), 
00187    fHistEdg(0), fHistDbg(0), fHistWt(0), 
00188    fMCvect(0), fMCwt(0), fRvec(0), 
00189    fRho(0), fMethodCall(0), fPseRan(0), 
00190    fNCalls(0), fNEffev(0), 
00191    fSumWt(0), fSumWt2(0), 
00192    fSumOve(0), fNevGen(0), 
00193    fWtMax(0), fWtMin(0), 
00194    fPrime(0), fMCresult(0), fMCerror(0), 
00195    fAlpha(0)
00196 {
00197 // User constructor, to be employed by the user
00198 
00199    if(strlen(Name)  >129) {
00200       Error("TFoam","Name too long %s \n",Name);
00201    }
00202    fName=Name;                                            // Class name
00203    fDate="  Release date:  2005.04.10";                   // Release date
00204    fVersion= "1.02M";                                     // Release version
00205    fMaskDiv  = 0;             // Dynamic Mask for  cell division, h-cubic
00206    fInhiDiv  = 0;             // Flag allowing to inhibit cell division in certain projection/edge
00207    fXdivPRD  = 0;             // Lists of division values encoded in one vector per direction
00208    fCells    = 0;
00209    fAlpha    = 0;
00210    fCellsAct = 0;
00211    fPrimAcu  = 0;
00212    fHistEdg  = 0;
00213    fHistWt   = 0;
00214    fHistDbg  = 0;
00215    fDim     = 0;                // dimension of hyp-cubical space
00216    fNCells   = 1000;             // Maximum number of Cells,    is usually re-set
00217    fNSampl   = 200;              // No of sampling when dividing cell
00218    fOptPRD   = 0;                // General Option switch for PRedefined Division, for quick check
00219    fOptDrive = 2;                // type of Drive =1,2 for TrueVol,Sigma,WtMax
00220    fChat     = 1;                // Chat=0,1,2 chat level in output, Chat=1 normal level
00221    fOptRej   = 1;                // OptRej=0, wted events; OptRej=1, wt=1 events
00222    //------------------------------------------------------
00223    fNBin     = 8;                // binning of edge-histogram in cell exploration
00224    fEvPerBin =25;                // maximum no. of EFFECTIVE event per bin, =0 option is inactive
00225    //------------------------------------------------------
00226    fNCalls = 0;                  // No of function calls
00227    fNEffev = 0;                  // Total no of eff. wt=1 events in build=up
00228    fLastCe =-1;                  // Index of the last cell
00229    fNoAct  = 0;                  // No of active cells (used in MC generation)
00230    fWtMin = gHigh;               // Minimal weight
00231    fWtMax = gVlow;               // Maximal weight
00232    fMaxWtRej =1.10;              // Maximum weight in rejection for getting wt=1 events
00233    fPseRan   = 0;                // Initialize private copy of random number generator
00234    fMCMonit  = 0;                // MC efficiency monitoring
00235    fRho = 0;                     // pointer to abstract class providing function to integrate
00236    fMCvect   = 0;
00237    fRvec     = 0;
00238    fPseRan   = 0;                // generator of pseudorandom numbers
00239    fMethodCall=0;                // ROOT's pointer to global distribution function
00240 }
00241 
00242 //_______________________________________________________________________________________________
00243 TFoam::~TFoam()
00244 {
00245 // Default destructor
00246 //  cout<<" DESTRUCTOR entered "<<endl;
00247    Int_t i;
00248 
00249    if(fCells!= 0) {
00250       for(i=0; i<fNCells; i++) delete fCells[i]; // TFoamCell*[]
00251       delete [] fCells;
00252    }
00253    if (fCellsAct) delete fCellsAct ; // WVE FIX LEAK
00254    if (fRvec)    delete [] fRvec;    //double[]
00255    if (fAlpha)   delete [] fAlpha;   //double[]
00256    if (fMCvect)  delete [] fMCvect;  //double[]
00257    if (fPrimAcu) delete [] fPrimAcu; //double[]
00258    if (fMaskDiv) delete [] fMaskDiv; //int[]
00259    if (fInhiDiv) delete [] fInhiDiv; //int[]
00260  
00261    if( fXdivPRD!= 0) {
00262       for(i=0; i<fDim; i++) delete fXdivPRD[i]; // TFoamVect*[]
00263       delete [] fXdivPRD;
00264    }
00265    if (fMCMonit) delete fMCMonit;
00266    if (fHistWt)  delete fHistWt;
00267 
00268    // delete histogram arrays
00269    if (fHistEdg) { 
00270       fHistEdg->Delete(); 
00271       delete fHistEdg; 
00272    }
00273    if (fHistDbg) { 
00274       fHistDbg->Delete(); 
00275       delete fHistDbg; 
00276    }
00277 }
00278 
00279 
00280 //_____________________________________________________________________________________________
00281 TFoam::TFoam(const TFoam &From): TObject(From)
00282 {
00283 // Copy Constructor  NOT IMPLEMENTED (NEVER USED)
00284    Error("TFoam", "COPY CONSTRUCTOR NOT IMPLEMENTED \n");
00285 }
00286 
00287 //_____________________________________________________________________________________________
00288 void TFoam::Initialize(TRandom *PseRan, TFoamIntegrand *fun )
00289 {
00290 // Basic initialization of FOAM invoked by the user. Mandatory!
00291 // ============================================================
00292 // This method starts the process of the cell build-up.
00293 // User must invoke Initialize with two arguments or Initialize without arguments.
00294 // This is done BEFORE generating first MC event and AFTER allocating FOAM object
00295 // and reseting (optionally) its internal parameters/switches.
00296 // The overall operational scheme of the FOAM is the following:
00297 //BEGIN_HTML <!--
00298 /* -->
00299 <img src="gif/foam_schema2.gif">
00300 <!--*/
00301 // -->END_HTML
00302 //
00303 // This method invokes several other methods:
00304 // ==========================================
00305 // InitCells initializes memory storage for cells and begins exploration process
00306 // from the root cell. The empty cells are allocated/filled using  CellFill.
00307 // The procedure Grow which loops over cells, picks up the cell with the biggest
00308 // ``driver integral'', see Comp. Phys. Commun. 152 152 (2003) 55 for explanations,
00309 // with the help of PeekMax procedure. The chosen cell is split using Divide.
00310 // Subsequently, the procedure Explore called by the Divide
00311 // (and by InitCells for the root cell) does the most important
00312 // job in the FOAM object build-up: it performs a small MC run for each
00313 // newly allocated daughter cell.
00314 // Explore calculates how profitable the future split of the cell will be
00315 // and defines the optimal cell division geometry with the help of Carver or Varedu
00316 // procedures, for maximum weight or variance optimization respectively.
00317 // All essential results of the exploration are written into
00318 // the explored cell object. At the very end of the foam build-up,
00319 // Finally, MakeActiveList is invoked to create a list of pointers to
00320 // all active cells, for the purpose of the quick access during the MC generation.
00321 // The procedure Explore employs MakeAlpha to generate random coordinates
00322 // inside a given cell with the uniform distribution.
00323 // The above sequence of the procedure calls is depicted in the following figure:
00324 //BEGIN_HTML <!--
00325 /* -->
00326 <img src="gif/foam_Initialize_schema.gif">
00327 <!--*/
00328 // -->END_HTML
00329 
00330    SetPseRan(PseRan);
00331    SetRho(fun);
00332    Initialize();
00333 }
00334 
00335 //_______________________________________________________________________________________
00336 void TFoam::Initialize()
00337 {
00338 // Basic initialization of FOAM invoked by the user.
00339 // IMPORTANT: Random number generator and the distribution object has to be
00340 // provided using SetPseRan and SetRho prior to invoking this initializator!
00341 
00342    Bool_t addStatus = TH1::AddDirectoryStatus();
00343    TH1::AddDirectory(kFALSE);
00344    Int_t i;
00345 
00346    if(fChat>0){
00347       BXOPE;
00348       BXTXT("****************************************");
00349       BXTXT("******      TFoam::Initialize    ******");
00350       BXTXT("****************************************");
00351       BXTXT(fName);
00352       BX1F("  Version",fVersion,  fDate);
00353       BX1I("     kDim",fDim,     " Dimension of the hyper-cubical space             ");
00354       BX1I("   nCells",fNCells,   " Requested number of Cells (half of them active)  ");
00355       BX1I("   nSampl",fNSampl,   " No of MC events in exploration of a cell         ");
00356       BX1I("     nBin",fNBin,     " No of bins in histograms, MC exploration of cell ");
00357       BX1I(" EvPerBin",fEvPerBin, " Maximum No effective_events/bin, MC exploration  ");
00358       BX1I(" OptDrive",fOptDrive, " Type of Driver   =1,2 for Sigma,WtMax            ");
00359       BX1I("   OptRej",fOptRej,   " MC rejection on/off for OptRej=0,1               ");
00360       BX1F(" MaxWtRej",fMaxWtRej, " Maximum wt in rejection for wt=1 evts");
00361       BXCLO;
00362    }
00363 
00364    if(fPseRan==0) Error("Initialize", "Random number generator not set \n");
00365    if(fRho==0 && fMethodCall==0 ) Error("Initialize", "Distribution function not set \n");
00366    if(fDim==0) Error("Initialize", "Zero dimension not allowed \n");
00367 
00368    /////////////////////////////////////////////////////////////////////////
00369    //                   ALLOCATE SMALL LISTS                              //
00370    //  it is done globally, not for each cell, to save on allocation time //
00371    /////////////////////////////////////////////////////////////////////////
00372    fRNmax= fDim+1;
00373    fRvec = new Double_t[fRNmax];   // Vector of random numbers
00374    if(fRvec==0)  Error("Initialize", "Cannot initialize buffer fRvec \n");
00375 
00376    if(fDim>0){
00377       fAlpha = new Double_t[fDim];    // sum<1 for internal parametrization of the simplex
00378       if(fAlpha==0)  Error("Initialize", "Cannot initialize buffer fAlpha \n" );
00379    }
00380    fMCvect = new Double_t[fDim]; // vector generated in the MC run
00381    if(fMCvect==0)  Error("Initialize", "Cannot initialize buffer fMCvect  \n" );
00382 
00383    //====== List of directions inhibited for division
00384    if(fInhiDiv == 0){
00385       fInhiDiv = new Int_t[fDim];
00386       for(i=0; i<fDim; i++) fInhiDiv[i]=0;
00387    }
00388    //====== Dynamic mask used in Explore for edge determination
00389    if(fMaskDiv == 0){
00390       fMaskDiv = new Int_t[fDim];
00391       for(i=0; i<fDim; i++) fMaskDiv[i]=1;
00392    }
00393    //====== List of predefined division values in all directions (initialized as empty)
00394    if(fXdivPRD == 0){
00395       fXdivPRD = new TFoamVect*[fDim];
00396       for(i=0; i<fDim; i++)  fXdivPRD[i]=0; // Artificially extended beyond fDim
00397    }
00398    //====== Initialize list of histograms
00399    fHistWt  = new TH1D("HistWt","Histogram of MC weight",100,0.0, 1.5*fMaxWtRej); // MC weight
00400    fHistEdg = new TObjArray(fDim);           // Initialize list of histograms
00401    TString hname;
00402    TString htitle;
00403    for(i=0;i<fDim;i++){
00404       hname=fName+TString("_HistEdge_");
00405       hname+=i;
00406       htitle=TString("Edge Histogram No. ");
00407       htitle+=i;
00408       //cout<<"i= "<<i<<"  hname= "<<hname<<"  htitle= "<<htitle<<endl;
00409       (*fHistEdg)[i] = new TH1D(hname.Data(),htitle.Data(),fNBin,0.0, 1.0); // Initialize histogram for each edge
00410       ((TH1D*)(*fHistEdg)[i])->Sumw2();
00411    }
00412    //======  extra histograms for debug purposes
00413    fHistDbg = new TObjArray(fDim);         // Initialize list of histograms
00414    for(i=0;i<fDim;i++){
00415       hname=fName+TString("_HistDebug_");
00416       hname+=i;
00417       htitle=TString("Debug Histogram ");
00418       htitle+=i;
00419       (*fHistDbg)[i] = new TH1D(hname.Data(),htitle.Data(),fNBin,0.0, 1.0); // Initialize histogram for each edge
00420    }
00421 
00422    // ||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||| //
00423    //                     BUILD-UP of the FOAM                            //
00424    // ||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||| //
00425    //
00426    //        Define and explore root cell(s)
00427    InitCells();
00428    //        PrintCells(); cout<<" ===== after InitCells ====="<<endl;
00429    Grow();
00430    //        PrintCells(); cout<<" ===== after Grow      ====="<<endl;
00431 
00432    MakeActiveList(); // Final Preperations for the M.C. generation
00433 
00434    // Preperations for the M.C. generation
00435    fSumWt  = 0.0;               // M.C. generation sum of Wt
00436    fSumWt2 = 0.0;               // M.C. generation sum of Wt**2
00437    fSumOve = 0.0;               // M.C. generation sum of overweighted
00438    fNevGen = 0.0;               // M.C. generation sum of 1d0
00439    fWtMax  = gVlow;               // M.C. generation maximum wt
00440    fWtMin  = gHigh;               // M.C. generation minimum wt
00441    fMCresult=fCells[0]->GetIntg(); // M.C. Value of INTEGRAL,temporary assignment
00442    fMCresult=fCells[0]->GetIntg(); // M.C. Value of INTEGRAL,temporary assignment
00443    fMCerror =fCells[0]->GetIntg(); // M.C. Value of ERROR   ,temporary assignment
00444    fMCMonit = new TFoamMaxwt(5.0,1000);  // monitoring M.C. efficiency
00445    //
00446    if(fChat>0){
00447       Double_t driver = fCells[0]->GetDriv();
00448       BXOPE;
00449       BXTXT("***  TFoam::Initialize FINISHED!!!  ***");
00450       BX1I("    nCalls",fNCalls,  "Total number of function calls         ");
00451       BX1F("    XPrime",fPrime,   "Primary total integral                 ");
00452       BX1F("    XDiver",driver,    "Driver  total integral                 ");
00453       BX1F("  mcResult",fMCresult,"Estimate of the true MC Integral       ");
00454       BXCLO;
00455    }
00456    if(fChat==2) PrintCells();
00457    TH1::AddDirectory(addStatus);
00458 } // Initialize
00459 
00460 //_______________________________________________________________________________________
00461 void TFoam::InitCells()
00462 {
00463 // Internal subprogram used by Initialize.
00464 // It initializes "root part" of the FOAM of the tree of cells.
00465 
00466    Int_t i;
00467 
00468    fLastCe =-1;                             // Index of the last cell
00469    if(fCells!= 0) {
00470       for(i=0; i<fNCells; i++) delete fCells[i];
00471       delete [] fCells;
00472    }
00473    //
00474    fCells = new TFoamCell*[fNCells];
00475    for(i=0;i<fNCells;i++){
00476       fCells[i]= new TFoamCell(fDim); // Allocate BIG list of cells
00477       fCells[i]->SetSerial(i);
00478    }
00479    if(fCells==0) Error("InitCells", "Cannot initialize CELLS \n"  );
00480 
00481    /////////////////////////////////////////////////////////////////////////////
00482    //              Single Root Hypercube                                      //
00483    /////////////////////////////////////////////////////////////////////////////
00484    CellFill(1,   0);  //  0-th cell ACTIVE
00485 
00486    // Exploration of the root cell(s)
00487    for(Long_t iCell=0; iCell<=fLastCe; iCell++){
00488       Explore( fCells[iCell] );               // Exploration of root cell(s)
00489    }
00490 }//InitCells
00491 
00492 //_______________________________________________________________________________________
00493 Int_t TFoam::CellFill(Int_t Status, TFoamCell *parent)
00494 {
00495 // Internal subprogram used by Initialize.
00496 // It initializes content of the newly allocated active cell.
00497 
00498    TFoamCell *cell;
00499    if (fLastCe==fNCells){
00500       Error( "CellFill", "Too many cells\n");
00501    }
00502    fLastCe++;   // 0-th cell is the first
00503    if (Status==1) fNoAct++;
00504 
00505    cell = fCells[fLastCe];
00506 
00507    cell->Fill(Status, parent, 0, 0);
00508 
00509    cell->SetBest( -1);         // pointer for planning division of the cell
00510    cell->SetXdiv(0.5);         // factor for division
00511    Double_t xInt2,xDri2;
00512    if(parent!=0){
00513       xInt2  = 0.5*parent->GetIntg();
00514       xDri2  = 0.5*parent->GetDriv();
00515       cell->SetIntg(xInt2);
00516       cell->SetDriv(xDri2);
00517    }else{
00518       cell->SetIntg(0.0);
00519       cell->SetDriv(0.0);
00520    }
00521    return fLastCe;
00522 }
00523 
00524 //______________________________________________________________________________________
00525 void TFoam::Explore(TFoamCell *cell)
00526 {
00527 // Internal subprogram used by Initialize.
00528 // It explores newly defined cell with help of special short MC sampling.
00529 // As a result, estimates of true and drive volume is defined/determined
00530 // Average and dispersion of the weight distribution will is found along
00531 // each edge and the best edge (minimum dispersion, best maximum weight)
00532 // is memorized for future use.
00533 // The optimal division point for eventual future cell division is
00534 // determined/recorded. Recorded are also minimum and maximum weight etc.
00535 // The volume estimate in all (inactive) parent cells is updated.
00536 // Note that links to parents and initial volume = 1/2 parent has to be
00537 // already defined prior to calling this routine.
00538 
00539    Double_t wt, dx, xBest=0, yBest=0;
00540    Double_t intOld, driOld;
00541 
00542    Long_t iev;
00543    Double_t nevMC;
00544    Int_t i, j, k;
00545    Int_t nProj, kBest;
00546    Double_t ceSum[5], xproj;
00547 
00548    TFoamVect  cellSize(fDim);
00549    TFoamVect  cellPosi(fDim);
00550 
00551    cell->GetHcub(cellPosi,cellSize);
00552 
00553    TFoamCell  *parent;
00554 
00555    Double_t *xRand = new Double_t[fDim];
00556 
00557    Double_t *volPart=0;
00558 
00559    cell->CalcVolume();
00560    dx = cell->GetVolume();
00561    intOld = cell->GetIntg(); //memorize old values,
00562    driOld = cell->GetDriv(); //will be needed for correcting parent cells
00563 
00564 
00565    /////////////////////////////////////////////////////
00566    //    Special Short MC sampling to probe cell      //
00567    /////////////////////////////////////////////////////
00568    ceSum[0]=0;
00569    ceSum[1]=0;
00570    ceSum[2]=0;
00571    ceSum[3]=gHigh;  //wtmin
00572    ceSum[4]=gVlow;  //wtmax
00573    //
00574    for(i=0;i<fDim;i++) ((TH1D *)(*fHistEdg)[i])->Reset(); // Reset histograms
00575    fHistWt->Reset();
00576    //
00577    // ||||||||||||||||||||||||||BEGIN MC LOOP|||||||||||||||||||||||||||||
00578    Double_t nevEff=0.;
00579    for(iev=0;iev<fNSampl;iev++){
00580       MakeAlpha();               // generate uniformly vector inside hypercube
00581 
00582       if(fDim>0){
00583       for(j=0; j<fDim; j++)
00584          xRand[j]= cellPosi[j] +fAlpha[j]*(cellSize[j]);
00585       }
00586 
00587       wt=dx*Eval(xRand);
00588 
00589       nProj = 0;
00590       if(fDim>0) {
00591          for(k=0; k<fDim; k++) {
00592             xproj =fAlpha[k];
00593             ((TH1D *)(*fHistEdg)[nProj])->Fill(xproj,wt);
00594             nProj++;
00595          }
00596       }
00597       //
00598       fNCalls++;
00599       ceSum[0] += wt;    // sum of weights
00600       ceSum[1] += wt*wt; // sum of weights squared
00601       ceSum[2]++;        // sum of 1
00602       if (ceSum[3]>wt) ceSum[3]=wt;  // minimum weight;
00603       if (ceSum[4]<wt) ceSum[4]=wt;  // maximum weight
00604       // test MC loop exit condition
00605       nevEff = ceSum[0]*ceSum[0]/ceSum[1];
00606       if( nevEff >= fNBin*fEvPerBin) break;
00607    }   // ||||||||||||||||||||||||||END MC LOOP|||||||||||||||||||||||||||||
00608    //------------------------------------------------------------------
00609    //---  predefine logics of searching for the best division edge ---
00610    for(k=0; k<fDim;k++){
00611       fMaskDiv[k] =1;                       // default is all
00612       if( fInhiDiv[k]==1) fMaskDiv[k] =0; // inhibit some...
00613    }
00614    // Note that predefined division below overrule inhibition above
00615    kBest=-1;
00616    Double_t rmin,rmax,rdiv;
00617    if(fOptPRD) {          // quick check
00618       for(k=0; k<fDim; k++) {
00619          rmin= cellPosi[k];
00620          rmax= cellPosi[k] +cellSize[k];
00621          if( fXdivPRD[k] != 0) {
00622             Int_t n= (fXdivPRD[k])->GetDim();
00623             for(j=0; j<n; j++) {
00624                rdiv=(*fXdivPRD[k])[j];
00625                // check predefined divisions is available in this cell
00626                if( (rmin +1e-99 <rdiv) && (rdiv< rmax -1e-99)) {
00627                   kBest=k;
00628                   xBest= (rdiv-cellPosi[k])/cellSize[k] ;
00629                   goto ee05;
00630                }
00631             }
00632          }
00633       }//k
00634    }
00635    ee05:
00636    //------------------------------------------------------------------
00637    fNEffev += (Long_t)nevEff;
00638    nevMC          = ceSum[2];
00639    Double_t intTrue = ceSum[0]/(nevMC+0.000001);
00640    Double_t intDriv=0.;
00641    Double_t intPrim=0.;
00642 
00643    switch(fOptDrive){
00644    case 1:                       // VARIANCE REDUCTION
00645       if(kBest == -1) Varedu(ceSum,kBest,xBest,yBest); // determine the best edge,
00646       //intDriv =sqrt( ceSum[1]/nevMC -intTrue*intTrue ); // Older ansatz, numerically not bad
00647       intDriv =sqrt(ceSum[1]/nevMC) -intTrue; // Foam build-up, sqrt(<w**2>) -<w>
00648       intPrim =sqrt(ceSum[1]/nevMC);          // MC gen. sqrt(<w**2>) =sqrt(<w>**2 +sigma**2)
00649       break;
00650    case 2:                       // WTMAX  REDUCTION
00651       if(kBest == -1) Carver(kBest,xBest,yBest);  // determine the best edge
00652       intDriv =ceSum[4] -intTrue; // Foam build-up, wtmax-<w>
00653       intPrim =ceSum[4];          // MC generation, wtmax!
00654       break;
00655    default:
00656       Error("Explore", "Wrong fOptDrive = \n" );
00657    }//switch
00658    //=================================================================================
00659    //hist_Neff_distrib.Fill( fLastCe/2.0+0.01, nevEff+0.01);  //
00660    //hist_kBest_distrib.Fill( kBest+0.50, 1.0 ); //  debug
00661    //hist_xBest_distrib.Fill( xBest+0.01, 1.0 ); //  debug
00662    //=================================================================================
00663    cell->SetBest(kBest);
00664    cell->SetXdiv(xBest);
00665    cell->SetIntg(intTrue);
00666    cell->SetDriv(intDriv);
00667    cell->SetPrim(intPrim);
00668    // correct/update integrals in all parent cells to the top of the tree
00669    Double_t  parIntg, parDriv;
00670    for(parent = cell->GetPare(); parent!=0; parent = parent->GetPare()){
00671       parIntg = parent->GetIntg();
00672       parDriv = parent->GetDriv();
00673       parent->SetIntg( parIntg   +intTrue -intOld );
00674       parent->SetDriv( parDriv   +intDriv -driOld );
00675    }
00676    delete [] volPart;
00677    delete [] xRand;
00678    //cell->Print();
00679 } // TFoam::Explore
00680 
00681 //______________________________________________________________________________________
00682 void TFoam::Varedu(Double_t ceSum[5], Int_t &kBest, Double_t &xBest, Double_t &yBest)
00683 {
00684 // Internal subrogram used by Initialize.
00685 // In determines the best edge candidate and the position of the cell division plane
00686 // in case of the variance reduction for future cell division,
00687 // using results of the MC exploration run stored in fHistEdg
00688 
00689    Double_t nent   = ceSum[2];
00690    Double_t swAll  = ceSum[0];
00691    Double_t sswAll = ceSum[1];
00692    Double_t ssw    = sqrt(sswAll)/sqrt(nent);
00693    //
00694    Double_t swIn,swOut,sswIn,sswOut,xLo,xUp;
00695    kBest =-1;
00696    xBest =0.5;
00697    yBest =1.0;
00698    Double_t maxGain=0.0;
00699    // Now go over all projections kProj
00700    for(Int_t kProj=0; kProj<fDim; kProj++) {
00701       if( fMaskDiv[kProj]) {
00702          // initialize search over bins
00703          Double_t sigmIn =0.0; Double_t sigmOut =0.0;
00704          Double_t sswtBest = gHigh;
00705          Double_t gain =0.0;
00706          Double_t xMin=0.0; Double_t xMax=0.0;
00707          // Double loop over all pairs jLo<jUp
00708          for(Int_t jLo=1; jLo<=fNBin; jLo++) {
00709             Double_t aswIn=0;  Double_t asswIn=0;
00710             for(Int_t jUp=jLo; jUp<=fNBin;jUp++) {
00711                aswIn  +=     ((TH1D *)(*fHistEdg)[kProj])->GetBinContent(jUp);
00712                asswIn += Sqr(((TH1D *)(*fHistEdg)[kProj])->GetBinError(  jUp));
00713                xLo=(jLo-1.0)/fNBin;
00714                xUp=(jUp*1.0)/fNBin;
00715                swIn  =        aswIn/nent;
00716                swOut = (swAll-aswIn)/nent;
00717                sswIn = sqrt(asswIn)       /sqrt(nent*(xUp-xLo))     *(xUp-xLo);
00718                sswOut= sqrt(sswAll-asswIn)/sqrt(nent*(1.0-xUp+xLo)) *(1.0-xUp+xLo);
00719                if( (sswIn+sswOut) < sswtBest) {
00720                   sswtBest = sswIn+sswOut;
00721                   gain     = ssw-sswtBest;
00722                   sigmIn   = sswIn -swIn;  // Debug
00723                   sigmOut  = sswOut-swOut; // Debug
00724                   xMin    = xLo;
00725                   xMax    = xUp;
00726                }
00727             }//jUp
00728          }//jLo
00729          Int_t iLo = (Int_t) (fNBin*xMin);
00730          Int_t iUp = (Int_t) (fNBin*xMax);
00731          //----------DEBUG printout
00732          //cout<<"@@@@@  xMin xMax = "<<xMin   <<" "<<xMax<<"  iLo= "<<iLo<<"  iUp= "<<iUp;
00733          //cout<<"  sswtBest/ssw= "<<sswtBest/ssw<<"  Gain/ssw= "<< Gain/ssw<<endl;
00734          //----------DEBUG auxilary Plot
00735          for(Int_t iBin=1;iBin<=fNBin;iBin++) {
00736             if( ((iBin-0.5)/fNBin > xMin) && ((iBin-0.5)/fNBin < xMax) ){
00737                ((TH1D *)(*fHistDbg)[kProj])->SetBinContent(iBin,sigmIn/(xMax-xMin));
00738             } else {
00739                ((TH1D *)(*fHistDbg)[kProj])->SetBinContent(iBin,sigmOut/(1-xMax+xMin));
00740             }
00741          }
00742          if(gain>=maxGain) {
00743             maxGain=gain;
00744             kBest=kProj; // <--- !!!!! The best edge
00745             xBest=xMin;
00746             yBest=xMax;
00747             if(iLo == 0)     xBest=yBest; // The best division point
00748             if(iUp == fNBin) yBest=xBest; // this is not really used
00749          }
00750       }
00751    } //kProj
00752    //----------DEBUG printout
00753    //cout<<"@@@@@@@>>>>> kBest= "<<kBest<<"  maxGain/ssw= "<< maxGain/ssw<<endl;
00754    if( (kBest >= fDim) || (kBest<0) ) Error("Varedu", "Something wrong with kBest \n" );
00755 }          //TFoam::Varedu
00756 
00757 //________________________________________________________________________________________
00758 void TFoam::Carver(Int_t &kBest, Double_t &xBest, Double_t &yBest)
00759 {
00760 // Internal subrogram used by Initialize.
00761 // Determines the best edge-candidate and the position of the division plane
00762 // for the future cell division, in the case of the optimization of the maximum weight.
00763 // It exploits results of the cell MC exploration run stored in fHistEdg.
00764 
00765    Int_t    kProj,iBin;
00766    Double_t carve,carvTot,carvMax,carvOne,binMax,binTot,primTot,primMax;
00767    Int_t    jLow,jUp,iLow,iUp;
00768    Double_t theBin;
00769    Int_t    jDivi; // TEST
00770    Int_t j;
00771 
00772    Double_t *bins  = new Double_t[fNBin];      // bins of histogram for single  PROJECTION
00773    if(bins==0)    Error("Carver", "Cannot initialize buffer Bins \n" );
00774 
00775    kBest =-1;
00776    xBest =0.5;
00777    yBest =1.0;
00778    carvMax = gVlow;
00779    primMax = gVlow;
00780    for(kProj=0; kProj<fDim; kProj++)
00781       if( fMaskDiv[kProj] ){
00782       //if( kProj==1 ){
00783       //cout<<"==================== Carver histogram: kProj ="<<kProj<<"==================="<<endl;
00784       //((TH1D *)(*fHistEdg)[kProj])->Print("all");
00785       binMax = gVlow;
00786       for(iBin=0; iBin<fNBin;iBin++){
00787          bins[iBin]= ((TH1D *)(*fHistEdg)[kProj])->GetBinContent(iBin+1);
00788          binMax = TMath::Max( binMax, bins[iBin]);       // Maximum content/bin
00789       }
00790       if(binMax < 0 ) {       //case of empty cell
00791          delete [] bins;
00792          return;
00793       }
00794       carvTot = 0.0;
00795       binTot  = 0.0;
00796       for(iBin=0;iBin<fNBin;iBin++){
00797          carvTot = carvTot + (binMax-bins[iBin]);     // Total Carve (more stable)
00798          binTot  +=bins[iBin];
00799       }
00800       primTot = binMax*fNBin;
00801       //cout <<"Carver:  CarvTot "<<CarvTot<< "    primTot "<<primTot<<endl;
00802       jLow =0;
00803       jUp  =fNBin-1;
00804       carvOne = gVlow;
00805       Double_t yLevel = gVlow;
00806       for(iBin=0; iBin<fNBin;iBin++) {
00807          theBin = bins[iBin];
00808          //-----  walk to the left and find first bin > theBin
00809          iLow = iBin;
00810          for(j=iBin; j>-1; j-- ) {
00811             if(theBin< bins[j]) break;
00812             iLow = j;
00813          }
00814          //iLow = iBin;
00815          //if(iLow>0)     while( (theBin >= bins[iLow-1])&&(iLow >0) ){iLow--;} // horror!!!
00816          //------ walk to the right and find first bin > theBin
00817          iUp  = iBin;
00818          for(j=iBin; j<fNBin; j++) {
00819             if(theBin< bins[j]) break;
00820             iUp = j;
00821          }
00822          //iUp  = iBin;
00823          //if(iUp<fNBin-1) while( (theBin >= bins[iUp+1])&&( iUp<fNBin-1 ) ){iUp++;} // horror!!!
00824          //
00825          carve = (iUp-iLow+1)*(binMax-theBin);
00826          if( carve > carvOne) {
00827             carvOne = carve;
00828             jLow = iLow;
00829             jUp  = iUp;
00830             yLevel = theBin;
00831          }
00832       }//iBin
00833       if( carvTot > carvMax) {
00834          carvMax   = carvTot;
00835          primMax   = primTot;
00836          //cout <<"Carver:   primMax "<<primMax<<endl;
00837          kBest = kProj;    // Best edge
00838          xBest = ((Double_t)(jLow))/fNBin;
00839          yBest = ((Double_t)(jUp+1))/fNBin;
00840          if(jLow == 0 )       xBest = yBest;
00841          if(jUp  == fNBin-1) yBest = xBest;
00842          // division ratio in units of 1/fNBin, testing
00843          jDivi = jLow;
00844          if(jLow == 0 )     jDivi=jUp+1;
00845       }
00846       //======  extra histograms for debug purposes
00847       //cout<<"kProj= "<<kProj<<" jLow= "<<jLow<<" jUp= "<<jUp<<endl;
00848       for(iBin=0;    iBin<fNBin;  iBin++)
00849          ((TH1D *)(*fHistDbg)[kProj])->SetBinContent(iBin+1,binMax);
00850       for(iBin=jLow; iBin<jUp+1;   iBin++)
00851          ((TH1D *)(*fHistDbg)[kProj])->SetBinContent(iBin+1,yLevel);
00852    }//kProj
00853    if( (kBest >= fDim) || (kBest<0) ) Error("Carver", "Something wrong with kBest \n" );
00854    delete [] bins;
00855 }          //TFoam::Carver
00856 
00857 //______________________________________________________________________________________________
00858 void TFoam::MakeAlpha()
00859 {
00860 // Internal subrogram used by Initialize.
00861 // Provides random vector Alpha  0< Alpha(i) < 1
00862    Int_t k;
00863    if(fDim<1) return;
00864 
00865    // simply generate and load kDim uniform random numbers
00866    fPseRan->RndmArray(fDim,fRvec);   // kDim random numbers needed
00867    for(k=0; k<fDim; k++) fAlpha[k] = fRvec[k];
00868 } //MakeAlpha
00869 
00870 
00871 //_____________________________________________________________________________________________
00872 void TFoam::Grow()
00873 {
00874 // Internal subrogram used by Initialize.
00875 // It grow new cells by the binary division process.
00876 
00877    Long_t iCell;
00878    TFoamCell* newCell;
00879 
00880    while ( (fLastCe+2) < fNCells ) {  // this condition also checked inside Divide
00881       iCell   = PeekMax();            // peek up cell with maximum driver integral
00882       if( (iCell<0) || (iCell>fLastCe) ) Error("Grow", "Wrong iCell \n");
00883       newCell = fCells[iCell];
00884 
00885       if(fLastCe !=0) {
00886          Int_t kEcho=10;
00887          if(fLastCe>=10000) kEcho=100;
00888          if( (fLastCe%kEcho)==0) {
00889            if (fChat>0) {
00890              if(fDim<10)
00891                cout<<fDim<<flush;
00892              else
00893                cout<<"."<<flush;
00894              if( (fLastCe%(100*kEcho))==0)  cout<<"|"<<fLastCe<<endl<<flush;
00895            }
00896          }
00897       }
00898       if( Divide( newCell )==0) break;  // and divide it into two
00899    }
00900    if (fChat>0) {
00901      cout<<endl<<flush;
00902    }
00903    CheckAll(0);   // set arg=1 for more info
00904 }// Grow
00905 
00906 //_____________________________________________________________________________________________
00907 Long_t  TFoam::PeekMax()
00908 {
00909 // Internal subprogram used by Initialize.
00910 // It finds cell with maximal driver integral for the purpose of the division.
00911 
00912    Long_t  i;
00913    Long_t iCell = -1;
00914    Double_t  drivMax, driv;
00915 
00916    drivMax = gVlow;
00917    for(i=0; i<=fLastCe; i++) {//without root
00918       if( fCells[i]->GetStat() == 1 ) {
00919          driv =  TMath::Abs( fCells[i]->GetDriv());
00920          //cout<<"PeekMax: Driv = "<<driv<<endl;
00921          if(driv > drivMax) {
00922             drivMax = driv;
00923             iCell = i;
00924          }
00925       }
00926    }
00927    //  cout << 'TFoam_PeekMax: iCell=' << iCell << endl;
00928    if (iCell == -1)
00929       cout << "STOP in TFoam::PeekMax: not found iCell=" <<  iCell << endl;
00930    return(iCell);
00931 }                 // TFoam_PeekMax
00932 
00933 //_____________________________________________________________________________________________
00934 Int_t TFoam::Divide(TFoamCell *cell)
00935 {
00936 // Internal subrogram used by Initialize.
00937 // It divides cell iCell into two daughter cells.
00938 // The iCell is retained and tagged as inactive, daughter cells are appended
00939 // at the end of the buffer.
00940 // New vertex is added to list of vertices.
00941 // List of active cells is updated, iCell removed, two daughters added
00942 // and their properties set with help of MC sampling (TFoam_Explore)
00943 // Returns Code RC=-1 of buffer limit is reached,  fLastCe=fnBuf.
00944 
00945    Double_t xdiv;
00946    Int_t   kBest;
00947 
00948    if(fLastCe+1 >= fNCells) Error("Divide", "Buffer limit is reached, fLastCe=fnBuf \n");
00949 
00950    cell->SetStat(0); // reset cell as inactive
00951    fNoAct--;
00952 
00953    xdiv  = cell->GetXdiv();
00954    kBest = cell->GetBest();
00955    if( kBest<0 || kBest>=fDim ) Error("Divide", "Wrong kBest \n");
00956 
00957    //////////////////////////////////////////////////////////////////
00958    //           define two daughter cells (active)                 //
00959    //////////////////////////////////////////////////////////////////
00960 
00961    Int_t d1 = CellFill(1,   cell);
00962    Int_t d2 = CellFill(1,   cell);
00963    cell->SetDau0((fCells[d1]));
00964    cell->SetDau1((fCells[d2]));
00965    Explore( (fCells[d1]) );
00966    Explore( (fCells[d2]) );
00967    return 1;
00968 } // TFoam_Divide
00969 
00970 
00971 //_________________________________________________________________________________________
00972 void TFoam::MakeActiveList()
00973 {
00974 // Internal subrogram used by Initialize.
00975 // It finds out number of active cells fNoAct,
00976 // creates list of active cell fCellsAct and primary cumulative fPrimAcu.
00977 // They are used during the MC generation to choose randomly an active cell.
00978 
00979    Long_t n, iCell;
00980    Double_t sum;
00981    // flush previous result
00982    if(fPrimAcu  != 0) delete [] fPrimAcu;
00983    if(fCellsAct != 0) delete fCellsAct;
00984 
00985    // Allocate tables of active cells
00986    fCellsAct = new TRefArray();
00987 
00988    // Count Active cells and find total Primary
00989    // Fill-in tables of active cells
00990 
00991    fPrime = 0.0; n = 0;
00992    for(iCell=0; iCell<=fLastCe; iCell++) { 
00993       if (fCells[iCell]->GetStat()==1) {
00994          fPrime += fCells[iCell]->GetPrim();
00995          fCellsAct->Add(fCells[iCell]);
00996          n++;
00997       }
00998    }
00999 
01000    if(fNoAct != n)  Error("MakeActiveList", "Wrong fNoAct               \n"  );
01001    if(fPrime == 0.) Error("MakeActiveList", "Integrand function is zero  \n"  );
01002 
01003    fPrimAcu  = new  Double_t[fNoAct]; // cumulative primary for MC generation
01004    if( fCellsAct==0 || fPrimAcu==0 ) Error("MakeActiveList", "Cant allocate fCellsAct or fPrimAcu \n");
01005 
01006    sum =0.0;
01007    for(iCell=0; iCell<fNoAct; iCell++) {
01008       sum = sum + ( (TFoamCell *) (fCellsAct->At(iCell)) )->GetPrim()/fPrime;
01009       fPrimAcu[iCell]=sum;
01010    }
01011 
01012 } //MakeActiveList
01013 
01014 //__________________________________________________________________________________________
01015 void TFoam::ResetPseRan(TRandom *PseRan)
01016 {
01017 // User may optionally reset random number generator using this method
01018 // Usually it is done when FOAM object is restored from the disk.
01019 // IMPORTANT: this method deletes existing  random number generator registered in the FOAM object.
01020 // In particular such an object is created by the streamer during the disk-read operation.
01021 
01022    if(fPseRan) {
01023       Info("ResetPseRan", "Resetting random number generator  \n");
01024       delete fPseRan;
01025    }
01026    SetPseRan(PseRan);
01027 }
01028 
01029 //__________________________________________________________________________________________
01030 void TFoam::SetRho(TFoamIntegrand *fun)
01031 {
01032 // User may use this method to set (register) random number generator used by
01033 // the given instance of the FOAM event generator. Note that single r.n. generator
01034 // may serve several FOAM objects.
01035 
01036    if (fun)
01037       fRho=fun;
01038    else
01039       Error("SetRho", "Bad function \n" );
01040 }
01041 
01042 //__________________________________________________________________________________________
01043 void TFoam::ResetRho(TFoamIntegrand *fun)
01044 {
01045 // User may optionally reset the distribution using this method
01046 // Usually it is done when FOAM object is restored from the disk.
01047 // IMPORTANT: this method deletes existing  distribution object registered in the FOAM object.
01048 // In particular such an object is created by the streamer diring the disk-read operation.
01049 // This method is used only in very special cases, because the distribution in most cases
01050 // should be "owned" by the FOAM object and should not be replaced by another one after initialization.
01051 
01052    if(fRho) {
01053       Info("ResetRho", "!!! Resetting distribution function  !!!\n");
01054       delete fRho;
01055    }
01056    SetRho(fun);
01057 }
01058 
01059 //__________________________________________________________________________________________
01060 void TFoam::SetRhoInt(void *fun)
01061 {
01062 // User may use this to set pointer to the global function (not descending
01063 // from TFoamIntegrand) serving as a distribution for FOAM.
01064 // It is useful for simple interactive applications.
01065 // Note that persistency for FOAM object will not work in the case of such
01066 // a distribution.
01067 
01068    const char *namefcn = gCint->Getp2f2funcname(fun); //name of integrand function
01069    if(namefcn) {
01070       fMethodCall=new TMethodCall();
01071       fMethodCall->InitWithPrototype(namefcn, "Int_t, Double_t *");
01072    }
01073    fRho=0;
01074 }
01075 
01076 //__________________________________________________________________________________________
01077 Double_t TFoam::Eval(Double_t *xRand)
01078 {
01079 // Internal subprogram.
01080 // Evaluates distribution to be generated.
01081 
01082    Double_t result;
01083 
01084    if(!fRho) {   //interactive mode
01085       Long_t paramArr[3];
01086       paramArr[0]=(Long_t)fDim;
01087       paramArr[1]=(Long_t)xRand;
01088       fMethodCall->SetParamPtrs(paramArr);
01089       fMethodCall->Execute(result);
01090    } else {       //compiled mode
01091       result=fRho->Density(fDim,xRand);
01092    }
01093 
01094    return result;
01095 }
01096 
01097 //___________________________________________________________________________________________
01098 void TFoam::GenerCel2(TFoamCell *&pCell)
01099 {
01100 // Internal subprogram.
01101 // Return randomly chosen active cell with probability equal to its
01102 // contribution into total driver integral using interpolation search.
01103 
01104    Long_t  lo, hi, hit;
01105    Double_t fhit, flo, fhi;
01106    Double_t random;
01107 
01108    random=fPseRan->Rndm();
01109    lo  = 0;              hi =fNoAct-1;
01110    flo = fPrimAcu[lo];  fhi=fPrimAcu[hi];
01111    while(lo+1<hi) {
01112       hit = lo + (Int_t)( (hi-lo)*(random-flo)/(fhi-flo)+0.5);
01113       if (hit<=lo)
01114          hit = lo+1;
01115       else if(hit>=hi)
01116          hit = hi-1;
01117       fhit=fPrimAcu[hit];
01118       if (fhit>random) {
01119          hi = hit;
01120          fhi = fhit;
01121       } else {
01122          lo = hit;
01123          flo = fhit;
01124       }
01125    }
01126    if (fPrimAcu[lo]>random)
01127       pCell = (TFoamCell *) fCellsAct->At(lo);
01128    else
01129       pCell = (TFoamCell *) fCellsAct->At(hi);
01130 }       // TFoam::GenerCel2
01131 
01132 
01133 //___________________________________________________________________________________________
01134 void TFoam::MakeEvent(void)
01135 {
01136 // User subprogram.
01137 // It generates randomly point/vector according to user-defined distribution.
01138 // Prior initialization with help of Initialize() is mandatory.
01139 // Generated MC point/vector is available using GetMCvect and the MC weight with GetMCwt.
01140 // MC point is generated with wt=1 or with variable weight, see OptRej switch.
01141 
01142    Int_t      j;
01143    Double_t   wt,dx,mcwt;
01144    TFoamCell *rCell;
01145    //
01146    //********************** MC LOOP STARS HERE **********************
01147 ee0:
01148    GenerCel2(rCell);   // choose randomly one cell
01149 
01150    MakeAlpha();
01151 
01152    TFoamVect  cellPosi(fDim); TFoamVect  cellSize(fDim);
01153    rCell->GetHcub(cellPosi,cellSize);
01154    for(j=0; j<fDim; j++)
01155       fMCvect[j]= cellPosi[j] +fAlpha[j]*cellSize[j];
01156    dx = rCell->GetVolume();      // Cartesian volume of the Cell
01157    //  weight average normalized to PRIMARY integral over the cell
01158 
01159    wt=dx*Eval(fMCvect);
01160 
01161    mcwt = wt / rCell->GetPrim();  // PRIMARY controls normalization
01162    fNCalls++;
01163    fMCwt   =  mcwt;
01164    // accumulation of statistics for the main MC weight
01165    fSumWt  += mcwt;           // sum of Wt
01166    fSumWt2 += mcwt*mcwt;      // sum of Wt**2
01167    fNevGen++;                 // sum of 1d0
01168    fWtMax  =  TMath::Max(fWtMax, mcwt);   // maximum wt
01169    fWtMin  =  TMath::Min(fWtMin, mcwt);   // minimum wt
01170    fMCMonit->Fill(mcwt);
01171    fHistWt->Fill(mcwt,1.0);          // histogram
01172    //*******  Optional rejection ******
01173    if(fOptRej == 1) {
01174       Double_t random;
01175       random=fPseRan->Rndm();
01176       if( fMaxWtRej*random > fMCwt) goto ee0;  // Wt=1 events, internal rejection
01177       if( fMCwt<fMaxWtRej ) {
01178          fMCwt = 1.0;                  // normal Wt=1 event
01179       } else {
01180          fMCwt = fMCwt/fMaxWtRej;    // weight for overweighted events! kept for debug
01181          fSumOve += fMCwt-fMaxWtRej; // contribution of overweighted
01182       }
01183    }
01184    //********************** MC LOOP ENDS HERE **********************
01185 } // MakeEvent
01186 
01187 //_________________________________________________________________________________
01188 void TFoam::GetMCvect(Double_t *MCvect)
01189 {
01190 // User may get generated MC point/vector with help of this method
01191 
01192    for ( Int_t k=0 ; k<fDim ; k++) *(MCvect +k) = fMCvect[k];
01193 }//GetMCvect
01194 
01195 //___________________________________________________________________________________
01196 Double_t TFoam::GetMCwt(void)
01197 {
01198 // User may get weight MC weight using this method
01199 
01200    return(fMCwt);
01201 }
01202 //___________________________________________________________________________________
01203 void TFoam::GetMCwt(Double_t &mcwt)
01204 {
01205 // User may get weight MC weight using this method
01206 
01207    mcwt=fMCwt;
01208 }
01209 
01210 //___________________________________________________________________________________
01211 Double_t TFoam::MCgenerate(Double_t *MCvect)
01212 {
01213 // User subprogram which generates MC event and returns MC weight
01214 
01215    MakeEvent();
01216    GetMCvect(MCvect);
01217    return(fMCwt);
01218 }//MCgenerate
01219 
01220 //___________________________________________________________________________________
01221 void TFoam::GetIntegMC(Double_t &mcResult, Double_t &mcError)
01222 {
01223 // User subprogram.
01224 // It provides the value of the integral calculated from the averages of the MC run
01225 // May be called after (or during) the MC run.
01226 
01227    Double_t mCerelat;
01228    mcResult = 0.0;
01229    mCerelat = 1.0;
01230    if (fNevGen>0) {
01231       mcResult = fPrime*fSumWt/fNevGen;
01232       mCerelat = sqrt( fSumWt2/(fSumWt*fSumWt) - 1/fNevGen);
01233    }
01234    mcError = mcResult *mCerelat;
01235 }//GetIntegMC
01236 
01237 //____________________________________________________________________________________
01238 void  TFoam::GetIntNorm(Double_t& IntNorm, Double_t& Errel )
01239 {
01240 // User subprogram.
01241 // It returns NORMALIZATION integral to be combined with the average weights
01242 // and content of the histograms in order to get proper absolute normalization
01243 // of the integrand and distributions.
01244 // It can be called after initialization, before or during the MC run.
01245 
01246    if(fOptRej == 1) {    // Wt=1 events, internal rejection
01247       Double_t intMC,errMC;
01248       GetIntegMC(intMC,errMC);
01249       IntNorm = intMC;
01250       Errel   = errMC;
01251    } else {                // Wted events, NO internal rejection
01252       IntNorm = fPrime;
01253       Errel   = 0;
01254    }
01255 }//GetIntNorm
01256 
01257 //______________________________________________________________________________________
01258 void  TFoam::GetWtParams(Double_t eps, Double_t &aveWt, Double_t &wtMax, Double_t &sigma)
01259 {
01260 // May be called optionally after the MC run.
01261 // Returns various parameters of the MC weight for efficiency evaluation
01262 
01263    Double_t mCeff, wtLim;
01264    fMCMonit->GetMCeff(eps, mCeff, wtLim);
01265    wtMax = wtLim;
01266    aveWt = fSumWt/fNevGen;
01267    sigma = sqrt( fSumWt2/fNevGen -aveWt*aveWt );
01268 }//GetmCeff
01269 
01270 //_______________________________________________________________________________________
01271 void TFoam::Finalize(Double_t& IntNorm, Double_t& Errel)
01272 {
01273 // May be called optionally by the user after the MC run.
01274 // It provides normalization and also prints some information/statistics on the MC run.
01275 
01276    GetIntNorm(IntNorm,Errel);
01277    Double_t mcResult,mcError;
01278    GetIntegMC(mcResult,mcError);
01279    Double_t mCerelat= mcError/mcResult;
01280    //
01281    if(fChat>0) {
01282       Double_t eps = 0.0005;
01283       Double_t mCeff, mcEf2, wtMax, aveWt, sigma;
01284       GetWtParams(eps, aveWt, wtMax, sigma);
01285       mCeff=0;
01286       if(wtMax>0.0) mCeff=aveWt/wtMax;
01287       mcEf2 = sigma/aveWt;
01288       Double_t driver = fCells[0]->GetDriv();
01289       //
01290       BXOPE;
01291       BXTXT("****************************************");
01292       BXTXT("******     TFoam::Finalize       ******");
01293       BXTXT("****************************************");
01294       BX1I("    NevGen",fNevGen, "Number of generated events in the MC generation   ");
01295       BX1I("    nCalls",fNCalls, "Total number of function calls                    ");
01296       BXTXT("----------------------------------------");
01297       BX1F("     AveWt",aveWt,    "Average MC weight                      ");
01298       BX1F("     WtMin",fWtMin,  "Minimum MC weight (absolute)           ");
01299       BX1F("     WtMax",fWtMax,  "Maximum MC weight (absolute)           ");
01300       BXTXT("----------------------------------------");
01301       BX1F("    XPrime",fPrime,  "Primary total integral, R_prime        ");
01302       BX1F("    XDiver",driver,   "Driver  total integral, R_loss         ");
01303       BXTXT("----------------------------------------");
01304       BX2F("    IntMC", mcResult,  mcError,      "Result of the MC Integral");
01305       BX1F(" mCerelat", mCerelat,  "Relative error of the MC integral      ");
01306       BX1F(" <w>/WtMax",mCeff,     "MC efficiency, acceptance rate");
01307       BX1F(" Sigma/<w>",mcEf2,     "MC efficiency, variance/ave_wt");
01308       BX1F("     WtMax",wtMax,     "WtMax(esp= 0.0005)            ");
01309       BX1F("     Sigma",sigma,     "variance of MC weight         ");
01310       if(fOptRej==1) {
01311          Double_t avOve=fSumOve/fSumWt;
01312          BX1F("<OveW>/<W>",avOve,     "Contrib. of events wt>MaxWtRej");
01313       }
01314       BXCLO;
01315    }
01316 }  // Finalize
01317 
01318 //_____________________________________________________________________________________
01319 void  TFoam::SetInhiDiv(Int_t iDim, Int_t InhiDiv)
01320 {
01321 // This can be called before Initialize, after setting kDim
01322 // It defines which variables are excluded in the process of the cell division.
01323 // For example 'FoamX->SetInhiDiv(1, 1);' inhibits division of y-variable.
01324 // The resulting map of cells in 2-dim. case will look as follows:
01325 //BEGIN_HTML <!--
01326 /* -->
01327 <img src="gif/foam_Map2.gif">
01328 <!--*/
01329 // -->END_HTML
01330 
01331    if(fDim==0) Error("TFoam","SetInhiDiv: fDim=0 \n");
01332    if(fInhiDiv == 0) {
01333       fInhiDiv = new Int_t[ fDim ];
01334       for(Int_t i=0; i<fDim; i++) fInhiDiv[i]=0;
01335    }
01336    //
01337    if( ( 0<=iDim) && (iDim<fDim)) {
01338       fInhiDiv[iDim] = InhiDiv;
01339    } else
01340       Error("SetInhiDiv:","Wrong iDim \n");
01341 }//SetInhiDiv
01342 
01343 //______________________________________________________________________________________
01344 void  TFoam::SetXdivPRD(Int_t iDim, Int_t len, Double_t xDiv[])
01345 {
01346 // This should be called before Initialize, after setting  kDim
01347 // It predefines values of the cell division for certain variable iDim.
01348 // For example setting 3 predefined division lines using:
01349 //     xDiv[0]=0.30; xDiv[1]=0.40; xDiv[2]=0.65;
01350 //     FoamX->SetXdivPRD(0,3,xDiv);
01351 // results in the following 2-dim. pattern of the cells:
01352 //BEGIN_HTML <!--
01353 /* -->
01354 <img src="gif/foam_Map3.gif">
01355 <!--*/
01356 // -->END_HTML
01357 
01358    Int_t i;
01359 
01360    if(fDim<=0)  Error("SetXdivPRD", "fDim=0 \n");
01361    if(   len<1 )  Error("SetXdivPRD", "len<1 \n");
01362    // allocate list of pointers, if it was not done before
01363    if(fXdivPRD == 0) {
01364       fXdivPRD = new TFoamVect*[fDim];
01365       for(i=0; i<fDim; i++)  fXdivPRD[i]=0;
01366    }
01367   // set division list for direction iDim in H-cubic space!!!
01368    if( ( 0<=iDim) && (iDim<fDim)) {
01369       fOptPRD =1;      // !!!!
01370       if( fXdivPRD[iDim] != 0)
01371          Error("SetXdivPRD", "Second allocation of XdivPRD not allowed \n");
01372       fXdivPRD[iDim] = new TFoamVect(len); // allocate list of division points
01373       for(i=0; i<len; i++) {
01374          (*fXdivPRD[iDim])[i]=xDiv[i]; // set list of division points
01375       }
01376    } else {
01377       Error("SetXdivPRD", "Wrong iDim  \n");
01378    }
01379    // Priting predefined division points
01380    cout<<" SetXdivPRD, idim= "<<iDim<<"  len= "<<len<<"   "<<endl;
01381    for(i=0; i<len; i++) {
01382       if (iDim < fDim) cout<< (*fXdivPRD[iDim])[i] <<"  ";
01383    }
01384    cout<<endl;
01385    for(i=0; i<len; i++)  cout<< xDiv[i] <<"   ";
01386    cout<<endl;
01387    //
01388 }//SetXdivPRD
01389 
01390 //_______________________________________________________________________________________
01391 void TFoam::CheckAll(Int_t level)
01392 {
01393 //  User utility, miscellaneous and debug.
01394 //  Checks all pointers in the tree of cells. This is useful autodiagnostic.
01395 //  level=0, no printout, failures causes STOP
01396 //  level=1, printout, failures lead to WARNINGS only
01397 
01398    Int_t errors, warnings;
01399    TFoamCell *cell;
01400    Long_t iCell;
01401 
01402    errors = 0; warnings = 0;
01403    if (level==1) cout << "///////////////////////////// FOAM_Checks /////////////////////////////////" << endl;
01404    for(iCell=1; iCell<=fLastCe; iCell++) {
01405       cell = fCells[iCell];
01406       //  checking general rules
01407       if( ((cell->GetDau0()==0) && (cell->GetDau1()!=0) ) ||
01408          ((cell->GetDau1()==0) && (cell->GetDau0()!=0) ) ) {
01409          errors++;
01410          if (level==1) Error("CheckAll","ERROR: Cell's no %ld has only one daughter \n",iCell);
01411       }
01412       if( (cell->GetDau0()==0) && (cell->GetDau1()==0) && (cell->GetStat()==0) ) {
01413          errors++;
01414          if (level==1) Error("CheckAll","ERROR: Cell's no %ld  has no daughter and is inactive \n",iCell);
01415       }
01416       if( (cell->GetDau0()!=0) && (cell->GetDau1()!=0) && (cell->GetStat()==1) ) {
01417          errors++;
01418          if (level==1) Error("CheckAll","ERROR: Cell's no %ld has two daughters and is active \n",iCell);
01419       }
01420 
01421       // checking parents
01422       if( (cell->GetPare())!=fCells[0] ) { // not child of the root
01423          if ( (cell != cell->GetPare()->GetDau0()) && (cell != cell->GetPare()->GetDau1()) ) {
01424             errors++;
01425             if (level==1) Error("CheckAll","ERROR: Cell's no %ld parent not pointing to this cell\n ",iCell);
01426          }
01427       }
01428 
01429       // checking daughters
01430       if(cell->GetDau0()!=0) {
01431          if(cell != (cell->GetDau0())->GetPare()) {
01432             errors++;
01433             if (level==1)  Error("CheckAll","ERROR: Cell's no %ld daughter 0 not pointing to this cell \n",iCell);
01434          }
01435       }
01436       if(cell->GetDau1()!=0) {
01437          if(cell != (cell->GetDau1())->GetPare()) {
01438             errors++;
01439             if (level==1) Error("CheckAll","ERROR: Cell's no %ld daughter 1 not pointing to this cell \n",iCell);
01440          }
01441       }
01442    }// loop after cells;
01443 
01444    // Check for empty cells
01445    for(iCell=0; iCell<=fLastCe; iCell++) {
01446       cell = fCells[iCell];
01447       if( (cell->GetStat()==1) && (cell->GetDriv()==0) ) {
01448          warnings++;
01449          if(level==1) Warning("CheckAll", "Warning: Cell no. %ld is active but empty \n", iCell);
01450       }
01451    }
01452    // summary
01453    if(level==1){
01454       Info("CheckAll","Check has found %d errors and %d warnings \n",errors, warnings);
01455    }
01456    if(errors>0){
01457       Info("CheckAll","Check - found total %d  errors \n",errors);
01458    }
01459 } // Check
01460 
01461 //________________________________________________________________________________________
01462 void TFoam::PrintCells(void)
01463 {
01464 // Prints geometry of ALL cells of the FOAM
01465 
01466    Long_t iCell;
01467 
01468    for(iCell=0; iCell<=fLastCe; iCell++) {
01469       cout<<"Cell["<<iCell<<"]={ ";
01470       //cout<<"  "<< fCells[iCell]<<"  ";  // extra DEBUG
01471       cout<<endl;
01472       fCells[iCell]->Print("1");
01473       cout<<"}"<<endl;
01474    }
01475 }
01476 
01477 //_________________________________________________________________________________________
01478 void TFoam::RootPlot2dim(Char_t *filename)
01479 {
01480 // Debugging tool which plots 2-dimensional cells as rectangles
01481 // in C++ format readable for root
01482 
01483    ofstream outfile(filename, ios::out);
01484    Double_t   x1,y1,x2,y2,x,y;
01485    Long_t    iCell;
01486    Double_t offs =0.1;
01487    Double_t lpag   =1-2*offs;
01488    outfile<<"{" << endl;
01489    outfile<<"cMap = new TCanvas(\"Map1\",\" Cell Map \",600,600);"<<endl;
01490    //
01491    outfile<<"TBox*a=new TBox();"<<endl;
01492    outfile<<"a->SetFillStyle(0);"<<endl;  // big frame
01493    outfile<<"a->SetLineWidth(4);"<<endl;
01494    outfile<<"a->SetLineColor(2);"<<endl;
01495    outfile<<"a->DrawBox("<<offs<<","<<offs<<","<<(offs+lpag)<<","<<(offs+lpag)<<");"<<endl;
01496    //
01497    outfile<<"TText*t=new TText();"<<endl;  // text for numbering
01498    outfile<<"t->SetTextColor(4);"<<endl;
01499    if(fLastCe<51)
01500       outfile<<"t->SetTextSize(0.025);"<<endl;  // text for numbering
01501    else if(fLastCe<251)
01502       outfile<<"t->SetTextSize(0.015);"<<endl;
01503    else
01504       outfile<<"t->SetTextSize(0.008);"<<endl;
01505    //
01506    outfile<<"TBox*b=new TBox();"<<endl;  // single cell
01507    outfile <<"b->SetFillStyle(0);"<<endl;
01508    //
01509    if(fDim==2 && fLastCe<=2000) {
01510       TFoamVect  cellPosi(fDim); TFoamVect  cellSize(fDim);
01511       outfile << "// =========== Rectangular cells  ==========="<< endl;
01512       for(iCell=1; iCell<=fLastCe; iCell++) {
01513          if( fCells[iCell]->GetStat() == 1) {
01514             fCells[iCell]->GetHcub(cellPosi,cellSize);
01515             x1 = offs+lpag*(        cellPosi[0]); y1 = offs+lpag*(        cellPosi[1]);
01516             x2 = offs+lpag*(cellPosi[0]+cellSize[0]); y2 = offs+lpag*(cellPosi[1]+cellSize[1]);
01517             //     cell rectangle
01518             if(fLastCe<=2000)
01519             outfile<<"b->DrawBox("<<x1<<","<<y1<<","<<x2<<","<<y2<<");"<<endl;
01520             //     cell number
01521             if(fLastCe<=250) {
01522                x = offs+lpag*(cellPosi[0]+0.5*cellSize[0]); y = offs+lpag*(cellPosi[1]+0.5*cellSize[1]);
01523                outfile<<"t->DrawText("<<x<<","<<y<<","<<"\""<<iCell<<"\""<<");"<<endl;
01524             }
01525          }
01526       }
01527       outfile<<"// ============== End Rectangles ==========="<< endl;
01528    }//kDim=2
01529    //
01530    //
01531    outfile << "}" << endl;
01532    outfile.close();
01533 }
01534 
01535 void TFoam::LinkCells()
01536 {
01537 // Void function for backward compatibility
01538 
01539    Info("LinkCells", "VOID function for backward compatibility \n");
01540    return;
01541 }
01542 
01543 ////////////////////////////////////////////////////////////////////////////////
01544 //       End of Class TFoam                                                   //
01545 ////////////////////////////////////////////////////////////////////////////////
01546 

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