00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046
00047
00048
00049
00050
00051
00052
00053
00054
00055
00056
00057
00058
00059
00060
00061
00062
00063
00064
00065
00066
00067
00068
00069
00070
00071
00072
00073
00074
00075
00076
00077
00078
00079
00080
00081
00082
00083
00084
00085
00086
00087
00088
00089
00090
00091
00092
00093
00094
00095
00096
00097
00098
00099
00100
00101
00102
00103
00104
00105
00106
00107
00108
00109
00110
00111
00112
00113
00114
00115
00116
00117
00118
00119
00120
00121 #include "TFoam.h"
00122 #include "TFoamIntegrand.h"
00123 #include "TFoamMaxwt.h"
00124 #include "TFoamVect.h"
00125 #include "TFoamCell.h"
00126 #include "Riostream.h"
00127 #include "TH1.h"
00128 #include "TRefArray.h"
00129 #include "TMethodCall.h"
00130 #include "TRandom.h"
00131 #include "TMath.h"
00132 #include "TInterpreter.h"
00133
00134 ClassImp(TFoam);
00135
00136
00137 #define BXOPE cout<<\
00138 "FFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFF"<<endl<<\
00139 "F F"<<endl
00140 #define BXTXT(text) cout<<\
00141 "F "<<setw(40)<< text <<" F"<<endl
00142 #define BX1I(name,numb,text) cout<<\
00143 "F "<<setw(10)<<name<<" = "<<setw(10)<<numb<<" = " <<setw(50)<<text<<" F"<<endl
00144 #define BX1F(name,numb,text) cout<<"F "<<setw(10)<<name<<\
00145 " = "<<setw(15)<<setprecision(8)<<numb<<" = "<<setw(40)<<text<<" F"<<endl
00146 #define BX2F(name,numb,err,text) cout<<"F "<<setw(10)<<name<<\
00147 " = "<<setw(15)<<setprecision(8)<<numb<<" +- "<<setw(15)<<setprecision(8)<<err<<\
00148 " = "<<setw(25)<<text<<" F"<<endl
00149 #define BXCLO cout<<\
00150 "F F"<<endl<<\
00151 "FFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFF"<<endl
00152
00153
00154 static const Double_t gHigh= 1.0e150;
00155 static const Double_t gVlow=-1.0e150;
00156
00157 #define SW2 setprecision(7) << setw(12)
00158
00159
00160 TFoam::TFoam() :
00161 fDim(0), fNCells(0), fRNmax(0),
00162 fOptDrive(0), fChat(0), fOptRej(0),
00163 fNBin(0), fNSampl(0), fEvPerBin(0),
00164 fMaskDiv(0), fInhiDiv(0), fOptPRD(0), fXdivPRD(0),
00165 fNoAct(0), fLastCe(0), fCells(0),
00166 fMCMonit(0), fMaxWtRej(0), fCellsAct(0), fPrimAcu(0),
00167 fHistEdg(0), fHistDbg(0), fHistWt(0),
00168 fMCvect(0), fMCwt(0), fRvec(0),
00169 fRho(0), fMethodCall(0), fPseRan(0),
00170 fNCalls(0), fNEffev(0),
00171 fSumWt(0), fSumWt2(0),
00172 fSumOve(0), fNevGen(0),
00173 fWtMax(0), fWtMin(0),
00174 fPrime(0), fMCresult(0), fMCerror(0),
00175 fAlpha(0)
00176 {
00177
00178 }
00179
00180 TFoam::TFoam(const Char_t* Name) :
00181 fDim(0), fNCells(0), fRNmax(0),
00182 fOptDrive(0), fChat(0), fOptRej(0),
00183 fNBin(0), fNSampl(0), fEvPerBin(0),
00184 fMaskDiv(0), fInhiDiv(0), fOptPRD(0), fXdivPRD(0),
00185 fNoAct(0), fLastCe(0), fCells(0),
00186 fMCMonit(0), fMaxWtRej(0), fCellsAct(0), fPrimAcu(0),
00187 fHistEdg(0), fHistDbg(0), fHistWt(0),
00188 fMCvect(0), fMCwt(0), fRvec(0),
00189 fRho(0), fMethodCall(0), fPseRan(0),
00190 fNCalls(0), fNEffev(0),
00191 fSumWt(0), fSumWt2(0),
00192 fSumOve(0), fNevGen(0),
00193 fWtMax(0), fWtMin(0),
00194 fPrime(0), fMCresult(0), fMCerror(0),
00195 fAlpha(0)
00196 {
00197
00198
00199 if(strlen(Name) >129) {
00200 Error("TFoam","Name too long %s \n",Name);
00201 }
00202 fName=Name;
00203 fDate=" Release date: 2005.04.10";
00204 fVersion= "1.02M";
00205 fMaskDiv = 0;
00206 fInhiDiv = 0;
00207 fXdivPRD = 0;
00208 fCells = 0;
00209 fAlpha = 0;
00210 fCellsAct = 0;
00211 fPrimAcu = 0;
00212 fHistEdg = 0;
00213 fHistWt = 0;
00214 fHistDbg = 0;
00215 fDim = 0;
00216 fNCells = 1000;
00217 fNSampl = 200;
00218 fOptPRD = 0;
00219 fOptDrive = 2;
00220 fChat = 1;
00221 fOptRej = 1;
00222
00223 fNBin = 8;
00224 fEvPerBin =25;
00225
00226 fNCalls = 0;
00227 fNEffev = 0;
00228 fLastCe =-1;
00229 fNoAct = 0;
00230 fWtMin = gHigh;
00231 fWtMax = gVlow;
00232 fMaxWtRej =1.10;
00233 fPseRan = 0;
00234 fMCMonit = 0;
00235 fRho = 0;
00236 fMCvect = 0;
00237 fRvec = 0;
00238 fPseRan = 0;
00239 fMethodCall=0;
00240 }
00241
00242
00243 TFoam::~TFoam()
00244 {
00245
00246
00247 Int_t i;
00248
00249 if(fCells!= 0) {
00250 for(i=0; i<fNCells; i++) delete fCells[i];
00251 delete [] fCells;
00252 }
00253 if (fCellsAct) delete fCellsAct ;
00254 if (fRvec) delete [] fRvec;
00255 if (fAlpha) delete [] fAlpha;
00256 if (fMCvect) delete [] fMCvect;
00257 if (fPrimAcu) delete [] fPrimAcu;
00258 if (fMaskDiv) delete [] fMaskDiv;
00259 if (fInhiDiv) delete [] fInhiDiv;
00260
00261 if( fXdivPRD!= 0) {
00262 for(i=0; i<fDim; i++) delete fXdivPRD[i];
00263 delete [] fXdivPRD;
00264 }
00265 if (fMCMonit) delete fMCMonit;
00266 if (fHistWt) delete fHistWt;
00267
00268
00269 if (fHistEdg) {
00270 fHistEdg->Delete();
00271 delete fHistEdg;
00272 }
00273 if (fHistDbg) {
00274 fHistDbg->Delete();
00275 delete fHistDbg;
00276 }
00277 }
00278
00279
00280
00281 TFoam::TFoam(const TFoam &From): TObject(From)
00282 {
00283
00284 Error("TFoam", "COPY CONSTRUCTOR NOT IMPLEMENTED \n");
00285 }
00286
00287
00288 void TFoam::Initialize(TRandom *PseRan, TFoamIntegrand *fun )
00289 {
00290
00291
00292
00293
00294
00295
00296
00297
00298
00299
00300
00301
00302
00303
00304
00305
00306
00307
00308
00309
00310
00311
00312
00313
00314
00315
00316
00317
00318
00319
00320
00321
00322
00323
00324
00325
00326
00327
00328
00329
00330 SetPseRan(PseRan);
00331 SetRho(fun);
00332 Initialize();
00333 }
00334
00335
00336 void TFoam::Initialize()
00337 {
00338
00339
00340
00341
00342 Bool_t addStatus = TH1::AddDirectoryStatus();
00343 TH1::AddDirectory(kFALSE);
00344 Int_t i;
00345
00346 if(fChat>0){
00347 BXOPE;
00348 BXTXT("****************************************");
00349 BXTXT("****** TFoam::Initialize ******");
00350 BXTXT("****************************************");
00351 BXTXT(fName);
00352 BX1F(" Version",fVersion, fDate);
00353 BX1I(" kDim",fDim, " Dimension of the hyper-cubical space ");
00354 BX1I(" nCells",fNCells, " Requested number of Cells (half of them active) ");
00355 BX1I(" nSampl",fNSampl, " No of MC events in exploration of a cell ");
00356 BX1I(" nBin",fNBin, " No of bins in histograms, MC exploration of cell ");
00357 BX1I(" EvPerBin",fEvPerBin, " Maximum No effective_events/bin, MC exploration ");
00358 BX1I(" OptDrive",fOptDrive, " Type of Driver =1,2 for Sigma,WtMax ");
00359 BX1I(" OptRej",fOptRej, " MC rejection on/off for OptRej=0,1 ");
00360 BX1F(" MaxWtRej",fMaxWtRej, " Maximum wt in rejection for wt=1 evts");
00361 BXCLO;
00362 }
00363
00364 if(fPseRan==0) Error("Initialize", "Random number generator not set \n");
00365 if(fRho==0 && fMethodCall==0 ) Error("Initialize", "Distribution function not set \n");
00366 if(fDim==0) Error("Initialize", "Zero dimension not allowed \n");
00367
00368
00369
00370
00371
00372 fRNmax= fDim+1;
00373 fRvec = new Double_t[fRNmax];
00374 if(fRvec==0) Error("Initialize", "Cannot initialize buffer fRvec \n");
00375
00376 if(fDim>0){
00377 fAlpha = new Double_t[fDim];
00378 if(fAlpha==0) Error("Initialize", "Cannot initialize buffer fAlpha \n" );
00379 }
00380 fMCvect = new Double_t[fDim];
00381 if(fMCvect==0) Error("Initialize", "Cannot initialize buffer fMCvect \n" );
00382
00383
00384 if(fInhiDiv == 0){
00385 fInhiDiv = new Int_t[fDim];
00386 for(i=0; i<fDim; i++) fInhiDiv[i]=0;
00387 }
00388
00389 if(fMaskDiv == 0){
00390 fMaskDiv = new Int_t[fDim];
00391 for(i=0; i<fDim; i++) fMaskDiv[i]=1;
00392 }
00393
00394 if(fXdivPRD == 0){
00395 fXdivPRD = new TFoamVect*[fDim];
00396 for(i=0; i<fDim; i++) fXdivPRD[i]=0;
00397 }
00398
00399 fHistWt = new TH1D("HistWt","Histogram of MC weight",100,0.0, 1.5*fMaxWtRej);
00400 fHistEdg = new TObjArray(fDim);
00401 TString hname;
00402 TString htitle;
00403 for(i=0;i<fDim;i++){
00404 hname=fName+TString("_HistEdge_");
00405 hname+=i;
00406 htitle=TString("Edge Histogram No. ");
00407 htitle+=i;
00408
00409 (*fHistEdg)[i] = new TH1D(hname.Data(),htitle.Data(),fNBin,0.0, 1.0);
00410 ((TH1D*)(*fHistEdg)[i])->Sumw2();
00411 }
00412
00413 fHistDbg = new TObjArray(fDim);
00414 for(i=0;i<fDim;i++){
00415 hname=fName+TString("_HistDebug_");
00416 hname+=i;
00417 htitle=TString("Debug Histogram ");
00418 htitle+=i;
00419 (*fHistDbg)[i] = new TH1D(hname.Data(),htitle.Data(),fNBin,0.0, 1.0);
00420 }
00421
00422
00423
00424
00425
00426
00427 InitCells();
00428
00429 Grow();
00430
00431
00432 MakeActiveList();
00433
00434
00435 fSumWt = 0.0;
00436 fSumWt2 = 0.0;
00437 fSumOve = 0.0;
00438 fNevGen = 0.0;
00439 fWtMax = gVlow;
00440 fWtMin = gHigh;
00441 fMCresult=fCells[0]->GetIntg();
00442 fMCresult=fCells[0]->GetIntg();
00443 fMCerror =fCells[0]->GetIntg();
00444 fMCMonit = new TFoamMaxwt(5.0,1000);
00445
00446 if(fChat>0){
00447 Double_t driver = fCells[0]->GetDriv();
00448 BXOPE;
00449 BXTXT("*** TFoam::Initialize FINISHED!!! ***");
00450 BX1I(" nCalls",fNCalls, "Total number of function calls ");
00451 BX1F(" XPrime",fPrime, "Primary total integral ");
00452 BX1F(" XDiver",driver, "Driver total integral ");
00453 BX1F(" mcResult",fMCresult,"Estimate of the true MC Integral ");
00454 BXCLO;
00455 }
00456 if(fChat==2) PrintCells();
00457 TH1::AddDirectory(addStatus);
00458 }
00459
00460
00461 void TFoam::InitCells()
00462 {
00463
00464
00465
00466 Int_t i;
00467
00468 fLastCe =-1;
00469 if(fCells!= 0) {
00470 for(i=0; i<fNCells; i++) delete fCells[i];
00471 delete [] fCells;
00472 }
00473
00474 fCells = new TFoamCell*[fNCells];
00475 for(i=0;i<fNCells;i++){
00476 fCells[i]= new TFoamCell(fDim);
00477 fCells[i]->SetSerial(i);
00478 }
00479 if(fCells==0) Error("InitCells", "Cannot initialize CELLS \n" );
00480
00481
00482
00483
00484 CellFill(1, 0);
00485
00486
00487 for(Long_t iCell=0; iCell<=fLastCe; iCell++){
00488 Explore( fCells[iCell] );
00489 }
00490 }
00491
00492
00493 Int_t TFoam::CellFill(Int_t Status, TFoamCell *parent)
00494 {
00495
00496
00497
00498 TFoamCell *cell;
00499 if (fLastCe==fNCells){
00500 Error( "CellFill", "Too many cells\n");
00501 }
00502 fLastCe++;
00503 if (Status==1) fNoAct++;
00504
00505 cell = fCells[fLastCe];
00506
00507 cell->Fill(Status, parent, 0, 0);
00508
00509 cell->SetBest( -1);
00510 cell->SetXdiv(0.5);
00511 Double_t xInt2,xDri2;
00512 if(parent!=0){
00513 xInt2 = 0.5*parent->GetIntg();
00514 xDri2 = 0.5*parent->GetDriv();
00515 cell->SetIntg(xInt2);
00516 cell->SetDriv(xDri2);
00517 }else{
00518 cell->SetIntg(0.0);
00519 cell->SetDriv(0.0);
00520 }
00521 return fLastCe;
00522 }
00523
00524
00525 void TFoam::Explore(TFoamCell *cell)
00526 {
00527
00528
00529
00530
00531
00532
00533
00534
00535
00536
00537
00538
00539 Double_t wt, dx, xBest=0, yBest=0;
00540 Double_t intOld, driOld;
00541
00542 Long_t iev;
00543 Double_t nevMC;
00544 Int_t i, j, k;
00545 Int_t nProj, kBest;
00546 Double_t ceSum[5], xproj;
00547
00548 TFoamVect cellSize(fDim);
00549 TFoamVect cellPosi(fDim);
00550
00551 cell->GetHcub(cellPosi,cellSize);
00552
00553 TFoamCell *parent;
00554
00555 Double_t *xRand = new Double_t[fDim];
00556
00557 Double_t *volPart=0;
00558
00559 cell->CalcVolume();
00560 dx = cell->GetVolume();
00561 intOld = cell->GetIntg();
00562 driOld = cell->GetDriv();
00563
00564
00565
00566
00567
00568 ceSum[0]=0;
00569 ceSum[1]=0;
00570 ceSum[2]=0;
00571 ceSum[3]=gHigh;
00572 ceSum[4]=gVlow;
00573
00574 for(i=0;i<fDim;i++) ((TH1D *)(*fHistEdg)[i])->Reset();
00575 fHistWt->Reset();
00576
00577
00578 Double_t nevEff=0.;
00579 for(iev=0;iev<fNSampl;iev++){
00580 MakeAlpha();
00581
00582 if(fDim>0){
00583 for(j=0; j<fDim; j++)
00584 xRand[j]= cellPosi[j] +fAlpha[j]*(cellSize[j]);
00585 }
00586
00587 wt=dx*Eval(xRand);
00588
00589 nProj = 0;
00590 if(fDim>0) {
00591 for(k=0; k<fDim; k++) {
00592 xproj =fAlpha[k];
00593 ((TH1D *)(*fHistEdg)[nProj])->Fill(xproj,wt);
00594 nProj++;
00595 }
00596 }
00597
00598 fNCalls++;
00599 ceSum[0] += wt;
00600 ceSum[1] += wt*wt;
00601 ceSum[2]++;
00602 if (ceSum[3]>wt) ceSum[3]=wt;
00603 if (ceSum[4]<wt) ceSum[4]=wt;
00604
00605 nevEff = ceSum[0]*ceSum[0]/ceSum[1];
00606 if( nevEff >= fNBin*fEvPerBin) break;
00607 }
00608
00609
00610 for(k=0; k<fDim;k++){
00611 fMaskDiv[k] =1;
00612 if( fInhiDiv[k]==1) fMaskDiv[k] =0;
00613 }
00614
00615 kBest=-1;
00616 Double_t rmin,rmax,rdiv;
00617 if(fOptPRD) {
00618 for(k=0; k<fDim; k++) {
00619 rmin= cellPosi[k];
00620 rmax= cellPosi[k] +cellSize[k];
00621 if( fXdivPRD[k] != 0) {
00622 Int_t n= (fXdivPRD[k])->GetDim();
00623 for(j=0; j<n; j++) {
00624 rdiv=(*fXdivPRD[k])[j];
00625
00626 if( (rmin +1e-99 <rdiv) && (rdiv< rmax -1e-99)) {
00627 kBest=k;
00628 xBest= (rdiv-cellPosi[k])/cellSize[k] ;
00629 goto ee05;
00630 }
00631 }
00632 }
00633 }
00634 }
00635 ee05:
00636
00637 fNEffev += (Long_t)nevEff;
00638 nevMC = ceSum[2];
00639 Double_t intTrue = ceSum[0]/(nevMC+0.000001);
00640 Double_t intDriv=0.;
00641 Double_t intPrim=0.;
00642
00643 switch(fOptDrive){
00644 case 1:
00645 if(kBest == -1) Varedu(ceSum,kBest,xBest,yBest);
00646
00647 intDriv =sqrt(ceSum[1]/nevMC) -intTrue;
00648 intPrim =sqrt(ceSum[1]/nevMC);
00649 break;
00650 case 2:
00651 if(kBest == -1) Carver(kBest,xBest,yBest);
00652 intDriv =ceSum[4] -intTrue;
00653 intPrim =ceSum[4];
00654 break;
00655 default:
00656 Error("Explore", "Wrong fOptDrive = \n" );
00657 }
00658
00659
00660
00661
00662
00663 cell->SetBest(kBest);
00664 cell->SetXdiv(xBest);
00665 cell->SetIntg(intTrue);
00666 cell->SetDriv(intDriv);
00667 cell->SetPrim(intPrim);
00668
00669 Double_t parIntg, parDriv;
00670 for(parent = cell->GetPare(); parent!=0; parent = parent->GetPare()){
00671 parIntg = parent->GetIntg();
00672 parDriv = parent->GetDriv();
00673 parent->SetIntg( parIntg +intTrue -intOld );
00674 parent->SetDriv( parDriv +intDriv -driOld );
00675 }
00676 delete [] volPart;
00677 delete [] xRand;
00678
00679 }
00680
00681
00682 void TFoam::Varedu(Double_t ceSum[5], Int_t &kBest, Double_t &xBest, Double_t &yBest)
00683 {
00684
00685
00686
00687
00688
00689 Double_t nent = ceSum[2];
00690 Double_t swAll = ceSum[0];
00691 Double_t sswAll = ceSum[1];
00692 Double_t ssw = sqrt(sswAll)/sqrt(nent);
00693
00694 Double_t swIn,swOut,sswIn,sswOut,xLo,xUp;
00695 kBest =-1;
00696 xBest =0.5;
00697 yBest =1.0;
00698 Double_t maxGain=0.0;
00699
00700 for(Int_t kProj=0; kProj<fDim; kProj++) {
00701 if( fMaskDiv[kProj]) {
00702
00703 Double_t sigmIn =0.0; Double_t sigmOut =0.0;
00704 Double_t sswtBest = gHigh;
00705 Double_t gain =0.0;
00706 Double_t xMin=0.0; Double_t xMax=0.0;
00707
00708 for(Int_t jLo=1; jLo<=fNBin; jLo++) {
00709 Double_t aswIn=0; Double_t asswIn=0;
00710 for(Int_t jUp=jLo; jUp<=fNBin;jUp++) {
00711 aswIn += ((TH1D *)(*fHistEdg)[kProj])->GetBinContent(jUp);
00712 asswIn += Sqr(((TH1D *)(*fHistEdg)[kProj])->GetBinError( jUp));
00713 xLo=(jLo-1.0)/fNBin;
00714 xUp=(jUp*1.0)/fNBin;
00715 swIn = aswIn/nent;
00716 swOut = (swAll-aswIn)/nent;
00717 sswIn = sqrt(asswIn) /sqrt(nent*(xUp-xLo)) *(xUp-xLo);
00718 sswOut= sqrt(sswAll-asswIn)/sqrt(nent*(1.0-xUp+xLo)) *(1.0-xUp+xLo);
00719 if( (sswIn+sswOut) < sswtBest) {
00720 sswtBest = sswIn+sswOut;
00721 gain = ssw-sswtBest;
00722 sigmIn = sswIn -swIn;
00723 sigmOut = sswOut-swOut;
00724 xMin = xLo;
00725 xMax = xUp;
00726 }
00727 }
00728 }
00729 Int_t iLo = (Int_t) (fNBin*xMin);
00730 Int_t iUp = (Int_t) (fNBin*xMax);
00731
00732
00733
00734
00735 for(Int_t iBin=1;iBin<=fNBin;iBin++) {
00736 if( ((iBin-0.5)/fNBin > xMin) && ((iBin-0.5)/fNBin < xMax) ){
00737 ((TH1D *)(*fHistDbg)[kProj])->SetBinContent(iBin,sigmIn/(xMax-xMin));
00738 } else {
00739 ((TH1D *)(*fHistDbg)[kProj])->SetBinContent(iBin,sigmOut/(1-xMax+xMin));
00740 }
00741 }
00742 if(gain>=maxGain) {
00743 maxGain=gain;
00744 kBest=kProj;
00745 xBest=xMin;
00746 yBest=xMax;
00747 if(iLo == 0) xBest=yBest;
00748 if(iUp == fNBin) yBest=xBest;
00749 }
00750 }
00751 }
00752
00753
00754 if( (kBest >= fDim) || (kBest<0) ) Error("Varedu", "Something wrong with kBest \n" );
00755 }
00756
00757
00758 void TFoam::Carver(Int_t &kBest, Double_t &xBest, Double_t &yBest)
00759 {
00760
00761
00762
00763
00764
00765 Int_t kProj,iBin;
00766 Double_t carve,carvTot,carvMax,carvOne,binMax,binTot,primTot,primMax;
00767 Int_t jLow,jUp,iLow,iUp;
00768 Double_t theBin;
00769 Int_t jDivi;
00770 Int_t j;
00771
00772 Double_t *bins = new Double_t[fNBin];
00773 if(bins==0) Error("Carver", "Cannot initialize buffer Bins \n" );
00774
00775 kBest =-1;
00776 xBest =0.5;
00777 yBest =1.0;
00778 carvMax = gVlow;
00779 primMax = gVlow;
00780 for(kProj=0; kProj<fDim; kProj++)
00781 if( fMaskDiv[kProj] ){
00782
00783
00784
00785 binMax = gVlow;
00786 for(iBin=0; iBin<fNBin;iBin++){
00787 bins[iBin]= ((TH1D *)(*fHistEdg)[kProj])->GetBinContent(iBin+1);
00788 binMax = TMath::Max( binMax, bins[iBin]);
00789 }
00790 if(binMax < 0 ) {
00791 delete [] bins;
00792 return;
00793 }
00794 carvTot = 0.0;
00795 binTot = 0.0;
00796 for(iBin=0;iBin<fNBin;iBin++){
00797 carvTot = carvTot + (binMax-bins[iBin]);
00798 binTot +=bins[iBin];
00799 }
00800 primTot = binMax*fNBin;
00801
00802 jLow =0;
00803 jUp =fNBin-1;
00804 carvOne = gVlow;
00805 Double_t yLevel = gVlow;
00806 for(iBin=0; iBin<fNBin;iBin++) {
00807 theBin = bins[iBin];
00808
00809 iLow = iBin;
00810 for(j=iBin; j>-1; j-- ) {
00811 if(theBin< bins[j]) break;
00812 iLow = j;
00813 }
00814
00815
00816
00817 iUp = iBin;
00818 for(j=iBin; j<fNBin; j++) {
00819 if(theBin< bins[j]) break;
00820 iUp = j;
00821 }
00822
00823
00824
00825 carve = (iUp-iLow+1)*(binMax-theBin);
00826 if( carve > carvOne) {
00827 carvOne = carve;
00828 jLow = iLow;
00829 jUp = iUp;
00830 yLevel = theBin;
00831 }
00832 }
00833 if( carvTot > carvMax) {
00834 carvMax = carvTot;
00835 primMax = primTot;
00836
00837 kBest = kProj;
00838 xBest = ((Double_t)(jLow))/fNBin;
00839 yBest = ((Double_t)(jUp+1))/fNBin;
00840 if(jLow == 0 ) xBest = yBest;
00841 if(jUp == fNBin-1) yBest = xBest;
00842
00843 jDivi = jLow;
00844 if(jLow == 0 ) jDivi=jUp+1;
00845 }
00846
00847
00848 for(iBin=0; iBin<fNBin; iBin++)
00849 ((TH1D *)(*fHistDbg)[kProj])->SetBinContent(iBin+1,binMax);
00850 for(iBin=jLow; iBin<jUp+1; iBin++)
00851 ((TH1D *)(*fHistDbg)[kProj])->SetBinContent(iBin+1,yLevel);
00852 }
00853 if( (kBest >= fDim) || (kBest<0) ) Error("Carver", "Something wrong with kBest \n" );
00854 delete [] bins;
00855 }
00856
00857
00858 void TFoam::MakeAlpha()
00859 {
00860
00861
00862 Int_t k;
00863 if(fDim<1) return;
00864
00865
00866 fPseRan->RndmArray(fDim,fRvec);
00867 for(k=0; k<fDim; k++) fAlpha[k] = fRvec[k];
00868 }
00869
00870
00871
00872 void TFoam::Grow()
00873 {
00874
00875
00876
00877 Long_t iCell;
00878 TFoamCell* newCell;
00879
00880 while ( (fLastCe+2) < fNCells ) {
00881 iCell = PeekMax();
00882 if( (iCell<0) || (iCell>fLastCe) ) Error("Grow", "Wrong iCell \n");
00883 newCell = fCells[iCell];
00884
00885 if(fLastCe !=0) {
00886 Int_t kEcho=10;
00887 if(fLastCe>=10000) kEcho=100;
00888 if( (fLastCe%kEcho)==0) {
00889 if (fChat>0) {
00890 if(fDim<10)
00891 cout<<fDim<<flush;
00892 else
00893 cout<<"."<<flush;
00894 if( (fLastCe%(100*kEcho))==0) cout<<"|"<<fLastCe<<endl<<flush;
00895 }
00896 }
00897 }
00898 if( Divide( newCell )==0) break;
00899 }
00900 if (fChat>0) {
00901 cout<<endl<<flush;
00902 }
00903 CheckAll(0);
00904 }
00905
00906
00907 Long_t TFoam::PeekMax()
00908 {
00909
00910
00911
00912 Long_t i;
00913 Long_t iCell = -1;
00914 Double_t drivMax, driv;
00915
00916 drivMax = gVlow;
00917 for(i=0; i<=fLastCe; i++) {
00918 if( fCells[i]->GetStat() == 1 ) {
00919 driv = TMath::Abs( fCells[i]->GetDriv());
00920
00921 if(driv > drivMax) {
00922 drivMax = driv;
00923 iCell = i;
00924 }
00925 }
00926 }
00927
00928 if (iCell == -1)
00929 cout << "STOP in TFoam::PeekMax: not found iCell=" << iCell << endl;
00930 return(iCell);
00931 }
00932
00933
00934 Int_t TFoam::Divide(TFoamCell *cell)
00935 {
00936
00937
00938
00939
00940
00941
00942
00943
00944
00945 Double_t xdiv;
00946 Int_t kBest;
00947
00948 if(fLastCe+1 >= fNCells) Error("Divide", "Buffer limit is reached, fLastCe=fnBuf \n");
00949
00950 cell->SetStat(0);
00951 fNoAct--;
00952
00953 xdiv = cell->GetXdiv();
00954 kBest = cell->GetBest();
00955 if( kBest<0 || kBest>=fDim ) Error("Divide", "Wrong kBest \n");
00956
00957
00958
00959
00960
00961 Int_t d1 = CellFill(1, cell);
00962 Int_t d2 = CellFill(1, cell);
00963 cell->SetDau0((fCells[d1]));
00964 cell->SetDau1((fCells[d2]));
00965 Explore( (fCells[d1]) );
00966 Explore( (fCells[d2]) );
00967 return 1;
00968 }
00969
00970
00971
00972 void TFoam::MakeActiveList()
00973 {
00974
00975
00976
00977
00978
00979 Long_t n, iCell;
00980 Double_t sum;
00981
00982 if(fPrimAcu != 0) delete [] fPrimAcu;
00983 if(fCellsAct != 0) delete fCellsAct;
00984
00985
00986 fCellsAct = new TRefArray();
00987
00988
00989
00990
00991 fPrime = 0.0; n = 0;
00992 for(iCell=0; iCell<=fLastCe; iCell++) {
00993 if (fCells[iCell]->GetStat()==1) {
00994 fPrime += fCells[iCell]->GetPrim();
00995 fCellsAct->Add(fCells[iCell]);
00996 n++;
00997 }
00998 }
00999
01000 if(fNoAct != n) Error("MakeActiveList", "Wrong fNoAct \n" );
01001 if(fPrime == 0.) Error("MakeActiveList", "Integrand function is zero \n" );
01002
01003 fPrimAcu = new Double_t[fNoAct];
01004 if( fCellsAct==0 || fPrimAcu==0 ) Error("MakeActiveList", "Cant allocate fCellsAct or fPrimAcu \n");
01005
01006 sum =0.0;
01007 for(iCell=0; iCell<fNoAct; iCell++) {
01008 sum = sum + ( (TFoamCell *) (fCellsAct->At(iCell)) )->GetPrim()/fPrime;
01009 fPrimAcu[iCell]=sum;
01010 }
01011
01012 }
01013
01014
01015 void TFoam::ResetPseRan(TRandom *PseRan)
01016 {
01017
01018
01019
01020
01021
01022 if(fPseRan) {
01023 Info("ResetPseRan", "Resetting random number generator \n");
01024 delete fPseRan;
01025 }
01026 SetPseRan(PseRan);
01027 }
01028
01029
01030 void TFoam::SetRho(TFoamIntegrand *fun)
01031 {
01032
01033
01034
01035
01036 if (fun)
01037 fRho=fun;
01038 else
01039 Error("SetRho", "Bad function \n" );
01040 }
01041
01042
01043 void TFoam::ResetRho(TFoamIntegrand *fun)
01044 {
01045
01046
01047
01048
01049
01050
01051
01052 if(fRho) {
01053 Info("ResetRho", "!!! Resetting distribution function !!!\n");
01054 delete fRho;
01055 }
01056 SetRho(fun);
01057 }
01058
01059
01060 void TFoam::SetRhoInt(void *fun)
01061 {
01062
01063
01064
01065
01066
01067
01068 const char *namefcn = gCint->Getp2f2funcname(fun);
01069 if(namefcn) {
01070 fMethodCall=new TMethodCall();
01071 fMethodCall->InitWithPrototype(namefcn, "Int_t, Double_t *");
01072 }
01073 fRho=0;
01074 }
01075
01076
01077 Double_t TFoam::Eval(Double_t *xRand)
01078 {
01079
01080
01081
01082 Double_t result;
01083
01084 if(!fRho) {
01085 Long_t paramArr[3];
01086 paramArr[0]=(Long_t)fDim;
01087 paramArr[1]=(Long_t)xRand;
01088 fMethodCall->SetParamPtrs(paramArr);
01089 fMethodCall->Execute(result);
01090 } else {
01091 result=fRho->Density(fDim,xRand);
01092 }
01093
01094 return result;
01095 }
01096
01097
01098 void TFoam::GenerCel2(TFoamCell *&pCell)
01099 {
01100
01101
01102
01103
01104 Long_t lo, hi, hit;
01105 Double_t fhit, flo, fhi;
01106 Double_t random;
01107
01108 random=fPseRan->Rndm();
01109 lo = 0; hi =fNoAct-1;
01110 flo = fPrimAcu[lo]; fhi=fPrimAcu[hi];
01111 while(lo+1<hi) {
01112 hit = lo + (Int_t)( (hi-lo)*(random-flo)/(fhi-flo)+0.5);
01113 if (hit<=lo)
01114 hit = lo+1;
01115 else if(hit>=hi)
01116 hit = hi-1;
01117 fhit=fPrimAcu[hit];
01118 if (fhit>random) {
01119 hi = hit;
01120 fhi = fhit;
01121 } else {
01122 lo = hit;
01123 flo = fhit;
01124 }
01125 }
01126 if (fPrimAcu[lo]>random)
01127 pCell = (TFoamCell *) fCellsAct->At(lo);
01128 else
01129 pCell = (TFoamCell *) fCellsAct->At(hi);
01130 }
01131
01132
01133
01134 void TFoam::MakeEvent(void)
01135 {
01136
01137
01138
01139
01140
01141
01142 Int_t j;
01143 Double_t wt,dx,mcwt;
01144 TFoamCell *rCell;
01145
01146
01147 ee0:
01148 GenerCel2(rCell);
01149
01150 MakeAlpha();
01151
01152 TFoamVect cellPosi(fDim); TFoamVect cellSize(fDim);
01153 rCell->GetHcub(cellPosi,cellSize);
01154 for(j=0; j<fDim; j++)
01155 fMCvect[j]= cellPosi[j] +fAlpha[j]*cellSize[j];
01156 dx = rCell->GetVolume();
01157
01158
01159 wt=dx*Eval(fMCvect);
01160
01161 mcwt = wt / rCell->GetPrim();
01162 fNCalls++;
01163 fMCwt = mcwt;
01164
01165 fSumWt += mcwt;
01166 fSumWt2 += mcwt*mcwt;
01167 fNevGen++;
01168 fWtMax = TMath::Max(fWtMax, mcwt);
01169 fWtMin = TMath::Min(fWtMin, mcwt);
01170 fMCMonit->Fill(mcwt);
01171 fHistWt->Fill(mcwt,1.0);
01172
01173 if(fOptRej == 1) {
01174 Double_t random;
01175 random=fPseRan->Rndm();
01176 if( fMaxWtRej*random > fMCwt) goto ee0;
01177 if( fMCwt<fMaxWtRej ) {
01178 fMCwt = 1.0;
01179 } else {
01180 fMCwt = fMCwt/fMaxWtRej;
01181 fSumOve += fMCwt-fMaxWtRej;
01182 }
01183 }
01184
01185 }
01186
01187
01188 void TFoam::GetMCvect(Double_t *MCvect)
01189 {
01190
01191
01192 for ( Int_t k=0 ; k<fDim ; k++) *(MCvect +k) = fMCvect[k];
01193 }
01194
01195
01196 Double_t TFoam::GetMCwt(void)
01197 {
01198
01199
01200 return(fMCwt);
01201 }
01202
01203 void TFoam::GetMCwt(Double_t &mcwt)
01204 {
01205
01206
01207 mcwt=fMCwt;
01208 }
01209
01210
01211 Double_t TFoam::MCgenerate(Double_t *MCvect)
01212 {
01213
01214
01215 MakeEvent();
01216 GetMCvect(MCvect);
01217 return(fMCwt);
01218 }
01219
01220
01221 void TFoam::GetIntegMC(Double_t &mcResult, Double_t &mcError)
01222 {
01223
01224
01225
01226
01227 Double_t mCerelat;
01228 mcResult = 0.0;
01229 mCerelat = 1.0;
01230 if (fNevGen>0) {
01231 mcResult = fPrime*fSumWt/fNevGen;
01232 mCerelat = sqrt( fSumWt2/(fSumWt*fSumWt) - 1/fNevGen);
01233 }
01234 mcError = mcResult *mCerelat;
01235 }
01236
01237
01238 void TFoam::GetIntNorm(Double_t& IntNorm, Double_t& Errel )
01239 {
01240
01241
01242
01243
01244
01245
01246 if(fOptRej == 1) {
01247 Double_t intMC,errMC;
01248 GetIntegMC(intMC,errMC);
01249 IntNorm = intMC;
01250 Errel = errMC;
01251 } else {
01252 IntNorm = fPrime;
01253 Errel = 0;
01254 }
01255 }
01256
01257
01258 void TFoam::GetWtParams(Double_t eps, Double_t &aveWt, Double_t &wtMax, Double_t &sigma)
01259 {
01260
01261
01262
01263 Double_t mCeff, wtLim;
01264 fMCMonit->GetMCeff(eps, mCeff, wtLim);
01265 wtMax = wtLim;
01266 aveWt = fSumWt/fNevGen;
01267 sigma = sqrt( fSumWt2/fNevGen -aveWt*aveWt );
01268 }
01269
01270
01271 void TFoam::Finalize(Double_t& IntNorm, Double_t& Errel)
01272 {
01273
01274
01275
01276 GetIntNorm(IntNorm,Errel);
01277 Double_t mcResult,mcError;
01278 GetIntegMC(mcResult,mcError);
01279 Double_t mCerelat= mcError/mcResult;
01280
01281 if(fChat>0) {
01282 Double_t eps = 0.0005;
01283 Double_t mCeff, mcEf2, wtMax, aveWt, sigma;
01284 GetWtParams(eps, aveWt, wtMax, sigma);
01285 mCeff=0;
01286 if(wtMax>0.0) mCeff=aveWt/wtMax;
01287 mcEf2 = sigma/aveWt;
01288 Double_t driver = fCells[0]->GetDriv();
01289
01290 BXOPE;
01291 BXTXT("****************************************");
01292 BXTXT("****** TFoam::Finalize ******");
01293 BXTXT("****************************************");
01294 BX1I(" NevGen",fNevGen, "Number of generated events in the MC generation ");
01295 BX1I(" nCalls",fNCalls, "Total number of function calls ");
01296 BXTXT("----------------------------------------");
01297 BX1F(" AveWt",aveWt, "Average MC weight ");
01298 BX1F(" WtMin",fWtMin, "Minimum MC weight (absolute) ");
01299 BX1F(" WtMax",fWtMax, "Maximum MC weight (absolute) ");
01300 BXTXT("----------------------------------------");
01301 BX1F(" XPrime",fPrime, "Primary total integral, R_prime ");
01302 BX1F(" XDiver",driver, "Driver total integral, R_loss ");
01303 BXTXT("----------------------------------------");
01304 BX2F(" IntMC", mcResult, mcError, "Result of the MC Integral");
01305 BX1F(" mCerelat", mCerelat, "Relative error of the MC integral ");
01306 BX1F(" <w>/WtMax",mCeff, "MC efficiency, acceptance rate");
01307 BX1F(" Sigma/<w>",mcEf2, "MC efficiency, variance/ave_wt");
01308 BX1F(" WtMax",wtMax, "WtMax(esp= 0.0005) ");
01309 BX1F(" Sigma",sigma, "variance of MC weight ");
01310 if(fOptRej==1) {
01311 Double_t avOve=fSumOve/fSumWt;
01312 BX1F("<OveW>/<W>",avOve, "Contrib. of events wt>MaxWtRej");
01313 }
01314 BXCLO;
01315 }
01316 }
01317
01318
01319 void TFoam::SetInhiDiv(Int_t iDim, Int_t InhiDiv)
01320 {
01321
01322
01323
01324
01325
01326
01327
01328
01329
01330
01331 if(fDim==0) Error("TFoam","SetInhiDiv: fDim=0 \n");
01332 if(fInhiDiv == 0) {
01333 fInhiDiv = new Int_t[ fDim ];
01334 for(Int_t i=0; i<fDim; i++) fInhiDiv[i]=0;
01335 }
01336
01337 if( ( 0<=iDim) && (iDim<fDim)) {
01338 fInhiDiv[iDim] = InhiDiv;
01339 } else
01340 Error("SetInhiDiv:","Wrong iDim \n");
01341 }
01342
01343
01344 void TFoam::SetXdivPRD(Int_t iDim, Int_t len, Double_t xDiv[])
01345 {
01346
01347
01348
01349
01350
01351
01352
01353
01354
01355
01356
01357
01358 Int_t i;
01359
01360 if(fDim<=0) Error("SetXdivPRD", "fDim=0 \n");
01361 if( len<1 ) Error("SetXdivPRD", "len<1 \n");
01362
01363 if(fXdivPRD == 0) {
01364 fXdivPRD = new TFoamVect*[fDim];
01365 for(i=0; i<fDim; i++) fXdivPRD[i]=0;
01366 }
01367
01368 if( ( 0<=iDim) && (iDim<fDim)) {
01369 fOptPRD =1;
01370 if( fXdivPRD[iDim] != 0)
01371 Error("SetXdivPRD", "Second allocation of XdivPRD not allowed \n");
01372 fXdivPRD[iDim] = new TFoamVect(len);
01373 for(i=0; i<len; i++) {
01374 (*fXdivPRD[iDim])[i]=xDiv[i];
01375 }
01376 } else {
01377 Error("SetXdivPRD", "Wrong iDim \n");
01378 }
01379
01380 cout<<" SetXdivPRD, idim= "<<iDim<<" len= "<<len<<" "<<endl;
01381 for(i=0; i<len; i++) {
01382 if (iDim < fDim) cout<< (*fXdivPRD[iDim])[i] <<" ";
01383 }
01384 cout<<endl;
01385 for(i=0; i<len; i++) cout<< xDiv[i] <<" ";
01386 cout<<endl;
01387
01388 }
01389
01390
01391 void TFoam::CheckAll(Int_t level)
01392 {
01393
01394
01395
01396
01397
01398 Int_t errors, warnings;
01399 TFoamCell *cell;
01400 Long_t iCell;
01401
01402 errors = 0; warnings = 0;
01403 if (level==1) cout << "///////////////////////////// FOAM_Checks /////////////////////////////////" << endl;
01404 for(iCell=1; iCell<=fLastCe; iCell++) {
01405 cell = fCells[iCell];
01406
01407 if( ((cell->GetDau0()==0) && (cell->GetDau1()!=0) ) ||
01408 ((cell->GetDau1()==0) && (cell->GetDau0()!=0) ) ) {
01409 errors++;
01410 if (level==1) Error("CheckAll","ERROR: Cell's no %ld has only one daughter \n",iCell);
01411 }
01412 if( (cell->GetDau0()==0) && (cell->GetDau1()==0) && (cell->GetStat()==0) ) {
01413 errors++;
01414 if (level==1) Error("CheckAll","ERROR: Cell's no %ld has no daughter and is inactive \n",iCell);
01415 }
01416 if( (cell->GetDau0()!=0) && (cell->GetDau1()!=0) && (cell->GetStat()==1) ) {
01417 errors++;
01418 if (level==1) Error("CheckAll","ERROR: Cell's no %ld has two daughters and is active \n",iCell);
01419 }
01420
01421
01422 if( (cell->GetPare())!=fCells[0] ) {
01423 if ( (cell != cell->GetPare()->GetDau0()) && (cell != cell->GetPare()->GetDau1()) ) {
01424 errors++;
01425 if (level==1) Error("CheckAll","ERROR: Cell's no %ld parent not pointing to this cell\n ",iCell);
01426 }
01427 }
01428
01429
01430 if(cell->GetDau0()!=0) {
01431 if(cell != (cell->GetDau0())->GetPare()) {
01432 errors++;
01433 if (level==1) Error("CheckAll","ERROR: Cell's no %ld daughter 0 not pointing to this cell \n",iCell);
01434 }
01435 }
01436 if(cell->GetDau1()!=0) {
01437 if(cell != (cell->GetDau1())->GetPare()) {
01438 errors++;
01439 if (level==1) Error("CheckAll","ERROR: Cell's no %ld daughter 1 not pointing to this cell \n",iCell);
01440 }
01441 }
01442 }
01443
01444
01445 for(iCell=0; iCell<=fLastCe; iCell++) {
01446 cell = fCells[iCell];
01447 if( (cell->GetStat()==1) && (cell->GetDriv()==0) ) {
01448 warnings++;
01449 if(level==1) Warning("CheckAll", "Warning: Cell no. %ld is active but empty \n", iCell);
01450 }
01451 }
01452
01453 if(level==1){
01454 Info("CheckAll","Check has found %d errors and %d warnings \n",errors, warnings);
01455 }
01456 if(errors>0){
01457 Info("CheckAll","Check - found total %d errors \n",errors);
01458 }
01459 }
01460
01461
01462 void TFoam::PrintCells(void)
01463 {
01464
01465
01466 Long_t iCell;
01467
01468 for(iCell=0; iCell<=fLastCe; iCell++) {
01469 cout<<"Cell["<<iCell<<"]={ ";
01470
01471 cout<<endl;
01472 fCells[iCell]->Print("1");
01473 cout<<"}"<<endl;
01474 }
01475 }
01476
01477
01478 void TFoam::RootPlot2dim(Char_t *filename)
01479 {
01480
01481
01482
01483 ofstream outfile(filename, ios::out);
01484 Double_t x1,y1,x2,y2,x,y;
01485 Long_t iCell;
01486 Double_t offs =0.1;
01487 Double_t lpag =1-2*offs;
01488 outfile<<"{" << endl;
01489 outfile<<"cMap = new TCanvas(\"Map1\",\" Cell Map \",600,600);"<<endl;
01490
01491 outfile<<"TBox*a=new TBox();"<<endl;
01492 outfile<<"a->SetFillStyle(0);"<<endl;
01493 outfile<<"a->SetLineWidth(4);"<<endl;
01494 outfile<<"a->SetLineColor(2);"<<endl;
01495 outfile<<"a->DrawBox("<<offs<<","<<offs<<","<<(offs+lpag)<<","<<(offs+lpag)<<");"<<endl;
01496
01497 outfile<<"TText*t=new TText();"<<endl;
01498 outfile<<"t->SetTextColor(4);"<<endl;
01499 if(fLastCe<51)
01500 outfile<<"t->SetTextSize(0.025);"<<endl;
01501 else if(fLastCe<251)
01502 outfile<<"t->SetTextSize(0.015);"<<endl;
01503 else
01504 outfile<<"t->SetTextSize(0.008);"<<endl;
01505
01506 outfile<<"TBox*b=new TBox();"<<endl;
01507 outfile <<"b->SetFillStyle(0);"<<endl;
01508
01509 if(fDim==2 && fLastCe<=2000) {
01510 TFoamVect cellPosi(fDim); TFoamVect cellSize(fDim);
01511 outfile << "// =========== Rectangular cells ==========="<< endl;
01512 for(iCell=1; iCell<=fLastCe; iCell++) {
01513 if( fCells[iCell]->GetStat() == 1) {
01514 fCells[iCell]->GetHcub(cellPosi,cellSize);
01515 x1 = offs+lpag*( cellPosi[0]); y1 = offs+lpag*( cellPosi[1]);
01516 x2 = offs+lpag*(cellPosi[0]+cellSize[0]); y2 = offs+lpag*(cellPosi[1]+cellSize[1]);
01517
01518 if(fLastCe<=2000)
01519 outfile<<"b->DrawBox("<<x1<<","<<y1<<","<<x2<<","<<y2<<");"<<endl;
01520
01521 if(fLastCe<=250) {
01522 x = offs+lpag*(cellPosi[0]+0.5*cellSize[0]); y = offs+lpag*(cellPosi[1]+0.5*cellSize[1]);
01523 outfile<<"t->DrawText("<<x<<","<<y<<","<<"\""<<iCell<<"\""<<");"<<endl;
01524 }
01525 }
01526 }
01527 outfile<<"// ============== End Rectangles ==========="<< endl;
01528 }
01529
01530
01531 outfile << "}" << endl;
01532 outfile.close();
01533 }
01534
01535 void TFoam::LinkCells()
01536 {
01537
01538
01539 Info("LinkCells", "VOID function for backward compatibility \n");
01540 return;
01541 }
01542
01543
01544
01545
01546