ROOT  6.06/08
Reference Guide
TFoam.cxx
Go to the documentation of this file.
1 // @(#)root/foam:$Id$
2 // Author: S. Jadach <mailto:Stanislaw.jadach@ifj.edu.pl>, P.Sawicki <mailto:Pawel.Sawicki@ifj.edu.pl>
3 
4 //______________________________________________________________________________
5 //
6 // FOAM Version 1.02M
7 // ===================
8 // Authors:
9 // S. Jadach and P.Sawicki
10 // Institute of Nuclear Physics, Cracow, Poland
11 // Stanislaw. Jadach@ifj.edu.pl, Pawel.Sawicki@ifj.edu.pl
12 //
13 // What is FOAM for?
14 // =================
15 // * Suppose you want to generate randomly points (vectors) according to
16 // an arbitrary probability distribution in n dimensions,
17 // for which you supply your own subprogram. FOAM can do it for you!
18 // Even if your distributions has quite strong peaks and is discontinuous!
19 // * FOAM generates random points with weight one or with variable weight.
20 // * FOAM is capable to integrate using efficient "adaptive" MC method.
21 // (The distribution does not need to be normalized to one.)
22 // How does it work?
23 // =================
24 // FOAM is the simplified version of the multi-dimensional general purpose
25 // Monte Carlo event generator (integrator) FOAM.
26 // It creates hyper-rectangular "foam of cells", which is more dense around its peaks.
27 // See the following 2-dim. example of the map of 1000 cells for doubly peaked distribution:
28 //
29 // <img src="gif/foam_MapCamel1000.gif">
30 //
31 // FOAM is now fully integrated with the ROOT package.
32 // The important bonus of the ROOT use is persistency of the FOAM objects!
33 //
34 // For more sophisticated problems full version of FOAM may be more appropriate:
35 // See [full version of FOAM](http://jadach.home.cern.ch/jadach/Foam/Index.html)
36 //
37 // Simple example of the use of FOAM:
38 // ==================================
39 // Int_t kanwa(){
40 // gSystem->Load("libFoam");
41 // TH2D *hst_xy = new TH2D("hst_xy" , "x-y plot", 50,0,1.0, 50,0,1.0);
42 // Double_t *MCvect =new Double_t[2]; // 2-dim vector generated in the MC run
43 // TRandom3 *PseRan = new TRandom3(); // Create random number generator
44 // PseRan->SetSeed(4357); // Set seed
45 // TFoam *FoamX = new TFoam("FoamX"); // Create Simulator
46 // FoamX->SetkDim(2); // No. of dimensions, obligatory!
47 // FoamX->SetnCells(500); // No. of cells, can be omitted, default=2000
48 // FoamX->SetRhoInt(Camel2); // Set 2-dim distribution, included below
49 // FoamX->SetPseRan(PseRan); // Set random number generator
50 // FoamX->Initialize(); // Initialize simulator, takes a few seconds...
51 // // From now on FoamX is ready to generate events according to Camel2(x,y)
52 // for(Long_t loop=0; loop<100000; loop++){
53 // FoamX->MakeEvent(); // generate MC event
54 // FoamX->GetMCvect( MCvect); // get generated vector (x,y)
55 // Double_t x=MCvect[0];
56 // Double_t y=MCvect[1];
57 // if(loop<10) std::cout<<"(x,y) = ( "<< x <<", "<< y <<" )"<<std::endl;
58 // hst_xy->Fill(x,y); // fill scattergram
59 // }// loop
60 // Double_t mcResult, mcError;
61 // FoamX->GetIntegMC( mcResult, mcError); // get MC integral, should be one
62 // std::cout << " mcResult= " << mcResult << " +- " << mcError <<std::endl;
63 // // now hst_xy will be plotted visualizing generated distribution
64 // TCanvas *cKanwa = new TCanvas("cKanwa","Canvas for plotting",600,600);
65 // cKanwa->cd();
66 // hst_xy->Draw("lego2");
67 // }//kanwa
68 // Double_t sqr(Double_t x){return x*x;};
69 // Double_t Camel2(Int_t nDim, Double_t *Xarg){
70 // // 2-dimensional distribution for FOAM, normalized to one (within 1e-5)
71 // Double_t x=Xarg[0];
72 // Double_t y=Xarg[1];
73 // Double_t GamSq= sqr(0.100e0);
74 // Double_t Dist=exp(-(sqr(x-1./3) +sqr(y-1./3))/GamSq)/GamSq/TMath::Pi();
75 // Dist +=exp(-(sqr(x-2./3) +sqr(y-2./3))/GamSq)/GamSq/TMath::Pi();
76 // return 0.5*Dist;
77 // }// Camel2
78 // Two-dim. histogram of the MC points generated with the above program looks as follows:
79 //
80 // <img src="gif/foam_cKanwa.gif">
81 //
82 // Canonical nine steering parameters of FOAM
83 // ===========================================
84 //------------------------------------------------------------------------------
85 // Name | default | Description
86 //------------------------------------------------------------------------------
87 // kDim | 0 | Dimension of the integration space. Must be redefined!
88 // nCells | 1000 | No of allocated number of cells,
89 // nSampl | 200 | No. of MC events in the cell MC exploration
90 // nBin | 8 | No. of bins in edge-histogram in cell exploration
91 // OptRej | 1 | OptRej = 0, weighted; OptRej=1, wt=1 MC events
92 // OptDrive | 2 | Maximum weight reduction, =1 for variance reduction
93 // EvPerBin | 25 | Maximum number of the effective wt=1 events/bin,
94 // | | EvPerBin=0 deactivates this option
95 // Chat | 1 | =0,1,2 is the ``chat level'' in the standard output
96 // MaxWtRej | 1.1 | Maximum weight used to get w=1 MC events
97 //------------------------------------------------------------------------------
98 // The above can be redefined before calling 'Initialize()' method,
99 // for instance FoamObject->SetkDim(15) sets dimension of the distribution to 15.
100 // Only kDim HAS TO BE redefined, the other parameters may be left at their defaults.
101 // nCell may be increased up to about million cells for wildly peaked distributions.
102 // Increasing nSampl sometimes helps, but it may cost CPU time.
103 // MaxWtRej may need to be increased for wild a distribution, while using OptRej=0.
104 //
105 // --------------------------------------------------------------------
106 // Past versions of FOAM: August 2003, v.1.00; September 2003 v.1.01
107 // Adopted starting from FOAM-2.06 by P. Sawicki
108 // --------------------------------------------------------------------
109 // Users of FOAM are kindly requested to cite the following work:
110 // S. Jadach, Computer Physics Communications 152 (2003) 55.
111 //
112 //______________________________________________________________________________
113 
114 #include "TFoam.h"
115 #include "TFoamIntegrand.h"
116 #include "TFoamMaxwt.h"
117 #include "TFoamVect.h"
118 #include "TFoamCell.h"
119 #include "Riostream.h"
120 #include "TH1.h"
121 #include "TRefArray.h"
122 #include "TMethodCall.h"
123 #include "TRandom.h"
124 #include "TMath.h"
125 #include "TInterpreter.h"
126 
127 ClassImp(TFoam);
128 
129 //FFFFFF BoX-FORMATs for nice and flexible outputs
130 #define BXOPE std::cout<<\
131 "FFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFF"<<std::endl<<\
132 "F F"<<std::endl
133 #define BXTXT(text) std::cout<<\
134 "F "<<std::setw(40)<< text <<" F"<<std::endl
135 #define BX1I(name,numb,text) std::cout<<\
136 "F "<<std::setw(10)<<name<<" = "<<std::setw(10)<<numb<<" = " <<std::setw(50)<<text<<" F"<<std::endl
137 #define BX1F(name,numb,text) std::cout<<"F "<<std::setw(10)<<name<<\
138  " = "<<std::setw(15)<<std::setprecision(8)<<numb<<" = "<<std::setw(40)<<text<<" F"<<std::endl
139 #define BX2F(name,numb,err,text) std::cout<<"F "<<std::setw(10)<<name<<\
140 " = "<<std::setw(15)<<std::setprecision(8)<<numb<<" +- "<<std::setw(15)<<std::setprecision(8)<<err<< \
141  " = "<<std::setw(25)<<text<<" F"<<std::endl
142 #define BXCLO std::cout<<\
143 "F F"<<std::endl<<\
144 "FFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFF"<<std::endl
145  //FFFFFF BoX-FORMATs ends here
146 
147 static const Double_t gHigh= 1.0e150;
148 static const Double_t gVlow=-1.0e150;
149 
150 #define SW2 setprecision(7) << std::setw(12)
151 
152 // class to wrap a global function in a TFoamIntegrand function
153 class FoamIntegrandFunction : public TFoamIntegrand {
154 
155 public:
156 
157  typedef Double_t (*FunctionPtr)(Int_t, Double_t*);
158 
159  FoamIntegrandFunction(FunctionPtr func) : fFunc(func) {}
160 
161  virtual ~FoamIntegrandFunction() {}
162 
163  // evaluate the density using the provided function pointer
164  Double_t Density (Int_t nDim, Double_t * x) {
165  return fFunc(nDim,x);
166  }
167 
168 private:
169 
170  FunctionPtr fFunc;
171 
172 };
173 
174 
175 ////////////////////////////////////////////////////////////////////////////////
176 /// Default constructor for streamer, user should not use it.
177 
179  fDim(0), fNCells(0), fRNmax(0),
180  fOptDrive(0), fChat(0), fOptRej(0),
181  fNBin(0), fNSampl(0), fEvPerBin(0),
182  fMaskDiv(0), fInhiDiv(0), fOptPRD(0), fXdivPRD(0),
183  fNoAct(0), fLastCe(0), fCells(0),
184  fMCMonit(0), fMaxWtRej(0), fCellsAct(0), fPrimAcu(0),
185  fHistEdg(0), fHistDbg(0), fHistWt(0),
186  fMCvect(0), fMCwt(0), fRvec(0),
187  fRho(0), fMethodCall(0), fPseRan(0),
188  fNCalls(0), fNEffev(0),
189  fSumWt(0), fSumWt2(0),
190  fSumOve(0), fNevGen(0),
191  fWtMax(0), fWtMin(0),
192  fPrime(0), fMCresult(0), fMCerror(0),
193  fAlpha(0)
194 {
195 }
196 ////////////////////////////////////////////////////////////////////////////////
197 /// User constructor, to be employed by the user
198 
200  fDim(0), fNCells(0), fRNmax(0),
201  fOptDrive(0), fChat(0), fOptRej(0),
202  fNBin(0), fNSampl(0), fEvPerBin(0),
203  fMaskDiv(0), fInhiDiv(0), fOptPRD(0), fXdivPRD(0),
204  fNoAct(0), fLastCe(0), fCells(0),
205  fMCMonit(0), fMaxWtRej(0), fCellsAct(0), fPrimAcu(0),
206  fHistEdg(0), fHistDbg(0), fHistWt(0),
207  fMCvect(0), fMCwt(0), fRvec(0),
208  fRho(0), fMethodCall(0), fPseRan(0),
209  fNCalls(0), fNEffev(0),
210  fSumWt(0), fSumWt2(0),
211  fSumOve(0), fNevGen(0),
212  fWtMax(0), fWtMin(0),
213  fPrime(0), fMCresult(0), fMCerror(0),
214  fAlpha(0)
215 {
216  if(strlen(Name) >129) {
217  Error("TFoam","Name too long %s \n",Name);
218  }
219  fName=Name; // Class name
220  fDate=" Release date: 2005.04.10"; // Release date
221  fVersion= "1.02M"; // Release version
222  fMaskDiv = 0; // Dynamic Mask for cell division, h-cubic
223  fInhiDiv = 0; // Flag allowing to inhibit cell division in certain projection/edge
224  fXdivPRD = 0; // Lists of division values encoded in one vector per direction
225  fCells = 0;
226  fAlpha = 0;
227  fCellsAct = 0;
228  fPrimAcu = 0;
229  fHistEdg = 0;
230  fHistWt = 0;
231  fHistDbg = 0;
232  fDim = 0; // dimension of hyp-cubical space
233  fNCells = 1000; // Maximum number of Cells, is usually re-set
234  fNSampl = 200; // No of sampling when dividing cell
235  fOptPRD = 0; // General Option switch for PRedefined Division, for quick check
236  fOptDrive = 2; // type of Drive =1,2 for TrueVol,Sigma,WtMax
237  fChat = 1; // Chat=0,1,2 chat level in output, Chat=1 normal level
238  fOptRej = 1; // OptRej=0, wted events; OptRej=1, wt=1 events
239  /////////////////////////////////////////////////////////////////////////////
240 
241  fNBin = 8; // binning of edge-histogram in cell exploration
242  fEvPerBin =25; // maximum no. of EFFECTIVE event per bin, =0 option is inactive
243  //------------------------------------------------------
244  fNCalls = 0; // No of function calls
245  fNEffev = 0; // Total no of eff. wt=1 events in build=up
246  fLastCe =-1; // Index of the last cell
247  fNoAct = 0; // No of active cells (used in MC generation)
248  fWtMin = gHigh; // Minimal weight
249  fWtMax = gVlow; // Maximal weight
250  fMaxWtRej =1.10; // Maximum weight in rejection for getting wt=1 events
251  fPseRan = 0; // Initialize private copy of random number generator
252  fMCMonit = 0; // MC efficiency monitoring
253  fRho = 0; // pointer to abstract class providing function to integrate
254  fMCvect = 0;
255  fRvec = 0;
256  fPseRan = 0; // generator of pseudorandom numbers
257  fMethodCall=0; // ROOT's pointer to global distribution function
258 }
259 
260 ////////////////////////////////////////////////////////////////////////////////
261 /// Default destructor
262 
264 {
265  //std::cout<<" DESTRUCTOR entered "<<std::endl;
266  Int_t i;
267 
268  if(fCells!= 0) {
269  for(i=0; i<fNCells; i++) delete fCells[i]; // TFoamCell*[]
270  delete [] fCells;
271  }
272  if (fCellsAct) delete fCellsAct ; // WVE FIX LEAK
273  if (fRvec) delete [] fRvec; //double[]
274  if (fAlpha) delete [] fAlpha; //double[]
275  if (fMCvect) delete [] fMCvect; //double[]
276  if (fPrimAcu) delete [] fPrimAcu; //double[]
277  if (fMaskDiv) delete [] fMaskDiv; //int[]
278  if (fInhiDiv) delete [] fInhiDiv; //int[]
279 
280  if( fXdivPRD!= 0) {
281  for(i=0; i<fDim; i++) delete fXdivPRD[i]; // TFoamVect*[]
282  delete [] fXdivPRD;
283  }
284  if (fMCMonit) delete fMCMonit;
285  if (fHistWt) delete fHistWt;
286 
287  // delete histogram arrays
288  if (fHistEdg) {
289  fHistEdg->Delete();
290  delete fHistEdg;
291  }
292  if (fHistDbg) {
293  fHistDbg->Delete();
294  delete fHistDbg;
295  }
296  // delete function object if it has been created here in SetRhoInt
297  if (fRho && dynamic_cast<FoamIntegrandFunction*>(fRho) ) delete fRho;
298 }
299 
300 
301 ////////////////////////////////////////////////////////////////////////////////
302 /// Copy Constructor NOT IMPLEMENTED (NEVER USED)
303 
304 TFoam::TFoam(const TFoam &From): TObject(From)
305 {
306  Error("TFoam", "COPY CONSTRUCTOR NOT IMPLEMENTED \n");
307 }
308 
309 ////////////////////////////////////////////////////////////////////////////////
310 /// Basic initialization of FOAM invoked by the user. Mandatory!
311 /// ============================================================
312 /// This method starts the process of the cell build-up.
313 /// User must invoke Initialize with two arguments or Initialize without arguments.
314 /// This is done BEFORE generating first MC event and AFTER allocating FOAM object
315 /// and reseting (optionally) its internal parameters/switches.
316 /// The overall operational scheme of the FOAM is the following:
317 ///
318 // <img src="gif/foam_schema2.gif">
319 ///
320 /// This method invokes several other methods:
321 /// ==========================================
322 /// InitCells initializes memory storage for cells and begins exploration process
323 /// from the root cell. The empty cells are allocated/filled using CellFill.
324 /// The procedure Grow which loops over cells, picks up the cell with the biggest
325 /// ``driver integral'', see Comp. Phys. Commun. 152 152 (2003) 55 for explanations,
326 /// with the help of PeekMax procedure. The chosen cell is split using Divide.
327 /// Subsequently, the procedure Explore called by the Divide
328 /// (and by InitCells for the root cell) does the most important
329 /// job in the FOAM object build-up: it performs a small MC run for each
330 /// newly allocated daughter cell.
331 /// Explore calculates how profitable the future split of the cell will be
332 /// and defines the optimal cell division geometry with the help of Carver or Varedu
333 /// procedures, for maximum weight or variance optimization respectively.
334 /// All essential results of the exploration are written into
335 /// the explored cell object. At the very end of the foam build-up,
336 /// Finally, MakeActiveList is invoked to create a list of pointers to
337 /// all active cells, for the purpose of the quick access during the MC generation.
338 /// The procedure Explore employs MakeAlpha to generate random coordinates
339 /// inside a given cell with the uniform distribution.
340 /// The above sequence of the procedure calls is depicted in the following figure:
341 ///
342 // <img src="gif/foam_Initialize_schema.gif">
343 
345 {
346  SetPseRan(PseRan);
347  SetRho(fun);
348  Initialize();
349 }
350 
351 ////////////////////////////////////////////////////////////////////////////////
352 /// Basic initialization of FOAM invoked by the user.
353 /// IMPORTANT: Random number generator and the distribution object has to be
354 /// provided using SetPseRan and SetRho prior to invoking this initializator!
355 
357 {
358  Bool_t addStatus = TH1::AddDirectoryStatus();
360  Int_t i;
361 
362  if(fChat>0){
363  BXOPE;
364  BXTXT("****************************************");
365  BXTXT("****** TFoam::Initialize ******");
366  BXTXT("****************************************");
367  BXTXT(fName);
368  BX1F(" Version",fVersion, fDate);
369  BX1I(" kDim",fDim, " Dimension of the hyper-cubical space ");
370  BX1I(" nCells",fNCells, " Requested number of Cells (half of them active) ");
371  BX1I(" nSampl",fNSampl, " No of MC events in exploration of a cell ");
372  BX1I(" nBin",fNBin, " No of bins in histograms, MC exploration of cell ");
373  BX1I(" EvPerBin",fEvPerBin, " Maximum No effective_events/bin, MC exploration ");
374  BX1I(" OptDrive",fOptDrive, " Type of Driver =1,2 for Sigma,WtMax ");
375  BX1I(" OptRej",fOptRej, " MC rejection on/off for OptRej=0,1 ");
376  BX1F(" MaxWtRej",fMaxWtRej, " Maximum wt in rejection for wt=1 evts");
377  BXCLO;
378  }
379 
380  if(fPseRan==0) Error("Initialize", "Random number generator not set \n");
381  if(fRho==0 && fMethodCall==0 ) Error("Initialize", "Distribution function not set \n");
382  if(fDim==0) Error("Initialize", "Zero dimension not allowed \n");
383 
384  /////////////////////////////////////////////////////////////////////////
385  // ALLOCATE SMALL LISTS //
386  // it is done globally, not for each cell, to save on allocation time //
387  /////////////////////////////////////////////////////////////////////////
388  fRNmax= fDim+1;
389  fRvec = new Double_t[fRNmax]; // Vector of random numbers
390  if(fRvec==0) Error("Initialize", "Cannot initialize buffer fRvec \n");
391 
392  if(fDim>0){
393  fAlpha = new Double_t[fDim]; // sum<1 for internal parametrization of the simplex
394  if(fAlpha==0) Error("Initialize", "Cannot initialize buffer fAlpha \n" );
395  }
396  fMCvect = new Double_t[fDim]; // vector generated in the MC run
397  if(fMCvect==0) Error("Initialize", "Cannot initialize buffer fMCvect \n" );
398 
399  //====== List of directions inhibited for division
400  if(fInhiDiv == 0){
401  fInhiDiv = new Int_t[fDim];
402  for(i=0; i<fDim; i++) fInhiDiv[i]=0;
403  }
404  //====== Dynamic mask used in Explore for edge determination
405  if(fMaskDiv == 0){
406  fMaskDiv = new Int_t[fDim];
407  for(i=0; i<fDim; i++) fMaskDiv[i]=1;
408  }
409  //====== List of predefined division values in all directions (initialized as empty)
410  if(fXdivPRD == 0){
411  fXdivPRD = new TFoamVect*[fDim];
412  for(i=0; i<fDim; i++) fXdivPRD[i]=0; // Artificially extended beyond fDim
413  }
414  //====== Initialize list of histograms
415  fHistWt = new TH1D("HistWt","Histogram of MC weight",100,0.0, 1.5*fMaxWtRej); // MC weight
416  fHistEdg = new TObjArray(fDim); // Initialize list of histograms
417  TString hname;
418  TString htitle;
419  for(i=0;i<fDim;i++){
420  hname=fName+TString("_HistEdge_");
421  hname+=i;
422  htitle=TString("Edge Histogram No. ");
423  htitle+=i;
424  //std::cout<<"i= "<<i<<" hname= "<<hname<<" htitle= "<<htitle<<std::endl;
425  (*fHistEdg)[i] = new TH1D(hname.Data(),htitle.Data(),fNBin,0.0, 1.0); // Initialize histogram for each edge
426  ((TH1D*)(*fHistEdg)[i])->Sumw2();
427  }
428  //====== extra histograms for debug purposes
429  fHistDbg = new TObjArray(fDim); // Initialize list of histograms
430  for(i=0;i<fDim;i++){
431  hname=fName+TString("_HistDebug_");
432  hname+=i;
433  htitle=TString("Debug Histogram ");
434  htitle+=i;
435  (*fHistDbg)[i] = new TH1D(hname.Data(),htitle.Data(),fNBin,0.0, 1.0); // Initialize histogram for each edge
436  }
437 
438  // ||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||| //
439  // BUILD-UP of the FOAM //
440  // ||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||| //
441  //
442  // Define and explore root cell(s)
443  InitCells();
444  // PrintCells(); std::cout<<" ===== after InitCells ====="<<std::endl;
445  Grow();
446  // PrintCells(); std::cout<<" ===== after Grow ====="<<std::endl;
447 
448  MakeActiveList(); // Final Preperations for the M.C. generation
449 
450  // Preperations for the M.C. generation
451  fSumWt = 0.0; // M.C. generation sum of Wt
452  fSumWt2 = 0.0; // M.C. generation sum of Wt**2
453  fSumOve = 0.0; // M.C. generation sum of overweighted
454  fNevGen = 0.0; // M.C. generation sum of 1d0
455  fWtMax = gVlow; // M.C. generation maximum wt
456  fWtMin = gHigh; // M.C. generation minimum wt
457  fMCresult=fCells[0]->GetIntg(); // M.C. Value of INTEGRAL,temporary assignment
458  fMCresult=fCells[0]->GetIntg(); // M.C. Value of INTEGRAL,temporary assignment
459  fMCerror =fCells[0]->GetIntg(); // M.C. Value of ERROR ,temporary assignment
460  fMCMonit = new TFoamMaxwt(5.0,1000); // monitoring M.C. efficiency
461  //
462  if(fChat>0){
463  Double_t driver = fCells[0]->GetDriv();
464  BXOPE;
465  BXTXT("*** TFoam::Initialize FINISHED!!! ***");
466  BX1I(" nCalls",fNCalls, "Total number of function calls ");
467  BX1F(" XPrime",fPrime, "Primary total integral ");
468  BX1F(" XDiver",driver, "Driver total integral ");
469  BX1F(" mcResult",fMCresult,"Estimate of the true MC Integral ");
470  BXCLO;
471  }
472  if(fChat==2) PrintCells();
473  TH1::AddDirectory(addStatus);
474 } // Initialize
475 
476 ////////////////////////////////////////////////////////////////////////////////
477 /// Internal subprogram used by Initialize.
478 /// It initializes "root part" of the FOAM of the tree of cells.
479 
481 {
482  Int_t i;
483 
484  fLastCe =-1; // Index of the last cell
485  if(fCells!= 0) {
486  for(i=0; i<fNCells; i++) delete fCells[i];
487  delete [] fCells;
488  }
489  //
490  fCells = new TFoamCell*[fNCells];
491  for(i=0;i<fNCells;i++){
492  fCells[i]= new TFoamCell(fDim); // Allocate BIG list of cells
493  fCells[i]->SetSerial(i);
494  }
495  if(fCells==0) Error("InitCells", "Cannot initialize CELLS \n" );
496 
497  /////////////////////////////////////////////////////////////////////////////
498  // Single Root Hypercube //
499  /////////////////////////////////////////////////////////////////////////////
500  CellFill(1, 0); // 0-th cell ACTIVE
501 
502  // Exploration of the root cell(s)
503  for(Long_t iCell=0; iCell<=fLastCe; iCell++){
504  Explore( fCells[iCell] ); // Exploration of root cell(s)
505  }
506 }//InitCells
507 
508 ////////////////////////////////////////////////////////////////////////////////
509 /// Internal subprogram used by Initialize.
510 /// It initializes content of the newly allocated active cell.
511 
513 {
514  TFoamCell *cell;
515  if (fLastCe==fNCells){
516  Error( "CellFill", "Too many cells\n");
517  }
518  fLastCe++; // 0-th cell is the first
519  if (Status==1) fNoAct++;
520 
521  cell = fCells[fLastCe];
522 
523  cell->Fill(Status, parent, 0, 0);
524 
525  cell->SetBest( -1); // pointer for planning division of the cell
526  cell->SetXdiv(0.5); // factor for division
527  Double_t xInt2,xDri2;
528  if(parent!=0){
529  xInt2 = 0.5*parent->GetIntg();
530  xDri2 = 0.5*parent->GetDriv();
531  cell->SetIntg(xInt2);
532  cell->SetDriv(xDri2);
533  }else{
534  cell->SetIntg(0.0);
535  cell->SetDriv(0.0);
536  }
537  return fLastCe;
538 }
539 
540 ////////////////////////////////////////////////////////////////////////////////
541 /// Internal subprogram used by Initialize.
542 /// It explores newly defined cell with help of special short MC sampling.
543 /// As a result, estimates of true and drive volume is defined/determined
544 /// Average and dispersion of the weight distribution will is found along
545 /// each edge and the best edge (minimum dispersion, best maximum weight)
546 /// is memorized for future use.
547 /// The optimal division point for eventual future cell division is
548 /// determined/recorded. Recorded are also minimum and maximum weight etc.
549 /// The volume estimate in all (inactive) parent cells is updated.
550 /// Note that links to parents and initial volume = 1/2 parent has to be
551 /// already defined prior to calling this routine.
552 
554 {
555  Double_t wt, dx, xBest=0, yBest=0;
556  Double_t intOld, driOld;
557 
558  Long_t iev;
559  Double_t nevMC;
560  Int_t i, j, k;
561  Int_t nProj, kBest;
562  Double_t ceSum[5], xproj;
563 
564  TFoamVect cellSize(fDim);
565  TFoamVect cellPosi(fDim);
566 
567  cell->GetHcub(cellPosi,cellSize);
568 
569  TFoamCell *parent;
570 
571  Double_t *xRand = new Double_t[fDim];
572 
573  Double_t *volPart=0;
574 
575  cell->CalcVolume();
576  dx = cell->GetVolume();
577  intOld = cell->GetIntg(); //memorize old values,
578  driOld = cell->GetDriv(); //will be needed for correcting parent cells
579 
580 
581  /////////////////////////////////////////////////////
582  // Special Short MC sampling to probe cell //
583  /////////////////////////////////////////////////////
584  ceSum[0]=0;
585  ceSum[1]=0;
586  ceSum[2]=0;
587  ceSum[3]=gHigh; //wtmin
588  ceSum[4]=gVlow; //wtmax
589  //
590  for(i=0;i<fDim;i++) ((TH1D *)(*fHistEdg)[i])->Reset(); // Reset histograms
591  fHistWt->Reset();
592  //
593  // ||||||||||||||||||||||||||BEGIN MC LOOP|||||||||||||||||||||||||||||
594  Double_t nevEff=0.;
595  for(iev=0;iev<fNSampl;iev++){
596  MakeAlpha(); // generate uniformly vector inside hypercube
597 
598  if(fDim>0){
599  for(j=0; j<fDim; j++)
600  xRand[j]= cellPosi[j] +fAlpha[j]*(cellSize[j]);
601  }
602 
603  wt=dx*Eval(xRand);
604 
605  nProj = 0;
606  if(fDim>0) {
607  for(k=0; k<fDim; k++) {
608  xproj =fAlpha[k];
609  ((TH1D *)(*fHistEdg)[nProj])->Fill(xproj,wt);
610  nProj++;
611  }
612  }
613  //
614  fNCalls++;
615  ceSum[0] += wt; // sum of weights
616  ceSum[1] += wt*wt; // sum of weights squared
617  ceSum[2]++; // sum of 1
618  if (ceSum[3]>wt) ceSum[3]=wt; // minimum weight;
619  if (ceSum[4]<wt) ceSum[4]=wt; // maximum weight
620  // test MC loop exit condition
621  nevEff = ceSum[0]*ceSum[0]/ceSum[1];
622  if( nevEff >= fNBin*fEvPerBin) break;
623  } // ||||||||||||||||||||||||||END MC LOOP|||||||||||||||||||||||||||||
624  //------------------------------------------------------------------
625  //--- predefine logics of searching for the best division edge ---
626  for(k=0; k<fDim;k++){
627  fMaskDiv[k] =1; // default is all
628  if( fInhiDiv[k]==1) fMaskDiv[k] =0; // inhibit some...
629  }
630  // Note that predefined division below overrule inhibition above
631  kBest=-1;
632  Double_t rmin,rmax,rdiv;
633  if(fOptPRD) { // quick check
634  for(k=0; k<fDim; k++) {
635  rmin= cellPosi[k];
636  rmax= cellPosi[k] +cellSize[k];
637  if( fXdivPRD[k] != 0) {
638  Int_t n= (fXdivPRD[k])->GetDim();
639  for(j=0; j<n; j++) {
640  rdiv=(*fXdivPRD[k])[j];
641  // check predefined divisions is available in this cell
642  if( (rmin +1e-99 <rdiv) && (rdiv< rmax -1e-99)) {
643  kBest=k;
644  xBest= (rdiv-cellPosi[k])/cellSize[k] ;
645  goto ee05;
646  }
647  }
648  }
649  }//k
650  }
651  ee05:
652  /////////////////////////////////////////////////////////////////////////////
653 
654  fNEffev += (Long_t)nevEff;
655  nevMC = ceSum[2];
656  Double_t intTrue = ceSum[0]/(nevMC+0.000001);
657  Double_t intDriv=0.;
658  Double_t intPrim=0.;
659 
660  switch(fOptDrive){
661  case 1: // VARIANCE REDUCTION
662  if(kBest == -1) Varedu(ceSum,kBest,xBest,yBest); // determine the best edge,
663  //intDriv =sqrt( ceSum[1]/nevMC -intTrue*intTrue ); // Older ansatz, numerically not bad
664  intDriv =sqrt(ceSum[1]/nevMC) -intTrue; // Foam build-up, sqrt(<w**2>) -<w>
665  intPrim =sqrt(ceSum[1]/nevMC); // MC gen. sqrt(<w**2>) =sqrt(<w>**2 +sigma**2)
666  break;
667  case 2: // WTMAX REDUCTION
668  if(kBest == -1) Carver(kBest,xBest,yBest); // determine the best edge
669  intDriv =ceSum[4] -intTrue; // Foam build-up, wtmax-<w>
670  intPrim =ceSum[4]; // MC generation, wtmax!
671  break;
672  default:
673  Error("Explore", "Wrong fOptDrive = \n" );
674  }//switch
675  //=================================================================================
676  //hist_Neff_distrib.Fill( fLastCe/2.0+0.01, nevEff+0.01); //
677  //hist_kBest_distrib.Fill( kBest+0.50, 1.0 ); // debug
678  //hist_xBest_distrib.Fill( xBest+0.01, 1.0 ); // debug
679  //=================================================================================
680  cell->SetBest(kBest);
681  cell->SetXdiv(xBest);
682  cell->SetIntg(intTrue);
683  cell->SetDriv(intDriv);
684  cell->SetPrim(intPrim);
685  // correct/update integrals in all parent cells to the top of the tree
686  Double_t parIntg, parDriv;
687  for(parent = cell->GetPare(); parent!=0; parent = parent->GetPare()){
688  parIntg = parent->GetIntg();
689  parDriv = parent->GetDriv();
690  parent->SetIntg( parIntg +intTrue -intOld );
691  parent->SetDriv( parDriv +intDriv -driOld );
692  }
693  delete [] volPart;
694  delete [] xRand;
695  //cell->Print();
696 } // TFoam::Explore
697 
698 ////////////////////////////////////////////////////////////////////////////////
699 /// Internal subrogram used by Initialize.
700 /// In determines the best edge candidate and the position of the cell division plane
701 /// in case of the variance reduction for future cell division,
702 /// using results of the MC exploration run stored in fHistEdg
703 
704 void TFoam::Varedu(Double_t ceSum[5], Int_t &kBest, Double_t &xBest, Double_t &yBest)
705 {
706  Double_t nent = ceSum[2];
707  Double_t swAll = ceSum[0];
708  Double_t sswAll = ceSum[1];
709  Double_t ssw = sqrt(sswAll)/sqrt(nent);
710  //
711  Double_t swIn,swOut,sswIn,sswOut,xLo,xUp;
712  kBest =-1;
713  xBest =0.5;
714  yBest =1.0;
715  Double_t maxGain=0.0;
716  // Now go over all projections kProj
717  for(Int_t kProj=0; kProj<fDim; kProj++) {
718  if( fMaskDiv[kProj]) {
719  // initialize search over bins
720  Double_t sigmIn =0.0; Double_t sigmOut =0.0;
721  Double_t sswtBest = gHigh;
722  Double_t gain =0.0;
723  Double_t xMin=0.0; Double_t xMax=0.0;
724  // Double loop over all pairs jLo<jUp
725  for(Int_t jLo=1; jLo<=fNBin; jLo++) {
726  Double_t aswIn=0; Double_t asswIn=0;
727  for(Int_t jUp=jLo; jUp<=fNBin;jUp++) {
728  aswIn += ((TH1D *)(*fHistEdg)[kProj])->GetBinContent(jUp);
729  asswIn += Sqr(((TH1D *)(*fHistEdg)[kProj])->GetBinError( jUp));
730  xLo=(jLo-1.0)/fNBin;
731  xUp=(jUp*1.0)/fNBin;
732  swIn = aswIn/nent;
733  swOut = (swAll-aswIn)/nent;
734  sswIn = sqrt(asswIn) /sqrt(nent*(xUp-xLo)) *(xUp-xLo);
735  sswOut= sqrt(sswAll-asswIn)/sqrt(nent*(1.0-xUp+xLo)) *(1.0-xUp+xLo);
736  if( (sswIn+sswOut) < sswtBest) {
737  sswtBest = sswIn+sswOut;
738  gain = ssw-sswtBest;
739  sigmIn = sswIn -swIn; // Debug
740  sigmOut = sswOut-swOut; // Debug
741  xMin = xLo;
742  xMax = xUp;
743  }
744  }//jUp
745  }//jLo
746  Int_t iLo = (Int_t) (fNBin*xMin);
747  Int_t iUp = (Int_t) (fNBin*xMax);
748  //----------DEBUG printout
749  //std::cout<<"@@@@@ xMin xMax = "<<xMin <<" "<<xMax<<" iLo= "<<iLo<<" iUp= "<<iUp;
750  //std::cout<<" sswtBest/ssw= "<<sswtBest/ssw<<" Gain/ssw= "<< Gain/ssw<<std::endl;
751  //----------DEBUG auxilary Plot
752  for(Int_t iBin=1;iBin<=fNBin;iBin++) {
753  if( ((iBin-0.5)/fNBin > xMin) && ((iBin-0.5)/fNBin < xMax) ){
754  ((TH1D *)(*fHistDbg)[kProj])->SetBinContent(iBin,sigmIn/(xMax-xMin));
755  } else {
756  ((TH1D *)(*fHistDbg)[kProj])->SetBinContent(iBin,sigmOut/(1-xMax+xMin));
757  }
758  }
759  if(gain>=maxGain) {
760  maxGain=gain;
761  kBest=kProj; // <--- !!!!! The best edge
762  xBest=xMin;
763  yBest=xMax;
764  if(iLo == 0) xBest=yBest; // The best division point
765  if(iUp == fNBin) yBest=xBest; // this is not really used
766  }
767  }
768  } //kProj
769  //----------DEBUG printout
770  //std::cout<<"@@@@@@@>>>>> kBest= "<<kBest<<" maxGain/ssw= "<< maxGain/ssw<<std::endl;
771  if( (kBest >= fDim) || (kBest<0) ) Error("Varedu", "Something wrong with kBest - kBest = %d dim = %d\n",kBest,fDim);
772 } //TFoam::Varedu
773 
774 ////////////////////////////////////////////////////////////////////////////////
775 /// Internal subrogram used by Initialize.
776 /// Determines the best edge-candidate and the position of the division plane
777 /// for the future cell division, in the case of the optimization of the maximum weight.
778 /// It exploits results of the cell MC exploration run stored in fHistEdg.
779 
780 void TFoam::Carver(Int_t &kBest, Double_t &xBest, Double_t &yBest)
781 {
782  Int_t kProj,iBin;
783  Double_t carve,carvTot,carvMax,carvOne,binMax,binTot;
784  Int_t jLow,jUp,iLow,iUp;
785  Double_t theBin;
786  // Int_t jDivi; // TEST
787  Int_t j;
788 
789  Double_t *bins = new Double_t[fNBin]; // bins of histogram for single PROJECTION
790  if(bins==0) Error("Carver", "Cannot initialize buffer Bins \n" );
791 
792  kBest =-1;
793  xBest =0.5;
794  yBest =1.0;
795  carvMax = gVlow;
796  // primMax = gVlow;
797  for(kProj=0; kProj<fDim; kProj++)
798  if( fMaskDiv[kProj] ) {
799  //if( kProj==1 ){
800  //std::cout<<"==================== Carver histogram: kProj ="<<kProj<<"==================="<<std::endl;
801  //((TH1D *)(*fHistEdg)[kProj])->Print("all");
802  binMax = gVlow;
803  for(iBin=0; iBin<fNBin;iBin++){
804  bins[iBin]= ((TH1D *)(*fHistEdg)[kProj])->GetBinContent(iBin+1);
805  binMax = TMath::Max( binMax, bins[iBin]); // Maximum content/bin
806  }
807  if(binMax < 0 ) { //case of empty cell
808  delete [] bins;
809  return;
810  }
811  carvTot = 0.0;
812  binTot = 0.0;
813  for(iBin=0;iBin<fNBin;iBin++){
814  carvTot = carvTot + (binMax-bins[iBin]); // Total Carve (more stable)
815  binTot +=bins[iBin];
816  }
817  // primTot = binMax*fNBin;
818  //std::cout <<"Carver: CarvTot "<<CarvTot<< " primTot "<<primTot<<std::endl;
819  jLow =0;
820  jUp =fNBin-1;
821  carvOne = gVlow;
822  Double_t yLevel = gVlow;
823  for(iBin=0; iBin<fNBin;iBin++) {
824  theBin = bins[iBin];
825  //----- walk to the left and find first bin > theBin
826  iLow = iBin;
827  for(j=iBin; j>-1; j-- ) {
828  if(theBin< bins[j]) break;
829  iLow = j;
830  }
831  //iLow = iBin;
832  //if(iLow>0) while( (theBin >= bins[iLow-1])&&(iLow >0) ){iLow--;} // horror!!!
833  //------ walk to the right and find first bin > theBin
834  iUp = iBin;
835  for(j=iBin; j<fNBin; j++) {
836  if(theBin< bins[j]) break;
837  iUp = j;
838  }
839  //iUp = iBin;
840  //if(iUp<fNBin-1) while( (theBin >= bins[iUp+1])&&( iUp<fNBin-1 ) ){iUp++;} // horror!!!
841  //
842  carve = (iUp-iLow+1)*(binMax-theBin);
843  if( carve > carvOne) {
844  carvOne = carve;
845  jLow = iLow;
846  jUp = iUp;
847  yLevel = theBin;
848  }
849  }//iBin
850  if( carvTot > carvMax) {
851  carvMax = carvTot;
852  //primMax = primTot;
853  //std::cout <<"Carver: primMax "<<primMax<<std::endl;
854  kBest = kProj; // Best edge
855  xBest = ((Double_t)(jLow))/fNBin;
856  yBest = ((Double_t)(jUp+1))/fNBin;
857  if(jLow == 0 ) xBest = yBest;
858  if(jUp == fNBin-1) yBest = xBest;
859  // division ratio in units of 1/fNBin, testing
860  // jDivi = jLow;
861  // if(jLow == 0 ) jDivi=jUp+1;
862  }
863  //====== extra histograms for debug purposes
864  //std::cout<<"kProj= "<<kProj<<" jLow= "<<jLow<<" jUp= "<<jUp<<std::endl;
865  for(iBin=0; iBin<fNBin; iBin++)
866  ((TH1D *)(*fHistDbg)[kProj])->SetBinContent(iBin+1,binMax);
867  for(iBin=jLow; iBin<jUp+1; iBin++)
868  ((TH1D *)(*fHistDbg)[kProj])->SetBinContent(iBin+1,yLevel);
869  }//kProj
870  if( (kBest >= fDim) || (kBest<0) ) Error("Carver", "Something wrong with kBest - kBest = %d dim = %d\n",kBest,fDim);
871  delete [] bins;
872 } //TFoam::Carver
873 
874 ////////////////////////////////////////////////////////////////////////////////
875 /// Internal subrogram used by Initialize.
876 /// Provides random vector Alpha 0< Alpha(i) < 1
877 
879 {
880  Int_t k;
881  if(fDim<1) return;
882 
883  // simply generate and load kDim uniform random numbers
884  fPseRan->RndmArray(fDim,fRvec); // kDim random numbers needed
885  for(k=0; k<fDim; k++) fAlpha[k] = fRvec[k];
886 } //MakeAlpha
887 
888 
889 ////////////////////////////////////////////////////////////////////////////////
890 /// Internal subrogram used by Initialize.
891 /// It grow new cells by the binary division process.
892 
894 {
895  Long_t iCell;
896  TFoamCell* newCell;
897 
898  while ( (fLastCe+2) < fNCells ) { // this condition also checked inside Divide
899  iCell = PeekMax(); // peek up cell with maximum driver integral
900  if( (iCell<0) || (iCell>fLastCe) ) Error("Grow", "Wrong iCell \n");
901  newCell = fCells[iCell];
902 
903  if(fLastCe !=0) {
904  Int_t kEcho=10;
905  if(fLastCe>=10000) kEcho=100;
906  if( (fLastCe%kEcho)==0) {
907  if (fChat>0) {
908  if(fDim<10)
909  std::cout<<fDim<<std::flush;
910  else
911  std::cout<<"."<<std::flush;
912  if( (fLastCe%(100*kEcho))==0) std::cout<<"|"<<fLastCe<<std::endl<<std::flush;
913  }
914  }
915  }
916  if( Divide( newCell )==0) break; // and divide it into two
917  }
918  if (fChat>0) {
919  std::cout<<std::endl<<std::flush;
920  }
921  CheckAll(0); // set arg=1 for more info
922 }// Grow
923 
924 ////////////////////////////////////////////////////////////////////////////////
925 /// Internal subprogram used by Initialize.
926 /// It finds cell with maximal driver integral for the purpose of the division.
927 
929 {
930  Long_t i;
931  Long_t iCell = -1;
932  Double_t drivMax, driv;
933 
934  drivMax = gVlow;
935  for(i=0; i<=fLastCe; i++) {//without root
936  if( fCells[i]->GetStat() == 1 ) {
937  driv = TMath::Abs( fCells[i]->GetDriv());
938  //std::cout<<"PeekMax: Driv = "<<driv<<std::endl;
939  if(driv > drivMax) {
940  drivMax = driv;
941  iCell = i;
942  }
943  }
944  }
945  // std::cout << 'TFoam_PeekMax: iCell=' << iCell << std::endl;
946  if (iCell == -1)
947  std::cout << "STOP in TFoam::PeekMax: not found iCell=" << iCell << std::endl;
948  return(iCell);
949 } // TFoam_PeekMax
950 
951 ////////////////////////////////////////////////////////////////////////////////
952 /// Internal subrogram used by Initialize.
953 /// It divides cell iCell into two daughter cells.
954 /// The iCell is retained and tagged as inactive, daughter cells are appended
955 /// at the end of the buffer.
956 /// New vertex is added to list of vertices.
957 /// List of active cells is updated, iCell removed, two daughters added
958 /// and their properties set with help of MC sampling (TFoam_Explore)
959 /// Returns Code RC=-1 of buffer limit is reached, fLastCe=fnBuf.
960 
962 {
963  Int_t kBest;
964 
965  if(fLastCe+1 >= fNCells) Error("Divide", "Buffer limit is reached, fLastCe=fnBuf \n");
966 
967  cell->SetStat(0); // reset cell as inactive
968  fNoAct--;
969 
970  kBest = cell->GetBest();
971  if( kBest<0 || kBest>=fDim ) Error("Divide", "Wrong kBest \n");
972 
973  //////////////////////////////////////////////////////////////////
974  // define two daughter cells (active) //
975  //////////////////////////////////////////////////////////////////
976 
977  Int_t d1 = CellFill(1, cell);
978  Int_t d2 = CellFill(1, cell);
979  cell->SetDau0((fCells[d1]));
980  cell->SetDau1((fCells[d2]));
981  Explore( (fCells[d1]) );
982  Explore( (fCells[d2]) );
983  return 1;
984 } // TFoam_Divide
985 
986 
987 ////////////////////////////////////////////////////////////////////////////////
988 /// Internal subrogram used by Initialize.
989 /// It finds out number of active cells fNoAct,
990 /// creates list of active cell fCellsAct and primary cumulative fPrimAcu.
991 /// They are used during the MC generation to choose randomly an active cell.
992 
994 {
995  Long_t n, iCell;
996  Double_t sum;
997  // flush previous result
998  if(fPrimAcu != 0) delete [] fPrimAcu;
999  if(fCellsAct != 0) delete fCellsAct;
1000 
1001  // Allocate tables of active cells
1002  fCellsAct = new TRefArray();
1003 
1004  // Count Active cells and find total Primary
1005  // Fill-in tables of active cells
1006 
1007  fPrime = 0.0; n = 0;
1008  for(iCell=0; iCell<=fLastCe; iCell++) {
1009  if (fCells[iCell]->GetStat()==1) {
1010  fPrime += fCells[iCell]->GetPrim();
1011  fCellsAct->Add(fCells[iCell]);
1012  n++;
1013  }
1014  }
1015 
1016  if(fNoAct != n) Error("MakeActiveList", "Wrong fNoAct \n" );
1017  if(fPrime == 0.) Error("MakeActiveList", "Integrand function is zero \n" );
1018 
1019  fPrimAcu = new Double_t[fNoAct]; // cumulative primary for MC generation
1020  if( fCellsAct==0 || fPrimAcu==0 ) Error("MakeActiveList", "Cant allocate fCellsAct or fPrimAcu \n");
1021 
1022  sum =0.0;
1023  for(iCell=0; iCell<fNoAct; iCell++) {
1024  sum = sum + ( (TFoamCell *) (fCellsAct->At(iCell)) )->GetPrim()/fPrime;
1025  fPrimAcu[iCell]=sum;
1026  }
1027 
1028 } //MakeActiveList
1029 
1030 ////////////////////////////////////////////////////////////////////////////////
1031 /// User may optionally reset random number generator using this method
1032 /// Usually it is done when FOAM object is restored from the disk.
1033 /// IMPORTANT: this method deletes existing random number generator registered in the FOAM object.
1034 /// In particular such an object is created by the streamer during the disk-read operation.
1035 
1037 {
1038  if(fPseRan) {
1039  Info("ResetPseRan", "Resetting random number generator \n");
1040  delete fPseRan;
1041  }
1042  SetPseRan(PseRan);
1043 }
1044 
1045 ////////////////////////////////////////////////////////////////////////////////
1046 /// User may use this method to set the distribution object
1047 
1049 {
1050  if (fun)
1051  fRho=fun;
1052  else
1053  Error("SetRho", "Bad function \n" );
1054 }
1055 ////////////////////////////////////////////////////////////////////////////////
1056 /// User may use this method to set the distribution object as a global function pointer
1057 /// (and not as an interpreted function).
1058 
1060 {
1061  // This is needed for both AClic and Cling
1062  if (fun) {
1063  // delete function object if it has been created here in SetRho
1064  if (fRho && dynamic_cast<FoamIntegrandFunction*>(fRho) ) delete fRho;
1065  fRho= new FoamIntegrandFunction(fun);
1066  } else
1067  Error("SetRho", "Bad function \n" );
1068 }
1069 
1070 ////////////////////////////////////////////////////////////////////////////////
1071 /// User may optionally reset the distribution using this method
1072 /// Usually it is done when FOAM object is restored from the disk.
1073 /// IMPORTANT: this method deletes existing distribution object registered in the FOAM object.
1074 /// In particular such an object is created by the streamer diring the disk-read operation.
1075 /// This method is used only in very special cases, because the distribution in most cases
1076 /// should be "owned" by the FOAM object and should not be replaced by another one after initialization.
1077 
1079 {
1080  if(fRho) {
1081  Info("ResetRho", "!!! Resetting distribution function !!!\n");
1082  delete fRho;
1083  }
1084  SetRho(fun);
1085 }
1086 
1087 ////////////////////////////////////////////////////////////////////////////////
1088 /// User may use this to set pointer to the global function (not descending
1089 /// from TFoamIntegrand) serving as a distribution for FOAM.
1090 /// It is useful for simple interactive applications.
1091 /// Note that persistency for FOAM object will not work in the case of such
1092 /// a distribution.
1093 
1094 void TFoam::SetRhoInt( void * fun)
1095 {
1096 
1097  const char *namefcn = gCling->Getp2f2funcname(fun); //name of integrand function
1098  if(namefcn) {
1099  fMethodCall=new TMethodCall();
1100  fMethodCall->InitWithPrototype(namefcn, "Int_t, Double_t *");
1101  }
1102  fRho=0;
1103 }
1104 
1105 ////////////////////////////////////////////////////////////////////////////////
1106 /// Internal subprogram.
1107 /// Evaluates distribution to be generated.
1108 
1110 {
1111  Double_t result;
1112 
1113  if(!fRho) { //interactive mode
1114  Long_t paramArr[3];
1115  paramArr[0]=(Long_t)fDim;
1116  paramArr[1]=(Long_t)xRand;
1117  fMethodCall->SetParamPtrs(paramArr);
1118  fMethodCall->Execute(result);
1119  } else { //compiled mode
1120  result=fRho->Density(fDim,xRand);
1121  }
1122 
1123  return result;
1124 }
1125 
1126 ////////////////////////////////////////////////////////////////////////////////
1127 /// Internal subprogram.
1128 /// Return randomly chosen active cell with probability equal to its
1129 /// contribution into total driver integral using interpolation search.
1130 
1132 {
1133  Long_t lo, hi, hit;
1134  Double_t fhit, flo, fhi;
1135  Double_t random;
1136 
1137  random=fPseRan->Rndm();
1138  lo = 0; hi =fNoAct-1;
1139  flo = fPrimAcu[lo]; fhi=fPrimAcu[hi];
1140  while(lo+1<hi) {
1141  hit = lo + (Int_t)( (hi-lo)*(random-flo)/(fhi-flo)+0.5);
1142  if (hit<=lo)
1143  hit = lo+1;
1144  else if(hit>=hi)
1145  hit = hi-1;
1146  fhit=fPrimAcu[hit];
1147  if (fhit>random) {
1148  hi = hit;
1149  fhi = fhit;
1150  } else {
1151  lo = hit;
1152  flo = fhit;
1153  }
1154  }
1155  if (fPrimAcu[lo]>random)
1156  pCell = (TFoamCell *) fCellsAct->At(lo);
1157  else
1158  pCell = (TFoamCell *) fCellsAct->At(hi);
1159 } // TFoam::GenerCel2
1160 
1161 
1162 ////////////////////////////////////////////////////////////////////////////////
1163 /// User subprogram.
1164 /// It generates randomly point/vector according to user-defined distribution.
1165 /// Prior initialization with help of Initialize() is mandatory.
1166 /// Generated MC point/vector is available using GetMCvect and the MC weight with GetMCwt.
1167 /// MC point is generated with wt=1 or with variable weight, see OptRej switch.
1168 
1170 {
1171  Int_t j;
1172  Double_t wt,dx,mcwt;
1173  TFoamCell *rCell;
1174  //
1175  //********************** MC LOOP STARS HERE **********************
1176 ee0:
1177  GenerCel2(rCell); // choose randomly one cell
1178 
1179  MakeAlpha();
1180 
1181  TFoamVect cellPosi(fDim); TFoamVect cellSize(fDim);
1182  rCell->GetHcub(cellPosi,cellSize);
1183  for(j=0; j<fDim; j++)
1184  fMCvect[j]= cellPosi[j] +fAlpha[j]*cellSize[j];
1185  dx = rCell->GetVolume(); // Cartesian volume of the Cell
1186  // weight average normalized to PRIMARY integral over the cell
1187 
1188  wt=dx*Eval(fMCvect);
1189 
1190  mcwt = wt / rCell->GetPrim(); // PRIMARY controls normalization
1191  fNCalls++;
1192  fMCwt = mcwt;
1193  // accumulation of statistics for the main MC weight
1194  fSumWt += mcwt; // sum of Wt
1195  fSumWt2 += mcwt*mcwt; // sum of Wt**2
1196  fNevGen++; // sum of 1d0
1197  fWtMax = TMath::Max(fWtMax, mcwt); // maximum wt
1198  fWtMin = TMath::Min(fWtMin, mcwt); // minimum wt
1199  fMCMonit->Fill(mcwt);
1200  fHistWt->Fill(mcwt,1.0); // histogram
1201  //******* Optional rejection ******
1202  if(fOptRej == 1) {
1203  Double_t random;
1204  random=fPseRan->Rndm();
1205  if( fMaxWtRej*random > fMCwt) goto ee0; // Wt=1 events, internal rejection
1206  if( fMCwt<fMaxWtRej ) {
1207  fMCwt = 1.0; // normal Wt=1 event
1208  } else {
1209  fMCwt = fMCwt/fMaxWtRej; // weight for overweighted events! kept for debug
1210  fSumOve += fMCwt-fMaxWtRej; // contribution of overweighted
1211  }
1212  }
1213  //********************** MC LOOP ENDS HERE **********************
1214 } // MakeEvent
1215 
1216 ////////////////////////////////////////////////////////////////////////////////
1217 /// User may get generated MC point/vector with help of this method
1218 
1220 {
1221  for ( Int_t k=0 ; k<fDim ; k++) *(MCvect +k) = fMCvect[k];
1222 }//GetMCvect
1223 
1224 ////////////////////////////////////////////////////////////////////////////////
1225 /// User may get weight MC weight using this method
1226 
1228 {
1229  return(fMCwt);
1230 }
1231 ////////////////////////////////////////////////////////////////////////////////
1232 /// User may get weight MC weight using this method
1233 
1235 {
1236  mcwt=fMCwt;
1237 }
1238 
1239 ////////////////////////////////////////////////////////////////////////////////
1240 /// User subprogram which generates MC event and returns MC weight
1241 
1243 {
1244  MakeEvent();
1245  GetMCvect(MCvect);
1246  return(fMCwt);
1247 }//MCgenerate
1248 
1249 ////////////////////////////////////////////////////////////////////////////////
1250 /// User subprogram.
1251 /// It provides the value of the integral calculated from the averages of the MC run
1252 /// May be called after (or during) the MC run.
1253 
1254 void TFoam::GetIntegMC(Double_t &mcResult, Double_t &mcError)
1255 {
1256  Double_t mCerelat;
1257  mcResult = 0.0;
1258  mCerelat = 1.0;
1259  if (fNevGen>0) {
1260  mcResult = fPrime*fSumWt/fNevGen;
1261  mCerelat = sqrt( fSumWt2/(fSumWt*fSumWt) - 1/fNevGen);
1262  }
1263  mcError = mcResult *mCerelat;
1264 }//GetIntegMC
1265 
1266 ////////////////////////////////////////////////////////////////////////////////
1267 /// User subprogram.
1268 /// It returns NORMALIZATION integral to be combined with the average weights
1269 /// and content of the histograms in order to get proper absolute normalization
1270 /// of the integrand and distributions.
1271 /// It can be called after initialization, before or during the MC run.
1272 
1273 void TFoam::GetIntNorm(Double_t& IntNorm, Double_t& Errel )
1274 {
1275  if(fOptRej == 1) { // Wt=1 events, internal rejection
1276  Double_t intMC,errMC;
1277  GetIntegMC(intMC,errMC);
1278  IntNorm = intMC;
1279  Errel = errMC;
1280  } else { // Wted events, NO internal rejection
1281  IntNorm = fPrime;
1282  Errel = 0;
1283  }
1284 }//GetIntNorm
1285 
1286 ////////////////////////////////////////////////////////////////////////////////
1287 /// May be called optionally after the MC run.
1288 /// Returns various parameters of the MC weight for efficiency evaluation
1289 
1290 void TFoam::GetWtParams(Double_t eps, Double_t &aveWt, Double_t &wtMax, Double_t &sigma)
1291 {
1292  Double_t mCeff, wtLim;
1293  fMCMonit->GetMCeff(eps, mCeff, wtLim);
1294  wtMax = wtLim;
1295  aveWt = fSumWt/fNevGen;
1296  sigma = sqrt( fSumWt2/fNevGen -aveWt*aveWt );
1297 }//GetmCeff
1298 
1299 ////////////////////////////////////////////////////////////////////////////////
1300 /// May be called optionally by the user after the MC run.
1301 /// It provides normalization and also prints some information/statistics on the MC run.
1302 
1303 void TFoam::Finalize(Double_t& IntNorm, Double_t& Errel)
1304 {
1305  GetIntNorm(IntNorm,Errel);
1306  Double_t mcResult,mcError;
1307  GetIntegMC(mcResult,mcError);
1308  Double_t mCerelat= mcError/mcResult;
1309  //
1310  if(fChat>0) {
1311  Double_t eps = 0.0005;
1312  Double_t mCeff, mcEf2, wtMax, aveWt, sigma;
1313  GetWtParams(eps, aveWt, wtMax, sigma);
1314  mCeff=0;
1315  if(wtMax>0.0) mCeff=aveWt/wtMax;
1316  mcEf2 = sigma/aveWt;
1317  Double_t driver = fCells[0]->GetDriv();
1318  //
1319  BXOPE;
1320  BXTXT("****************************************");
1321  BXTXT("****** TFoam::Finalize ******");
1322  BXTXT("****************************************");
1323  BX1I(" NevGen",fNevGen, "Number of generated events in the MC generation ");
1324  BX1I(" nCalls",fNCalls, "Total number of function calls ");
1325  BXTXT("----------------------------------------");
1326  BX1F(" AveWt",aveWt, "Average MC weight ");
1327  BX1F(" WtMin",fWtMin, "Minimum MC weight (absolute) ");
1328  BX1F(" WtMax",fWtMax, "Maximum MC weight (absolute) ");
1329  BXTXT("----------------------------------------");
1330  BX1F(" XPrime",fPrime, "Primary total integral, R_prime ");
1331  BX1F(" XDiver",driver, "Driver total integral, R_loss ");
1332  BXTXT("----------------------------------------");
1333  BX2F(" IntMC", mcResult, mcError, "Result of the MC Integral");
1334  BX1F(" mCerelat", mCerelat, "Relative error of the MC integral ");
1335  BX1F(" <w>/WtMax",mCeff, "MC efficiency, acceptance rate");
1336  BX1F(" Sigma/<w>",mcEf2, "MC efficiency, variance/ave_wt");
1337  BX1F(" WtMax",wtMax, "WtMax(esp= 0.0005) ");
1338  BX1F(" Sigma",sigma, "variance of MC weight ");
1339  if(fOptRej==1) {
1340  Double_t avOve=fSumOve/fSumWt;
1341  BX1F("<OveW>/<W>",avOve, "Contrib. of events wt>MaxWtRej");
1342  }
1343  BXCLO;
1344  }
1345 } // Finalize
1346 
1347 ////////////////////////////////////////////////////////////////////////////////
1348 /// This can be called before Initialize, after setting kDim
1349 /// It defines which variables are excluded in the process of the cell division.
1350 /// For example 'FoamX->SetInhiDiv(1, 1);' inhibits division of y-variable.
1351 /// The resulting map of cells in 2-dim. case will look as follows:
1352 ///
1353 // <img src="gif/foam_Map2.gif">
1354 
1355 void TFoam::SetInhiDiv(Int_t iDim, Int_t InhiDiv)
1356 {
1357  if(fDim==0) Error("TFoam","SetInhiDiv: fDim=0 \n");
1358  if(fInhiDiv == 0) {
1359  fInhiDiv = new Int_t[ fDim ];
1360  for(Int_t i=0; i<fDim; i++) fInhiDiv[i]=0;
1361  }
1362  //
1363  if( ( 0<=iDim) && (iDim<fDim)) {
1364  fInhiDiv[iDim] = InhiDiv;
1365  } else
1366  Error("SetInhiDiv:","Wrong iDim \n");
1367 }//SetInhiDiv
1368 
1369 ////////////////////////////////////////////////////////////////////////////////
1370 /// This should be called before Initialize, after setting kDim
1371 /// It predefines values of the cell division for certain variable iDim.
1372 /// For example setting 3 predefined division lines using:
1373 /// xDiv[0]=0.30; xDiv[1]=0.40; xDiv[2]=0.65;
1374 /// FoamX->SetXdivPRD(0,3,xDiv);
1375 /// results in the following 2-dim. pattern of the cells:
1376 ///
1377 // <img src="gif/foam_Map3.gif">
1378 
1379 void TFoam::SetXdivPRD(Int_t iDim, Int_t len, Double_t xDiv[])
1380 {
1381  Int_t i;
1382 
1383  if(fDim<=0) Error("SetXdivPRD", "fDim=0 \n");
1384  if( len<1 ) Error("SetXdivPRD", "len<1 \n");
1385  // allocate list of pointers, if it was not done before
1386  if(fXdivPRD == 0) {
1387  fXdivPRD = new TFoamVect*[fDim];
1388  for(i=0; i<fDim; i++) fXdivPRD[i]=0;
1389  }
1390  // set division list for direction iDim in H-cubic space!!!
1391  if( ( 0<=iDim) && (iDim<fDim)) {
1392  fOptPRD =1; // !!!!
1393  if( fXdivPRD[iDim] != 0)
1394  Error("SetXdivPRD", "Second allocation of XdivPRD not allowed \n");
1395  fXdivPRD[iDim] = new TFoamVect(len); // allocate list of division points
1396  for(i=0; i<len; i++) {
1397  (*fXdivPRD[iDim])[i]=xDiv[i]; // set list of division points
1398  }
1399  } else {
1400  Error("SetXdivPRD", "Wrong iDim \n");
1401  }
1402  // Priting predefined division points
1403  std::cout<<" SetXdivPRD, idim= "<<iDim<<" len= "<<len<<" "<<std::endl;
1404  for(i=0; i<len; i++) {
1405  if (iDim < fDim) std::cout<< (*fXdivPRD[iDim])[i] <<" ";
1406  }
1407  std::cout<<std::endl;
1408  for(i=0; i<len; i++) std::cout<< xDiv[i] <<" ";
1409  std::cout<<std::endl;
1410  //
1411 }//SetXdivPRD
1412 
1413 ////////////////////////////////////////////////////////////////////////////////
1414 /// User utility, miscellaneous and debug.
1415 /// Checks all pointers in the tree of cells. This is useful autodiagnostic.
1416 /// level=0, no printout, failures causes STOP
1417 /// level=1, printout, failures lead to WARNINGS only
1418 
1420 {
1421  Int_t errors, warnings;
1422  TFoamCell *cell;
1423  Long_t iCell;
1424 
1425  errors = 0; warnings = 0;
1426  if (level==1) std::cout << "///////////////////////////// FOAM_Checks /////////////////////////////////" << std::endl;
1427  for(iCell=1; iCell<=fLastCe; iCell++) {
1428  cell = fCells[iCell];
1429  // checking general rules
1430  if( ((cell->GetDau0()==0) && (cell->GetDau1()!=0) ) ||
1431  ((cell->GetDau1()==0) && (cell->GetDau0()!=0) ) ) {
1432  errors++;
1433  if (level==1) Error("CheckAll","ERROR: Cell's no %ld has only one daughter \n",iCell);
1434  }
1435  if( (cell->GetDau0()==0) && (cell->GetDau1()==0) && (cell->GetStat()==0) ) {
1436  errors++;
1437  if (level==1) Error("CheckAll","ERROR: Cell's no %ld has no daughter and is inactive \n",iCell);
1438  }
1439  if( (cell->GetDau0()!=0) && (cell->GetDau1()!=0) && (cell->GetStat()==1) ) {
1440  errors++;
1441  if (level==1) Error("CheckAll","ERROR: Cell's no %ld has two daughters and is active \n",iCell);
1442  }
1443 
1444  // checking parents
1445  if( (cell->GetPare())!=fCells[0] ) { // not child of the root
1446  if ( (cell != cell->GetPare()->GetDau0()) && (cell != cell->GetPare()->GetDau1()) ) {
1447  errors++;
1448  if (level==1) Error("CheckAll","ERROR: Cell's no %ld parent not pointing to this cell\n ",iCell);
1449  }
1450  }
1451 
1452  // checking daughters
1453  if(cell->GetDau0()!=0) {
1454  if(cell != (cell->GetDau0())->GetPare()) {
1455  errors++;
1456  if (level==1) Error("CheckAll","ERROR: Cell's no %ld daughter 0 not pointing to this cell \n",iCell);
1457  }
1458  }
1459  if(cell->GetDau1()!=0) {
1460  if(cell != (cell->GetDau1())->GetPare()) {
1461  errors++;
1462  if (level==1) Error("CheckAll","ERROR: Cell's no %ld daughter 1 not pointing to this cell \n",iCell);
1463  }
1464  }
1465  }// loop after cells;
1466 
1467  // Check for empty cells
1468  for(iCell=0; iCell<=fLastCe; iCell++) {
1469  cell = fCells[iCell];
1470  if( (cell->GetStat()==1) && (cell->GetDriv()==0) ) {
1471  warnings++;
1472  if(level==1) Warning("CheckAll", "Warning: Cell no. %ld is active but empty \n", iCell);
1473  }
1474  }
1475  // summary
1476  if(level==1){
1477  Info("CheckAll","Check has found %d errors and %d warnings \n",errors, warnings);
1478  }
1479  if(errors>0){
1480  Info("CheckAll","Check - found total %d errors \n",errors);
1481  }
1482 } // Check
1483 
1484 ////////////////////////////////////////////////////////////////////////////////
1485 /// Prints geometry of ALL cells of the FOAM
1486 
1488 {
1489  Long_t iCell;
1490 
1491  for(iCell=0; iCell<=fLastCe; iCell++) {
1492  std::cout<<"Cell["<<iCell<<"]={ ";
1493  //std::cout<<" "<< fCells[iCell]<<" "; // extra DEBUG
1494  std::cout<<std::endl;
1495  fCells[iCell]->Print("1");
1496  std::cout<<"}"<<std::endl;
1497  }
1498 }
1499 
1500 ////////////////////////////////////////////////////////////////////////////////
1501 /// Debugging tool which plots 2-dimensional cells as rectangles
1502 /// in C++ format readable for root
1503 
1505 {
1506  std::ofstream outfile(filename, std::ios::out);
1507  Double_t x1,y1,x2,y2,x,y;
1508  Long_t iCell;
1509  Double_t offs =0.1;
1510  Double_t lpag =1-2*offs;
1511  outfile<<"{" << std::endl;
1512  outfile<<"cMap = new TCanvas(\"Map1\",\" Cell Map \",600,600);"<<std::endl;
1513  //
1514  outfile<<"TBox*a=new TBox();"<<std::endl;
1515  outfile<<"a->SetFillStyle(0);"<<std::endl; // big frame
1516  outfile<<"a->SetLineWidth(4);"<<std::endl;
1517  outfile<<"a->SetLineColor(2);"<<std::endl;
1518  outfile<<"a->DrawBox("<<offs<<","<<offs<<","<<(offs+lpag)<<","<<(offs+lpag)<<");"<<std::endl;
1519  //
1520  outfile<<"TText*t=new TText();"<<std::endl; // text for numbering
1521  outfile<<"t->SetTextColor(4);"<<std::endl;
1522  if(fLastCe<51)
1523  outfile<<"t->SetTextSize(0.025);"<<std::endl; // text for numbering
1524  else if(fLastCe<251)
1525  outfile<<"t->SetTextSize(0.015);"<<std::endl;
1526  else
1527  outfile<<"t->SetTextSize(0.008);"<<std::endl;
1528  //
1529  outfile<<"TBox*b=new TBox();"<<std::endl; // single cell
1530  outfile <<"b->SetFillStyle(0);"<<std::endl;
1531  //
1532  if(fDim==2 && fLastCe<=2000) {
1533  TFoamVect cellPosi(fDim); TFoamVect cellSize(fDim);
1534  outfile << "// =========== Rectangular cells ==========="<< std::endl;
1535  for(iCell=1; iCell<=fLastCe; iCell++) {
1536  if( fCells[iCell]->GetStat() == 1) {
1537  fCells[iCell]->GetHcub(cellPosi,cellSize);
1538  x1 = offs+lpag*( cellPosi[0]); y1 = offs+lpag*( cellPosi[1]);
1539  x2 = offs+lpag*(cellPosi[0]+cellSize[0]); y2 = offs+lpag*(cellPosi[1]+cellSize[1]);
1540  // cell rectangle
1541  if(fLastCe<=2000)
1542  outfile<<"b->DrawBox("<<x1<<","<<y1<<","<<x2<<","<<y2<<");"<<std::endl;
1543  // cell number
1544  if(fLastCe<=250) {
1545  x = offs+lpag*(cellPosi[0]+0.5*cellSize[0]); y = offs+lpag*(cellPosi[1]+0.5*cellSize[1]);
1546  outfile<<"t->DrawText("<<x<<","<<y<<","<<"\""<<iCell<<"\""<<");"<<std::endl;
1547  }
1548  }
1549  }
1550  outfile<<"// ============== End Rectangles ==========="<< std::endl;
1551  }//kDim=2
1552  //
1553  //
1554  outfile << "}" << std::endl;
1555  outfile.close();
1556 }
1557 
1559 {
1560 // Void function for backward compatibility
1561 
1562  Info("LinkCells", "VOID function for backward compatibility \n");
1563  return;
1564 }
1565 
1566 ////////////////////////////////////////////////////////////////////////////////
1567 // End of Class TFoam //
1568 ////////////////////////////////////////////////////////////////////////////////
1569 
TRefArray * fCellsAct
Definition: TFoam.h:59
Int_t fNCells
Definition: TFoam.h:37
Int_t fNoAct
Lists of division values encoded in one vector per direction.
Definition: TFoam.h:53
virtual void CheckAll(Int_t)
User utility, miscellaneous and debug.
Definition: TFoam.cxx:1419
virtual void RootPlot2dim(Char_t *)
Debugging tool which plots 2-dimensional cells as rectangles in C++ format readable for root...
Definition: TFoam.cxx:1504
Int_t fLastCe
Definition: TFoam.h:54
Int_t GetStat() const
Definition: TFoamCell.h:70
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
Definition: TH1.cxx:3165
Double_t * fMCvect
Definition: TFoam.h:65
An array of TObjects.
Definition: TObjArray.h:39
Double_t fSumWt
Definition: TFoam.h:75
virtual void Info(const char *method, const char *msgfmt,...) const
Issue info message.
Definition: TObject.cxx:892
#define BX2F(name, numb, err, text)
Definition: TFoam.cxx:139
void Add(TObject *obj)
Definition: TRefArray.h:84
virtual void LinkCells(void)
Definition: TFoam.cxx:1558
virtual Double_t Rndm(Int_t i=0)
Machine independent random number generator.
Definition: TRandom.cxx:512
virtual Long_t PeekMax()
Internal subprogram used by Initialize.
Definition: TFoam.cxx:928
virtual void Delete(Option_t *option="")
Remove all objects from the array AND delete all heap based objects.
Definition: TObjArray.cxx:328
virtual void GetWtParams(Double_t, Double_t &, Double_t &, Double_t &)
May be called optionally after the MC run.
Definition: TFoam.cxx:1290
TFoamCell ** fCells
Definition: TFoam.h:55
Int_t fOptPRD
[fDim] Flags for inhibiting cell division
Definition: TFoam.h:50
void SetDau1(TFoamCell *Daug)
Definition: TFoamCell.h:76
virtual void RndmArray(Int_t n, Float_t *array)
Return an array of n random numbers uniformly distributed in ]0,1].
Definition: TRandom.cxx:548
virtual void Varedu(Double_t [], Int_t &, Double_t &, Double_t &)
Internal subrogram used by Initialize.
Definition: TFoam.cxx:704
virtual Double_t MCgenerate(Double_t *MCvect)
User subprogram which generates MC event and returns MC weight.
Definition: TFoam.cxx:1242
Double_t fMCwt
Definition: TFoam.h:66
static const char * filename()
Double_t * fRvec
Definition: TFoam.h:67
Int_t fNBin
Definition: TFoam.h:44
Double_t fMaxWtRej
Definition: TFoam.h:58
Double_t fNevGen
Definition: TFoam.h:77
TString fVersion
Definition: TFoam.h:34
Basic string class.
Definition: TString.h:137
static Bool_t AddDirectoryStatus()
static function: cannot be inlined on Windows/NT
Definition: TH1.cxx:714
void Print(Option_t *option) const
Printout of the cell geometry parameters for the debug purpose.
Definition: TFoamCell.cxx:189
Short_t Min(Short_t a, Short_t b)
Definition: TMathBase.h:170
int Int_t
Definition: RtypesCore.h:41
virtual void MakeEvent()
User subprogram.
Definition: TFoam.cxx:1169
bool Bool_t
Definition: RtypesCore.h:59
TFoamCell * GetDau1() const
Definition: TFoamCell.h:74
Int_t fChat
Definition: TFoam.h:41
const Bool_t kFALSE
Definition: Rtypes.h:92
virtual ~TFoam()
Default destructor.
Definition: TFoam.cxx:263
Double_t fSumWt2
Definition: TFoam.h:75
TFoamIntegrand * fRho
Definition: TFoam.h:69
An array of references to TObjects.
Definition: TRefArray.h:43
void GetMCeff(Double_t, Double_t &, Double_t &)
Calculates Efficiency= aveWt/wtLim for a given tolerance level epsilon<<1 using information stored in...
Definition: TFoamMaxwt.cxx:124
TFoamVect ** fXdivPRD
Definition: TFoam.h:51
virtual void GenerCel2(TFoamCell *&)
Internal subprogram.
Definition: TFoam.cxx:1131
Double_t fMCerror
Definition: TFoam.h:81
virtual void MakeActiveList()
Internal subrogram used by Initialize.
Definition: TFoam.cxx:993
Double_t fWtMin
Definition: TFoam.h:78
Double_t fMCresult
Definition: TFoam.h:80
Short_t Abs(Short_t d)
Definition: TMathBase.h:110
THist< 1, double > TH1D
Definition: THist.h:314
virtual void Grow()
Internal subrogram used by Initialize.
Definition: TFoam.cxx:893
static void AddDirectory(Bool_t add=kTRUE)
Sets the flag controlling the automatic add of histograms in memory.
Definition: TH1.cxx:1236
const char * Name
Definition: TXMLSetup.cxx:67
virtual void SetXdivPRD(Int_t, Int_t, Double_t[])
This should be called before Initialize, after setting kDim It predefines values of the cell division...
Definition: TFoam.cxx:1379
void SetParamPtrs(void *paramArr, Int_t nparam=-1)
ParamArr is an array containing the function argument values.
virtual Double_t GetMCwt()
User may get weight MC weight using this method.
Definition: TFoam.cxx:1227
double sqrt(double)
static const double x2[5]
virtual void ResetRho(TFoamIntegrand *Rho)
User may optionally reset the distribution using this method Usually it is done when FOAM object is r...
Definition: TFoam.cxx:1078
Double_t x[n]
Definition: legend1.C:17
TMethodCall * fMethodCall
Pointer to the user-defined integrand function/distribution.
Definition: TFoam.h:70
This is the base class for the ROOT Random number generators.
Definition: TRandom.h:29
virtual void PrintCells()
Prints geometry of ALL cells of the FOAM.
Definition: TFoam.cxx:1487
Long_t fNEffev
Definition: TFoam.h:74
Int_t fOptRej
Definition: TFoam.h:42
TObject * At(Int_t idx) const
Definition: TRefArray.h:184
#define BX1I(name, numb, text)
Definition: TFoam.cxx:135
Double_t Sqr(Double_t x) const
Definition: TFoam.h:152
TObjArray * fHistEdg
Definition: TFoam.h:61
Method or function calling interface.
Definition: TMethodCall.h:41
virtual void Reset(Option_t *option="")
Reset.
Definition: TH1.cxx:9704
void CalcVolume()
Calculates volume of the cell using size params which are calculated.
Definition: TFoamCell.cxx:174
TFoamCell * GetDau0() const
Definition: TFoamCell.h:73
Double_t GetPrim() const
Definition: TFoamCell.h:65
virtual void Finalize(Double_t &, Double_t &)
May be called optionally by the user after the MC run.
Definition: TFoam.cxx:1303
void SetXdiv(Double_t Xdiv)
Definition: TFoamCell.h:57
#define BXTXT(text)
Definition: TFoam.cxx:133
Double_t * fPrimAcu
Definition: TFoam.h:60
virtual void MakeAlpha()
Internal subrogram used by Initialize.
Definition: TFoam.cxx:878
void SetDriv(Double_t Driv)
Definition: TFoamCell.h:67
virtual void GetMCvect(Double_t *)
User may get generated MC point/vector with help of this method.
Definition: TFoam.cxx:1219
Double_t * fAlpha
Definition: TFoam.h:83
Int_t fEvPerBin
Definition: TFoam.h:46
virtual void SetInhiDiv(Int_t, Int_t)
This can be called before Initialize, after setting kDim It defines which variables are excluded in t...
Definition: TFoam.cxx:1355
virtual Double_t Density(Int_t ndim, Double_t *)=0
void SetIntg(Double_t Intg)
Definition: TFoamCell.h:66
Double_t GetIntg() const
Definition: TFoamCell.h:63
Double_t GetVolume() const
Definition: TFoamCell.h:62
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition: TObject.cxx:918
virtual void Initialize()
Basic initialization of FOAM invoked by the user.
Definition: TFoam.cxx:356
Int_t * fMaskDiv
Definition: TFoam.h:48
Double_t GetDriv() const
Definition: TFoamCell.h:64
virtual void SetRhoInt(void *Rho)
User may use this to set pointer to the global function (not descending from TFoamIntegrand) serving ...
Definition: TFoam.cxx:1094
TFoamCell * GetPare() const
Definition: TFoamCell.h:72
TString fName
Definition: TFoam.h:33
1-D histogram with a double per channel (see TH1 documentation)}
Definition: TH1.h:613
void SetStat(Int_t Stat)
Definition: TFoamCell.h:71
Double_t fWtMax
Definition: TFoam.h:78
Int_t fDim
Definition: TFoam.h:36
static const Double_t gVlow
Definition: TFoam.cxx:148
#define BXCLO
Definition: TFoam.cxx:142
void InitWithPrototype(TClass *cl, const char *method, const char *proto, Bool_t objectIsConst=kFALSE, ROOT::EFunctionMatchMode mode=ROOT::kConversionMatch)
Initialize the method invocation environment.
void Fill(Int_t, TFoamCell *, TFoamCell *, TFoamCell *)
Fills in certain data into newly allocated cell.
Definition: TFoamCell.cxx:103
long Long_t
Definition: RtypesCore.h:50
Int_t fRNmax
Definition: TFoam.h:38
Int_t GetBest() const
Definition: TFoamCell.h:55
virtual void Carver(Int_t &, Double_t &, Double_t &)
Internal subrogram used by Initialize.
Definition: TFoam.cxx:780
static const double x1[5]
#define ClassImp(name)
Definition: Rtypes.h:279
double Double_t
Definition: RtypesCore.h:55
virtual const char * Getp2f2funcname(void *) const
Definition: TInterpreter.h:213
void GetHcub(TFoamVect &, TFoamVect &) const
Provides size and position of the cell These parameter are calculated by analyzing information in all...
Definition: TFoamCell.cxx:120
Double_t fSumOve
Definition: TFoam.h:76
Double_t y[n]
Definition: legend1.C:17
double func(double *x, double *p)
Definition: stressTF1.cxx:213
#define BXOPE
Definition: TFoam.cxx:130
virtual void SetPseRan(TRandom *PseRan)
Definition: TFoam.h:124
virtual void InitCells()
Internal subprogram used by Initialize.
Definition: TFoam.cxx:480
Int_t fOptDrive
Definition: TFoam.h:40
void SetSerial(Int_t Serial)
Definition: TFoamCell.h:77
Mother of all ROOT objects.
Definition: TObject.h:58
char Char_t
Definition: RtypesCore.h:29
virtual Int_t CellFill(Int_t, TFoamCell *)
Internal subprogram used by Initialize.
Definition: TFoam.cxx:512
virtual void GetIntNorm(Double_t &, Double_t &)
User subprogram.
Definition: TFoam.cxx:1273
virtual void Explore(TFoamCell *Cell)
Internal subprogram used by Initialize.
Definition: TFoam.cxx:553
TString fDate
Definition: TFoam.h:35
TObjArray * fHistDbg
Definition: TFoam.h:62
void Execute(const char *, const char *, int *=0)
Execute method on this object with the given parameter string, e.g.
Definition: TMethodCall.h:68
Short_t Max(Short_t a, Short_t b)
Definition: TMathBase.h:202
Double_t fPrime
Definition: TFoam.h:79
void SetBest(Int_t Best)
Definition: TFoamCell.h:56
void SetPrim(Double_t Prim)
Definition: TFoamCell.h:68
virtual void ResetPseRan(TRandom *PseRan)
User may optionally reset random number generator using this method Usually it is done when FOAM obje...
Definition: TFoam.cxx:1036
void Fill(Double_t)
Filling analyzed weight.
Definition: TFoamMaxwt.cxx:97
#define BX1F(name, numb, text)
Definition: TFoam.cxx:137
TFoamMaxwt * fMCMonit
Definition: TFoam.h:57
float type_of_call hi(const int &, const int &)
double result[121]
virtual Double_t Eval(Double_t *)
Internal subprogram.
Definition: TFoam.cxx:1109
virtual Int_t Divide(TFoamCell *)
Internal subrogram used by Initialize.
Definition: TFoam.cxx:961
static const Double_t gHigh
Definition: TFoam.cxx:147
TH1D * fHistWt
Definition: TFoam.h:63
R__EXTERN TInterpreter * gCling
Definition: TInterpreter.h:504
Int_t * fInhiDiv
[fDim] Dynamic Mask for cell division
Definition: TFoam.h:49
void SetDau0(TFoamCell *Daug)
Definition: TFoamCell.h:75
const Int_t n
Definition: legend1.C:16
virtual void SetRho(TFoamIntegrand *Rho)
User may use this method to set the distribution object.
Definition: TFoam.cxx:1048
Long_t fNCalls
Definition: TFoam.h:73
virtual void GetIntegMC(Double_t &, Double_t &)
User subprogram.
Definition: TFoam.cxx:1254
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition: TObject.cxx:904
Definition: TFoam.h:29
Int_t fNSampl
Definition: TFoam.h:45
const char * Data() const
Definition: TString.h:349
TFoam()
Default constructor for streamer, user should not use it.
Definition: TFoam.cxx:178
TRandom * fPseRan
ROOT&#39;s pointer to user-defined global distribution function.
Definition: TFoam.h:71