From $ROOTSYS/tutorials/foam/foam_kanwa.C

// This program can be execute from the command line as folows:
//
//      root -l foam_kanwa.C
//_____________________________________________________________________________

#include "Riostream.h"
#include "TFoam.h"
#include "TCanvas.h"
#include "TH2.h"
#include "TMath.h"
#include "TFoamIntegrand.h"
#include "TRandom3.h"

//_____________________________________________________________________________
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= 0;
  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
//_____________________________________________________________________________

Int_t foam_kanwa(){
  cout<<"--- kanwa started ---"<<endl;
  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
  //
  TRandom     *PseRan   = new TRandom3();  // Create random number generator
  PseRan->SetSeed(4357);
  TFoam   *FoamX    = new TFoam("FoamX");   // Create Simulator
  FoamX->SetkDim(2);         // No. of dimensions, obligatory!
  FoamX->SetnCells(500);     // Optionally No. of cells, default=2000
  FoamX->SetRhoInt(Camel2);  // Set 2-dim distribution, included below
  FoamX->SetPseRan(PseRan);  // Set random number generator
  FoamX->Initialize();       // Initialize simulator, may take time...
  //
  // visualising generated distribution
  TCanvas *cKanwa = new TCanvas("cKanwa","Canvas for plotting",600,600);
  cKanwa->cd();
  // From now on FoamX is ready to generate events
  int nshow=5000;
  for(long 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);
    // live plot
    if(loop == nshow){
      nshow += 5000;
      hst_xy->Draw("lego2");
      cKanwa->Update();
    }
  }// loop
  //
  hst_xy->Draw("lego2");  // final plot
  cKanwa->Update();
  //
  Double_t MCresult, MCerror;
  FoamX->GetIntegMC( MCresult, MCerror);  // get MC integral, should be one
  cout << " MCresult= " << MCresult << " +- " << MCerror <<endl;
  cout<<"--- kanwa ended ---"<<endl;

  return 0;
}//kanwa

 foam_kanwa.C:1
 foam_kanwa.C:2
 foam_kanwa.C:3
 foam_kanwa.C:4
 foam_kanwa.C:5
 foam_kanwa.C:6
 foam_kanwa.C:7
 foam_kanwa.C:8
 foam_kanwa.C:9
 foam_kanwa.C:10
 foam_kanwa.C:11
 foam_kanwa.C:12
 foam_kanwa.C:13
 foam_kanwa.C:14
 foam_kanwa.C:15
 foam_kanwa.C:16
 foam_kanwa.C:17
 foam_kanwa.C:18
 foam_kanwa.C:19
 foam_kanwa.C:20
 foam_kanwa.C:21
 foam_kanwa.C:22
 foam_kanwa.C:23
 foam_kanwa.C:24
 foam_kanwa.C:25
 foam_kanwa.C:26
 foam_kanwa.C:27
 foam_kanwa.C:28
 foam_kanwa.C:29
 foam_kanwa.C:30
 foam_kanwa.C:31
 foam_kanwa.C:32
 foam_kanwa.C:33
 foam_kanwa.C:34
 foam_kanwa.C:35
 foam_kanwa.C:36
 foam_kanwa.C:37
 foam_kanwa.C:38
 foam_kanwa.C:39
 foam_kanwa.C:40
 foam_kanwa.C:41
 foam_kanwa.C:42
 foam_kanwa.C:43
 foam_kanwa.C:44
 foam_kanwa.C:45
 foam_kanwa.C:46
 foam_kanwa.C:47
 foam_kanwa.C:48
 foam_kanwa.C:49
 foam_kanwa.C:50
 foam_kanwa.C:51
 foam_kanwa.C:52
 foam_kanwa.C:53
 foam_kanwa.C:54
 foam_kanwa.C:55
 foam_kanwa.C:56
 foam_kanwa.C:57
 foam_kanwa.C:58
 foam_kanwa.C:59
 foam_kanwa.C:60
 foam_kanwa.C:61
 foam_kanwa.C:62
 foam_kanwa.C:63
 foam_kanwa.C:64
 foam_kanwa.C:65
 foam_kanwa.C:66
 foam_kanwa.C:67
 foam_kanwa.C:68
 foam_kanwa.C:69
 foam_kanwa.C:70
 foam_kanwa.C:71
 foam_kanwa.C:72
 foam_kanwa.C:73
 foam_kanwa.C:74