foam_demo.C

Go to the documentation of this file.
00001 // Demonstrate the TFoam class.
00002 //
00003 //  To run this macro type from CINT command line
00004 //
00005 //  root [0] gSystem->Load("libFoam.so")
00006 //  root [1] .x foam_demo.C+
00007 //Author: Stascek Jadach
00008 //____________________________________________________________________________
00009 
00010 #include "Riostream.h"
00011 #include "TFile.h"
00012 #include "TFoam.h"
00013 #include "TH1.h"
00014 #include "TMath.h"
00015 #include "TFoamIntegrand.h"
00016 #include "TRandom3.h"
00017 
00018 class TFDISTR: public TFoamIntegrand {
00019 public:
00020   TFDISTR(){};
00021   Double_t Density(int nDim, Double_t *Xarg){
00022   // Integrand for mFOAM
00023   Double_t Fun1,Fun2,R1,R2;
00024   Double_t pos1=1e0/3e0;
00025   Double_t pos2=2e0/3e0;
00026   Double_t Gam1= 0.100e0;  // as in JPC
00027   Double_t Gam2= 0.100e0;  // as in JPC
00028   Double_t sPi = sqrt(TMath::Pi());
00029   Double_t xn1=1e0;
00030   Double_t xn2=1e0;
00031   int i;
00032   R1=0;
00033   R2=0;
00034   for(i = 0 ; i<nDim ; i++){
00035     R1=R1+(Xarg[i] -pos1)*(Xarg[i] -pos1);
00036     R2=R2+(Xarg[i] -pos2)*(Xarg[i] -pos2);
00037     xn1=xn1*Gam1*sPi;
00038     xn2=xn2*Gam2*sPi;      
00039   }
00040   R1   = sqrt(R1);
00041   R2   = sqrt(R2);
00042   Fun1 = exp(-(R1*R1)/(Gam1*Gam1))/xn1;  // Gaussian delta-like profile
00043   Fun2 = exp(-(R2*R2)/(Gam2*Gam2))/xn2;  // Gaussian delta-like profile
00044   return 0.5e0*(Fun1+ Fun2);
00045 }
00046   ClassDef(TFDISTR,1) //Class of testing functions for FOAM
00047 };
00048 ClassImp(TFDISTR)
00049 
00050 Int_t foam_demo()
00051 {
00052   TFile RootFile("foam_demo.root","RECREATE","histograms");
00053   long   loop;
00054   Double_t MCresult,MCerror,MCwt;
00055 //=========================================================
00056   long NevTot   =     50000;   // Total MC statistics
00057   Int_t  kDim   =         2;   // total dimension
00058   Int_t  nCells   =     500;   // Number of Cells
00059   Int_t  nSampl   =     200;   // Number of MC events per cell in build-up
00060   Int_t  nBin     =       8;   // Number of bins in build-up
00061   Int_t  OptRej   =       1;   // Wted events for OptRej=0; wt=1 for OptRej=1 (default)
00062   Int_t  OptDrive =       2;   // (D=2) Option, type of Drive =0,1,2 for TrueVol,Sigma,WtMax
00063   Int_t  EvPerBin =      25;   // Maximum events (equiv.) per bin in buid-up
00064   Int_t  Chat     =       1;   // Chat level
00065 //=========================================================
00066   TRandom *PseRan   = new TRandom3();  // Create random number generator
00067   TFoam   *FoamX    = new TFoam("FoamX");   // Create Simulator
00068   TFoamIntegrand    *rho= new TFDISTR(); 
00069   PseRan->SetSeed(4357);
00070 //=========================================================
00071   cout<<"*****   Demonstration Program for Foam version "<<FoamX->GetVersion()<<"    *****"<<endl;
00072   FoamX->SetkDim(        kDim);      // Mandatory!!!
00073   FoamX->SetnCells(      nCells);    // optional
00074   FoamX->SetnSampl(      nSampl);    // optional
00075   FoamX->SetnBin(        nBin);      // optional
00076   FoamX->SetOptRej(      OptRej);    // optional
00077   FoamX->SetOptDrive(    OptDrive);  // optional
00078   FoamX->SetEvPerBin(    EvPerBin);  // optional
00079   FoamX->SetChat(        Chat);      // optional
00080 //===============================
00081   FoamX->SetRho(rho);
00082   FoamX->SetPseRan(PseRan);
00083   FoamX->Initialize(); // Initialize simulator
00084   FoamX->Write("FoamX");     // Writing Foam on the disk, TESTING PERSISTENCY!!!
00085 //===============================
00086   long nCalls=FoamX->GetnCalls();
00087   cout << "====== Initialization done, entering MC loop" << endl;
00088 //======================================================================
00089   //cout<<" About to start MC loop: ";  cin.getline(question,20);
00090   Double_t *MCvect =new Double_t[kDim]; // vector generated in the MC run
00091 //======================================================================
00092   TH1D  *hst_Wt = new TH1D("hst_Wt" ,  "Main weight of Foam",25,0,1.25);
00093   hst_Wt->Sumw2();
00094 //======================================================================
00095   for(loop=0; loop<NevTot; loop++){
00096 //===============================
00097     FoamX->MakeEvent();           // generate MC event
00098 //===============================
00099     FoamX->GetMCvect( MCvect);
00100     MCwt=FoamX->GetMCwt();
00101     hst_Wt->Fill(MCwt,1.0);    
00102     if(loop<15){
00103       cout<<"MCwt= "<<MCwt<<",  ";
00104       cout<<"MCvect= ";
00105       for ( Int_t k=0 ; k<kDim ; k++) cout<<MCvect[k]<<" "; cout<< endl;
00106     }
00107     if( ((loop)%100000)==0 ){
00108       cout<<"   loop= "<<loop<<endl;
00109     }
00110   }
00111 //======================================================================
00112 //======================================================================
00113   cout << "====== Events generated, entering Finalize" << endl;
00114 
00115   hst_Wt->Print("all");
00116   Double_t eps = 0.0005;
00117   Double_t Effic, WtMax, AveWt, Sigma;
00118   Double_t IntNorm, Errel;
00119   FoamX->Finalize(   IntNorm, Errel);     // final printout
00120   FoamX->GetIntegMC( MCresult, MCerror);  // get MC intnegral
00121   FoamX->GetWtParams(eps, AveWt, WtMax, Sigma); // get MC wt parameters
00122   Effic=0; if(WtMax>0) Effic=AveWt/WtMax;
00123   cout << "================================================================" << endl;
00124   cout << " MCresult= " << MCresult << " +- " << MCerror << " RelErr= "<< MCerror/MCresult << endl;
00125   cout << " Dispersion/<wt>= " << Sigma/AveWt << endl;
00126   cout << "      <wt>/WtMax= " << Effic <<",    for epsilon = "<<eps << endl;
00127   cout << " nCalls (initialization only) =   " << nCalls << endl;
00128   cout << "================================================================" << endl;
00129 
00130   delete [] MCvect;
00131   //
00132   RootFile.ls();
00133   RootFile.Write();
00134   RootFile.Close();
00135   cout << "***** End of Demonstration Program  *****" << endl;
00136 
00137   return 0;
00138 } // end of Demo
00139 
00140 

Generated on Tue Jul 5 15:44:11 2011 for ROOT_528-00b_version by  doxygen 1.5.1