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