stress.cxx

Go to the documentation of this file.
00001 // @(#)root/test:$Id: stress.cxx 37531 2010-12-10 20:38:06Z pcanal $
00002 // Author: Rene Brun   05/11/98
00003 
00004 /////////////////////////////////////////////////////////////////
00005 //
00006 //    R O O T   T E S T  S U I T E  and  B E N C H M A R K S
00007 //    ======================================================
00008 //
00009 // The suite of programs below test the essential parts of Root.
00010 // In particular, there is an extensive test of the I/O and Trees.
00011 // The test can be run in batch or with the interpreter.
00012 // You must run
00013 //   gmake  on Unix systems
00014 //   nmake  on Windows
00015 // To run in batch, do
00016 //   stress -b 1000  : with no parameters, run standard test with 1000 events
00017 //   stress -b 30:     run test with 30 events only
00018 //
00019 // To run interactively, do
00020 // root -b
00021 //  Root > .L stress.cxx
00022 //  Root > stress(1000)  run standard test with 1000 events
00023 //  Root > stress(30)    run with 30 events only
00024 //
00025 // The standard test with 1000 events will create several files.
00026 // The size of all files is around 100 Mbytes.
00027 // The test with 30 events only require around  20 Mbytes
00028 // NB: The test must be run with more than 10 events
00029 //
00030 // The tests runs sequentially 16 tests. Each test will produce
00031 // one line (Test OK or Test failed) with some result parameters.
00032 // At the end of the test a table is printed showing the global results
00033 // with the amount of I/O, Real Time and Cpu Time.
00034 // One single number (ROOTMARKS) is also calculated showing the relative
00035 // performance of your machine compared to a reference machine
00036 // a Pentium IV 2.4 Ghz) with 512 MBytes of memory
00037 // and 120 GBytes IDE disk.
00038 //
00039 // An example of output when all the tests run OK is shown below:
00040 // ******************************************************************
00041 // *  Starting  R O O T - S T R E S S test suite with 1000 events
00042 // ******************************************************************
00043 // Test  1 : Functions, Random Numbers, Histogram Fits............. OK
00044 // Test  2 : Check size & compression factor of a Root file........ OK
00045 // Test  3 : Purge, Reuse of gaps in TFile......................... OK
00046 // Test  4 : Test of 2-d histograms, functions, 2-d fits........... OK
00047 // Test  5 : Test graphics & Postscript............................ OK
00048 // Test  6 : Test subdirectories in a Root file.................... OK
00049 // Test  7 : TNtuple, selections, TCut, TCutG, TEventList.......... OK
00050 // Test  8 : Trees split and compression modes..................... OK
00051 // Test  9 : Analyze Event.root file of stress 8................... OK
00052 // Test 10 : Create 10 files starting from Event.root.............. OK
00053 // Test 11 : Test chains of Trees using the 10 files............... OK
00054 // Test 12 : Compare histograms of test 9 and 11................... OK
00055 // Test 13 : Test merging files of a chain......................... OK
00056 // Test 14 : Check correct rebuilt of Event.root in test 13........ OK
00057 // Test 15 : Divert Tree branches to separate files................ OK
00058 // Test 16 : CINT test (3 nested loops) with LHCb trigger.......... OK
00059 // ******************************************************************
00060 //*  Linux pcbrun.cern.ch 2.4.20 #1 Thu Jan 9 12:21:02 MET 2003
00061 //******************************************************************
00062 //stress    : Total I/O =  703.7 Mbytes, I =  535.2, O = 168.5
00063 //stress    : Compr I/O =  557.0 Mbytes, I =  425.1, O = 131.9
00064 //stress    : Real Time =  64.84 seconds Cpu Time =  61.00 seconds
00065 //******************************************************************
00066 //*  ROOTMARKS = 600.1   *  Root4.02/00   20041217/1146
00067 //******************************************************************
00068 // 
00069 //_____________________________batch only_____________________
00070 #ifndef __CINT__
00071 
00072 #include <stdlib.h>
00073 #include <TROOT.h>
00074 #include <TSystem.h>
00075 #include <TH1.h>
00076 #include <TH2.h>
00077 #include <TFile.h>
00078 #include <TMath.h>
00079 #include <TF1.h>
00080 #include <TF2.h>
00081 #include <TProfile.h>
00082 #include <TKey.h>
00083 #include <TCanvas.h>
00084 #include <TGraph.h>
00085 #include <TRandom.h>
00086 #include <TPostScript.h>
00087 #include <TNtuple.h>
00088 #include <TTreeCache.h>
00089 #include <TChain.h>
00090 #include <TCut.h>
00091 #include <TCutG.h>
00092 #include <TEventList.h>
00093 #include <TBenchmark.h>
00094 #include <TSystem.h>
00095 #include <TApplication.h>
00096 #include <TClassTable.h>
00097 #include "Event.h"
00098 
00099 void stress(Int_t nevent, Int_t style, Int_t printSubBenchmark, UInt_t portion );
00100 void stress1();
00101 void stress2();
00102 void stress3();
00103 void stress4();
00104 void stress5();
00105 void stress6();
00106 void stress7();
00107 void stress8(Int_t nevent);
00108 void stress9tree(TTree *tree, Int_t realTestNum);
00109 void stress9();
00110 void stress10();
00111 void stress11();
00112 void stress12(Int_t testid);
00113 void stress13();
00114 void stress14();
00115 void stress15();
00116 void stress16();
00117 void cleanup();
00118 
00119 
00120 int main(int argc, char **argv)
00121 {
00122    TApplication theApp("App", &argc, argv);
00123    gBenchmark = new TBenchmark();
00124    Int_t nevent = 1000;      // by default create 1000 events
00125    if (argc > 1)  nevent = atoi(argv[1]);
00126    Int_t style  = 1;        // by default the new branch style
00127    if (argc > 2) style  = atoi(argv[2]);
00128    Int_t printSubBench = kFALSE;
00129    if (argc > 3) printSubBench = atoi(argv[3]);
00130    Int_t portion = 65535;
00131    if (argc > 4) portion  = atoi(argv[4]);
00132    stress(nevent, style, printSubBench, portion);
00133    return 0;
00134 }
00135 
00136 #endif
00137 
00138 class TH1;
00139 class TTree;
00140 
00141 int gPrintSubBench = 0;
00142 
00143 //_______________________common part_________________________
00144 
00145 Double_t ntotin=0, ntotout=0;
00146 
00147 void stress(Int_t nevent, Int_t style = 1, 
00148             Int_t printSubBenchmark = kFALSE, UInt_t portion = 65535)
00149 {
00150    //Main control function invoking all test programs
00151    
00152    gPrintSubBench = printSubBenchmark;
00153    
00154    if (nevent < 11) nevent = 11; // must have at least 10 events
00155    //Delete all possible objects in memory (to execute stress several times)
00156    gROOT->GetListOfFunctions()->Delete();
00157    gROOT->GetList()->Delete();
00158 
00159    printf("******************************************************************\n");
00160    printf("*  Starting  R O O T - S T R E S S test suite with %d events\n",nevent);
00161    printf("******************************************************************\n");
00162    // select the branch style
00163    TTree::SetBranchStyle(style);
00164 
00165    //Run the standard test suite
00166    gBenchmark->Start("stress");
00167    if (portion&1) stress1();
00168    if (portion&2) stress2();
00169    if (portion&4) stress3();
00170    if (portion&8) stress4();
00171    if (portion&16) stress5();
00172    if (portion&32) stress6();
00173    if (portion&64) stress7();
00174    if (portion&128) stress8(nevent);
00175    if (portion&256) stress9();
00176    if (portion&512) stress10();
00177    if (portion&1024) stress11();
00178    if (portion&2048) stress12(12);
00179    if (portion&4096) stress13();
00180    if (portion&8192) stress14();
00181    if (portion&16384) stress15();
00182    if (portion&32768) stress16();
00183    gBenchmark->Stop("stress");
00184 
00185    cleanup();
00186 
00187    //Print table with results
00188    Bool_t UNIX = strcmp(gSystem->GetName(), "Unix") == 0;
00189    printf("******************************************************************\n");
00190    if (UNIX) {
00191       TString sp = gSystem->GetFromPipe("uname -a");
00192       sp.Resize(60);
00193       printf("*  SYS: %s\n",sp.Data());
00194       if (strstr(gSystem->GetBuildNode(),"Linux")) {
00195          sp = gSystem->GetFromPipe("lsb_release -d -s");
00196          printf("*  SYS: %s\n",sp.Data());
00197       }
00198       if (strstr(gSystem->GetBuildNode(),"Darwin")) {
00199          sp  = gSystem->GetFromPipe("sw_vers -productVersion");
00200          sp += " Mac OS X ";
00201          printf("*  SYS: %s\n",sp.Data());
00202       }
00203    } else {
00204       const char *os = gSystem->Getenv("OS");
00205       if (!os) printf("*  SYS: Windows 95\n");
00206       else     printf("*  SYS: %s %s \n",os,gSystem->Getenv("PROCESSOR_IDENTIFIER"));
00207    }
00208 
00209    printf("******************************************************************\n");
00210    Float_t mbtot = (Float_t)(ntotin+ntotout)/1000000.;
00211    Float_t mbin  = (Float_t)ntotin/1000000.;
00212    Float_t mbout = (Float_t)ntotout/1000000.;
00213    printf("stress    : Total I/O =%7.1f Mbytes, I =%7.1f, O =%6.1f\n",mbtot,mbin,mbout);
00214    Float_t mbin1  = (Float_t)(TFile::GetFileBytesRead()/1000000.);
00215    Float_t mbout1 = (Float_t)(TFile::GetFileBytesWritten()/1000000.);
00216    Float_t mbtot1 = mbin1+mbout1;
00217    printf("stress    : Compr I/O =%7.1f Mbytes, I =%7.1f, O =%6.1f\n",mbtot1,mbin1,mbout1);
00218    gBenchmark->Print("stress");
00219 #ifndef __CINT__
00220    Float_t cp_brun_30   = 12.73;
00221    Float_t cp_brun_1000 = 61.88;
00222 #else
00223    Float_t cp_brun_30   = 31.03;  //The difference is essentially coming from stress16
00224    Float_t cp_brun_1000 = 84.30;
00225 #endif
00226    Float_t cp_brun = cp_brun_1000 - (cp_brun_1000 - cp_brun_30)*(1000-nevent)/(1000-30);
00227    Float_t ct = gBenchmark->GetCpuTime("stress");
00228    Float_t rootmarks = 600*cp_brun/ct;
00229    printf("******************************************************************\n");
00230    printf("*  ROOTMARKS =%6.1f   *  Root%-8s  %d/%d\n",rootmarks,gROOT->GetVersion(),gROOT->GetVersionDate(),gROOT->GetVersionTime());
00231    printf("******************************************************************\n");
00232    
00233    delete gBenchmark;
00234 }
00235 
00236 //_______________________________________________________________
00237 Double_t f1int(Double_t *x, Double_t *p)
00238 {
00239    //Compute a function sum of 3 gaussians
00240    Double_t e1 = (x[0]-p[1])/p[2];
00241    Double_t e2 = (x[0]-p[4])/p[5];
00242    Double_t e3 = (x[0]-p[7])/p[8];
00243    Double_t f  = p[0]*TMath::Exp(-0.5*e1*e1)
00244                 +p[3]*TMath::Exp(-0.5*e2*e2)
00245                 +p[6]*TMath::Exp(-0.5*e3*e3);
00246    return f;
00247 }
00248 
00249 //_______________________________________________________________
00250 void Bprint(Int_t id, const char *title)
00251 {
00252   // Print test program number and its title
00253    const Int_t kMAX = 65;
00254    char header[80];
00255    sprintf(header,"Test %2d : %s",id,title);
00256    Int_t nch = strlen(header);
00257    for (Int_t i=nch;i<kMAX;i++) header[i] = '.';
00258    header[kMAX] = 0;
00259    header[kMAX-1] = ' ';
00260    printf("%s",header);
00261 }
00262 
00263 //_______________________________________________________________
00264 void stress1()
00265 {
00266    //Generate two functions supposed to produce the same result
00267    //One function "f1form" will be computed by the TFormula class
00268    //The second function "f1int" will be
00269    //   - compiled when running in batch mode
00270    //   - interpreted by CINT when running in interactive mode
00271 
00272    Bprint(1,"Functions, Random Numbers, Histogram Fits");
00273 
00274    //Start with a function inline expression (managed by TFormula)
00275    Double_t f1params[9] = {100,-3,3,60,0,0.5,40,4,0.7};
00276    TF1 *f1form = new TF1("f1form","gaus(0)+gaus(3)+gaus(6)",-10,10);
00277    f1form->SetParameters(f1params);
00278 
00279    //Create an histogram and fill it randomly with f1form
00280    gRandom->SetSeed(65539);
00281    TH1F *h1form = new TH1F("h1form","distribution from f1form",100,-10,10);
00282    TH1F *h1diff = (TH1F*)h1form->Clone();
00283    h1diff->SetName("h1diff");
00284    h1form->FillRandom("f1form",10000);
00285 
00286    //Fit h1form with original function f1form
00287    h1form->Fit("f1form","q0");
00288 
00289    //same operation with an interpreted function f1int
00290    TF1 *f1 = new TF1("f1int",f1int,-10,10,9);
00291    f1->SetParameters(f1params);
00292 
00293    //Create an histogram and fill it randomly with f1int
00294    gRandom->SetSeed(65539); //make sure we start with the same random numbers
00295    TH1F *h1int = new TH1F("h1int","distribution from f1int",100,-10,10);
00296    h1int->FillRandom("f1int",10000);
00297 
00298    //Fit h1int with original function f1int
00299    h1int->Fit("f1int","q0");
00300 
00301    //The difference between the two histograms must be null
00302    h1diff->Add(h1form, h1int, 1, -1);
00303    Double_t hdiff = h1diff->Integral(0,101);
00304 
00305    //Compare fitted parameters and value of integral of f1form in [-8,6]
00306    Int_t npar = f1form->GetNpar();
00307    Double_t pdiff, pdifftot = 0;
00308    for (Int_t i=0;i<npar;i++) {
00309       pdiff = (f1form->GetParameter(i) - f1->GetParameter(i))/f1form->GetParameter(i);
00310       pdifftot += TMath::Abs(pdiff);
00311    }
00312    // The integral in the range [-8,6] must be = 1923.74578
00313    Double_t rint = TMath::Abs(f1form->Integral(-8,6) - 1923.74578);
00314 
00315    //Some slight differences are authorized to take into account
00316    //different math libraries used by the compiler, CINT and TFormula
00317    Bool_t OK = kTRUE;
00318    if (hdiff > 0.1 || pdifftot > 2.e-3 || rint > 10) OK = kFALSE;
00319    if (OK) printf("OK\n");
00320    else    {
00321       printf("failed\n");
00322       printf("%-8s hdiff=%g, pdifftot=%g, rint=%g\n"," ",hdiff,pdifftot,rint);
00323    }
00324    if (gPrintSubBench) { printf("Test  1 : "); gBenchmark->Show("stress");gBenchmark->Start("stress"); }
00325    //Save all objects in a Root file (will be checked by stress2)
00326    TFile local("stress.root","recreate");
00327    f1form->Write();
00328    f1->Write();
00329    h1form->Write();
00330    h1int->Write();
00331    ntotout += local.GetBytesWritten();
00332    //do not close the file. should be done by the destructor automatically
00333    delete h1int;
00334    delete h1form;
00335    delete h1diff;
00336 }
00337 
00338 //_______________________________________________________________
00339 void stress2()
00340 {
00341    //check length and compression factor in stress.root
00342    Bprint(2,"Check size & compression factor of a Root file");
00343    TFile f("stress.root");
00344    Long64_t last = f.GetEND();
00345    Float_t comp = f.GetCompressionFactor();
00346 
00347    Bool_t OK = kTRUE;
00348    Long64_t lastgood = 9428;
00349    if (last <lastgood-200 || last > lastgood+200 || comp <2.0 || comp > 2.4) OK = kFALSE;
00350    if (OK) printf("OK\n");
00351    else    {
00352       printf("failed\n");
00353       printf("%-8s last =%lld, comp=%f\n"," ",last,comp);
00354    }
00355    if (gPrintSubBench) { printf("Test  2 : "); gBenchmark->Show("stress");gBenchmark->Start("stress"); }
00356 }
00357 
00358 //_______________________________________________________________
00359 void stress3()
00360 {
00361    //Open stress.root, read all objects, save 10 times and purge
00362    //This function tests the generation and reuse of gaps in files
00363 
00364    Bprint(3,"Purge, Reuse of gaps in TFile");
00365    TFile f("stress.root","update");
00366    f.ReadAll();
00367    for (Int_t i=0;i<10;i++) {
00368       f.Write();
00369    }
00370    f.Purge();
00371    f.Write();
00372 
00373    //check length and compression level in stress.root
00374    ntotin  += f.GetBytesRead();
00375    ntotout += f.GetBytesWritten();
00376    Long64_t last = f.GetEND();
00377    Float_t comp = f.GetCompressionFactor();
00378    Bool_t OK = kTRUE;
00379    Long64_t lastgood = 49203;
00380    if (last <lastgood-900 || last > lastgood+900 || comp <1.8 || comp > 2.4) OK = kFALSE;
00381    if (OK) printf("OK\n");
00382    else    {
00383       printf("failed\n");
00384       printf("%-8s last =%lld, comp=%f\n"," ",last,comp);
00385    }
00386    if (gPrintSubBench) { printf("Test  3 : "); gBenchmark->Show("stress");gBenchmark->Start("stress"); }
00387 }
00388 
00389 //_______________________________________________________________
00390 void stress4()
00391 {
00392 // Test of 2-d histograms, functions, 2-d fits
00393 
00394    Bprint(4,"Test of 2-d histograms, functions, 2-d fits");
00395 
00396    Double_t f2params[15] = {100,-3,3,-3,3,160,0,0.8,0,0.9,40,4,0.7,4,0.7};
00397    TF2 *f2form = new TF2("f2form","xygaus(0)+xygaus(5)+xygaus(10)",-10,10,-10,10);
00398    f2form->SetParameters(f2params);
00399 
00400    //Create an histogram and fill it randomly with f2form
00401    gRandom->SetSeed(65539);
00402    TH2F *h2form = new TH2F("h2form","distribution from f2form",40,-10,10,40,-10,10);
00403    Int_t nentries = 100000;
00404    h2form->FillRandom("f2form",nentries);
00405    //Fit h2form with original function f2form
00406    Float_t ratio = 4*nentries/100000;
00407    f2params[ 0] *= ratio;
00408    f2params[ 5] *= ratio;
00409    f2params[10] *= ratio;
00410    f2form->SetParameters(f2params);
00411    h2form->Fit("f2form","q0");
00412    //Update stress.root
00413    TFile f("stress.root","update");
00414    h2form->Write();
00415    f2form->Write();
00416 
00417    ntotin  += f.GetBytesRead();
00418    ntotout += f.GetBytesWritten();
00419 
00420    //Compare results of fit with expected parameters
00421    Double_t th0=395.085, th5=620.957, th10=151.308;
00422    Double_t dp0  = TMath::Abs((f2form->GetParameter(0) -th0)/th0);
00423    Double_t dp5  = TMath::Abs((f2form->GetParameter(5) -th5)/th5);
00424    Double_t dp10 = TMath::Abs((f2form->GetParameter(10)-th10)/th10);
00425    Bool_t OK = kTRUE;
00426    if (dp0 > 1.e-2 || dp5 > 1.e-2 || dp10 > 1.e-2) OK = kFALSE;
00427    if (OK) printf("OK\n");
00428    else    {
00429       printf("failed\n");
00430       printf("%-8s dp0=%g, dp5=%g, dp10=%g\n"," ",dp0,dp5,dp10);
00431       printf("%-8s  p0=%g,  p5=%g,  p10=%g\n"," ",f2form->GetParameter(0),f2form->GetParameter(5),f2form->GetParameter(10));
00432    }
00433    if (gPrintSubBench) { printf("Test  4 : "); gBenchmark->Show("stress");gBenchmark->Start("stress"); }
00434 }
00435 
00436 //_______________________________________________________________
00437 void stress5()
00438 {
00439 // Test of Postscript.
00440 // Make a complex picture. Verify number of lines on ps file
00441 // Testing automatically the graphics package is a complex problem.
00442 // The best way we have found is to generate a Postscript image
00443 // of a complex canvas containing many objects.
00444 // The number of lines in the ps file is compared with a reference run.
00445 // A few lines (up to 2 or 3) of difference may be expected because
00446 // Postscript works with floats. The date and time of the run are also
00447 // different.
00448 // You can also inspect visually the ps file with a ps viewer.
00449 
00450    Bprint(5,"Test graphics & Postscript");
00451 
00452    TCanvas *c1 = new TCanvas("c1","stress canvas",800,600);
00453    gROOT->LoadClass("TPostScript","Postscript");
00454    TPostScript ps("stress.ps",112);
00455 
00456    //Get objects generated in previous test
00457    TFile f("stress.root");
00458    TF1  *f1form = (TF1*)f.Get("f1form");
00459    TF2  *f2form = (TF2*)f.Get("f2form");
00460    TH1F *h1form = (TH1F*)f.Get("h1form");
00461    TH2F *h2form = (TH2F*)f.Get("h2form");
00462 
00463    //Divide the canvas in subpads. Plot with different options
00464    c1->Divide(2,2);
00465    c1->cd(1);
00466    f1form->Draw();
00467    c1->cd(2);
00468    h1form->Draw();
00469    c1->cd(3);
00470    h2form->Draw("box");
00471    f2form->Draw("cont1same");
00472    c1->cd(4);
00473    f2form->Draw("surf");
00474 
00475    ps.Close();
00476 
00477    //count number of lines in ps file
00478    FILE *fp = fopen("stress.ps","r");
00479    char line[260];
00480    Int_t nlines = 0;
00481    Int_t nlinesGood = 632;
00482    while (fgets(line,255,fp)) {
00483       nlines++;
00484    }
00485    fclose(fp);
00486    ntotin  += f.GetBytesRead();
00487    ntotout += f.GetBytesWritten();
00488    Bool_t OK = kTRUE;
00489    if (nlines < nlinesGood-40 || nlines > nlinesGood+40) OK = kFALSE;
00490    if (OK) printf("OK\n");
00491    else    {
00492       printf("failed\n");
00493       printf("%-8s nlines in stress.ps file = %d\n"," ",nlines);
00494    }
00495    delete c1;
00496    if (gPrintSubBench) { printf("Test  5 : "); gBenchmark->Show("stress");gBenchmark->Start("stress"); }
00497    
00498 }
00499 
00500 //_______________________________________________________________
00501 void stress6()
00502 {
00503 // Test subdirectories in a Root file
00504 // Create many TH1S histograms, make operations between them
00505 
00506    Bprint(6,"Test subdirectories in a Root file");
00507 
00508    TFile f("stress.root","update");
00509    // create a new subdirectory for each plane
00510    gRandom->SetSeed(65539);
00511    const Int_t nplanes = 10;
00512    const Int_t ncounters = 100;
00513    char dirname[50];
00514    char hname[20];
00515    char htitle[80];
00516    TH1S *hn[ncounters];
00517    TH1S *hs[ncounters];
00518    Int_t i,j,k,id;
00519    TH1F *hsumPlanes = new TH1F("hsumPlanes","Sum of all planes",100,0,100);
00520    //Create a subdirectory per detector plane
00521    for (i=0;i<nplanes;i++) {
00522       sprintf(dirname,"plane%d",i);
00523       TDirectory *cdplane = f.mkdir(dirname);
00524       if (cdplane == 0) continue;
00525       cdplane->cd();
00526       // create counter histograms
00527       for (j=0;j<ncounters;j++) {
00528          sprintf(hname,"h%d_%dN",i,j);
00529          sprintf(htitle,"hist for counter:%d in plane:%d North",j,i);
00530          hn[j] = new TH1S(hname,htitle,100,0,100);
00531          sprintf(hname,"h%d_%dS",i,j);
00532          sprintf(htitle,"hist for counter:%d in plane:%d South",j,i);
00533          hs[j] = new TH1S(hname,htitle,100,0,100);
00534       }
00535       // fill counter histograms randomly
00536       for (k=0;k<10000;k++) {
00537          id = Int_t(ncounters*gRandom->Rndm());
00538          hn[id]->Fill(gRandom->Gaus(60,10));
00539          hs[id]->Fill(gRandom->Gaus(40,5));
00540       }
00541       // Write all objects in directory in memory to disk
00542       cdplane->Write();
00543       // Delete all objects from memory
00544       cdplane->GetList()->Delete();
00545       f.cd();
00546    }
00547    // Now read back all objects from all subdirectories
00548    // Add North and south histograms in hsumPlanes
00549    for (i=0;i<nplanes;i++) {
00550       sprintf(dirname,"plane%d",i);
00551       f.cd(dirname);
00552       for (j=0;j<ncounters;j++) {
00553          sprintf(hname,"h%d_%dN",i,j);
00554          TH1S *hnorth; gDirectory->GetObject(hname,hnorth);
00555          sprintf(hname,"h%d_%dS",i,j);
00556          TH1S *hsouth; gDirectory->GetObject(hname,hsouth);
00557          if (hnorth == 0 || hsouth == 0) continue;
00558          hsumPlanes->Add(hnorth);
00559          hsumPlanes->Add(hsouth);
00560          delete hnorth; delete hsouth;
00561       }
00562       f.cd();    // change current directory to top
00563    }
00564    // Verify number of entries, rms and mean value
00565    ntotin  += f.GetBytesRead();
00566    ntotout += f.GetBytesWritten();
00567    Int_t nentries = (Int_t)hsumPlanes->GetEntries();
00568    Double_t rms   = hsumPlanes->GetRMS();
00569    Double_t mean  = hsumPlanes->GetMean();
00570    Int_t nentriesGood = 200000;
00571    Double_t rmsGood  = 12.745;
00572    Double_t meanGood = 50.01;
00573    Double_t diffrms  = TMath::Abs(rmsGood -rms)/rmsGood;
00574    Double_t diffmean = TMath::Abs(meanGood -mean)/meanGood;
00575    Bool_t OK = kTRUE;
00576    if (nentriesGood != nentries || diffrms > 1.e-2 || diffmean > 1.e-2) OK = kFALSE;
00577    if (OK) printf("OK\n");
00578    else    {
00579       printf("failed\n");
00580       printf("%-8s nentries=%d, diffmean=%g, diffrms=%g\n"," ",nentries,diffmean,diffrms);
00581    }
00582    if (gPrintSubBench) { printf("Test  6 : "); gBenchmark->Show("stress");gBenchmark->Start("stress"); }
00583 }
00584 
00585 //_______________________________________________________________
00586 void stress7()
00587 {
00588 // Test TNtuple class with several selection mechanisms
00589 // Test expression cuts
00590 // Test graphical cuts
00591 // Test event lists and operations on event lists
00592 // Compare results of TTree::Draw with results of an explict loop
00593 
00594    Bprint(7,"TNtuple, selections, TCut, TCutG, TEventList");
00595 
00596    TFile f("stress.root","update");
00597    // Create and fill a TNtuple
00598    gRandom->SetSeed(65539);
00599    TNtuple *ntuple = new TNtuple("ntuple","Demo ntuple","px:py:pz:random:i");
00600    Float_t px, py, pz;
00601    Int_t nall = 50000;
00602    Int_t i;
00603    for (i = 0; i < nall; i++) {
00604       gRandom->Rannor(px,py);
00605       pz = px*px + py*py;
00606       Float_t random = gRandom->Rndm(1);
00607       ntuple->Fill(px,py,pz,random,i);
00608    }
00609    ntuple->Write();
00610 
00611    // Create a graphical cut. Select only events in cut
00612    TCutG *cutg = new TCutG("cutg",9);
00613    cutg->SetVarX("py");
00614    cutg->SetVarY("px");
00615    cutg->SetPoint(0,-1.75713,2.46193);
00616    cutg->SetPoint(1,-2.58656,-0.786802);
00617    cutg->SetPoint(2,-0.179195,-0.101523);
00618    cutg->SetPoint(3,2.12702,-1.49746);
00619    cutg->SetPoint(4,2.2484,1.95431);
00620    cutg->SetPoint(5,0.630004,0.583756);
00621    cutg->SetPoint(6,-0.381495,2.28426);
00622    cutg->SetPoint(7,-1.27161,1.01523);
00623    cutg->SetPoint(8,-1.75713,2.46193);
00624    TH2F *hpxpy = new TH2F("hpxpy","px vx py with cutg",40,-4,4,40,-4,4);
00625    ntuple->Draw("px:py>>hpxpy","cutg","goff");
00626    Int_t npxpy = (Int_t)hpxpy->GetEntries();
00627    Int_t npxpyGood = 27918;
00628    hpxpy->Write();
00629    cutg->Write();
00630    delete cutg;
00631 
00632    // Fill a TEventList using the standard cut
00633    ntuple->Draw(">>elist","py<0 && pz>4 && random<0.5","goff");
00634    TEventList *elist; gDirectory->GetObject("elist",elist);
00635    // Fill hist htemp using the standard cut
00636    ntuple->Draw("px>>htemp0","py<0 && pz>4 && random<0.5","goff");
00637    TH1F *htemp0;  gDirectory->GetObject("htemp0",htemp0);
00638    Double_t pxmean0 = htemp0->GetMean();
00639    Double_t pxrms0  = htemp0->GetRMS();
00640 
00641    // Fill hist hcut using a TCut = the standard cut
00642    TCut cut1 = "py<0 && pz>4 && random<0.5";
00643    TCut vcut = "px>>hcut";
00644    ntuple->Draw(vcut,cut1,"goff");
00645    // Fill hist helist looping on the eventlist in TTree::Draw
00646    ntuple->SetEventList(elist);
00647    ntuple->Draw("px>>helist","","goff");
00648    ntuple->SetEventList(0);
00649    TH1F *hcut;   gDirectory->GetObject("hcut",hcut);
00650    TH1F *helist; gDirectory->GetObject("helist",helist);
00651    Int_t n1 = (Int_t)hcut->GetEntries();
00652    Int_t n2 = (Int_t)helist->GetEntries();
00653    htemp0->Write();
00654    cut1.Write();
00655    helist->Write();
00656    hcut->Write();
00657 
00658    // now loop on eventlist explicitely and fill helist again
00659    Float_t pxr;
00660    ntuple->SetBranchAddress("px",&pxr);
00661    TH1F *helistc = (TH1F*)helist->Clone();
00662    helistc->Reset();
00663    helistc->SetName("helistc");
00664    Int_t nlist = elist->GetN();
00665    for (i=0;i<nlist;i++) {
00666       Long64_t event = elist->GetEntry(i);
00667       ntuple->GetEntry(event);
00668       helistc->Fill(pxr);
00669    }
00670    Int_t n3 = (Int_t)helistc->GetEntries();
00671    Double_t pxmean2 = helistc->GetMean();
00672    Double_t pxrms2  = helistc->GetRMS();
00673    helistc->Write();
00674    elist->Write();
00675 
00676    // Generate several TEventlist objects + total and save them
00677    char elistname[20];
00678    char cutname[20];
00679    TEventList *el[10];
00680    TEventList *elistall = new TEventList("elistall","Sum of all cuts");
00681    for (i=0;i<10;i++) {
00682       sprintf(elistname,">>elist%d",i);
00683       sprintf(cutname,"i 10 == %d",i); cutname[1] ='%';
00684       ntuple->Draw(elistname,cutname,"goff");
00685       gDirectory->GetObject(&elistname[2],el[i]);
00686       el[i]->Write();
00687       elistall->Add(el[i]);
00688    }
00689    elistall->Write();
00690 
00691    // Read big list from file and check that the distribution with the list
00692    // correspond to all events (no cuts)
00693    delete ntuple;
00694    TNtuple *nt; gDirectory->GetObject("ntuple",nt);
00695    nt->SetBranchAddress("px",&pxr);
00696    TH1F *hpx = new TH1F("hpx","hpx",100,-3,3);
00697    nt->Draw("px>>hpx","","goff");
00698    TEventList *all; gDirectory->GetObject("elistall",all);
00699    nt->SetEstimate(nall); //must be done because the order in eventlist is different
00700    nt->SetEventList(all);
00701    TH1F *hall = (TH1F*)hpx->Clone();
00702    hall->SetName("hall");
00703    nt->Draw("px>>hall","","goff");
00704    // Take the difference between the two histograms. Must be empty
00705    //TH1F hcomp = (*hall) - (*hpx);
00706    //Double_t compsum = hcomp.GetSum();
00707    hall->Add(hpx,-1);
00708    Double_t compsum = hall->GetSum();
00709    ntotin  += f.GetBytesRead();
00710    ntotout += f.GetBytesWritten();
00711 
00712    // We can compare entries, means and rms
00713    Bool_t OK = kTRUE;
00714    if (n1 != n2 || n1 != n3 || n3 != nlist || nall !=elistall->GetN()
00715                 || npxpy != npxpyGood
00716                 || compsum != 0
00717                 || TMath::Abs(pxmean0-pxmean2) > 0.1
00718                 || TMath::Abs(pxrms0-pxrms2) > 0.01) OK = kFALSE;
00719    if (OK) printf("OK\n");
00720    else    {
00721       printf("failed\n");
00722       printf("%-8s n1=%d, n2=%d, n3=%d, elistallN=%d\n"," ",n1,n2,n3,elistall->GetN());
00723       printf("%-8s pxmean0=%g, pxmean2=%g, pxrms0=%g\n"," ",pxmean0,pxmean2,pxrms0);
00724       printf("%-8s pxrms2=%g, compsum=%g, npxpy=%d\n"," ",pxrms2,compsum,npxpy);
00725    }
00726    if (gPrintSubBench) { printf("Test  7 : "); gBenchmark->Show("stress");gBenchmark->Start("stress"); }
00727 }
00728 
00729 //_______________________________________________________________
00730 Int_t stress8read(Int_t nevent)
00731 {
00732 //  Read the event file
00733 //  Loop on all events in the file (reading everything).
00734 //  Count number of bytes read
00735 
00736    TFile *hfile = new TFile("Event.root");
00737    TTree *tree; hfile->GetObject("T",tree);
00738    Event *event = 0;
00739    tree->SetBranchAddress("event",&event);
00740    Int_t nentries = (Int_t)tree->GetEntries();
00741    Int_t nev = TMath::Max(nevent,nentries);
00742    //activate the treeCache
00743    Int_t cachesize = 10000000; //this is the default value: 10 MBytes
00744    tree->SetCacheSize(cachesize);
00745    TTreeCache::SetLearnEntries(1); //one entry is sufficient to learn
00746    TTreeCache *tc = (TTreeCache*)hfile->GetCacheRead();
00747    tc->SetEntryRange(0,nevent);
00748    Int_t nb = 0;
00749    for (Int_t ev = 0; ev < nev; ev++) {
00750       nb += tree->GetEntry(ev);        //read complete event in memory
00751    }
00752    ntotin  += hfile->GetBytesRead();
00753 
00754    delete event;
00755    delete hfile;
00756    return nb;
00757 }
00758 
00759 
00760 //_______________________________________________________________
00761 Int_t stress8write(Int_t nevent, Int_t comp, Int_t split)
00762 {
00763 //  Create the Event file in various modes
00764    // comp = compression level
00765    // split = 1 split mode, 0 = no split
00766 
00767    // Create the Event file, the Tree and the branches
00768    TFile *hfile = new TFile("Event.root","RECREATE","TTree benchmark ROOT file");
00769    hfile->SetCompressionLevel(comp);
00770 
00771    // Create one event
00772    Event *event = new Event();
00773 
00774    // Create a ROOT Tree and one superbranch
00775    TTree *tree = new TTree("T","An example of a ROOT tree");
00776    tree->SetAutoSave(100000000);  // autosave when 100 Mbytes written
00777    Int_t bufsize = 64000;
00778    if (split)  bufsize /= 4;
00779    tree->Branch("event", &event, bufsize,split);
00780 
00781    //Fill the Tree
00782    Int_t ev, nb=0, meanTracks=600;
00783    Float_t ptmin = 1;
00784    for (ev = 0; ev < nevent; ev++) {
00785       event->Build(ev,meanTracks,ptmin);
00786 
00787       nb += tree->Fill();  //fill the tree
00788    }
00789    hfile->Write();
00790    ntotout += hfile->GetBytesWritten();
00791    delete event;
00792    delete hfile;
00793    return nb;
00794 }
00795 
00796 
00797 //_______________________________________________________________
00798 void stress8(Int_t nevent)
00799 {
00800 //  Run the $ROOTSYS/test/Event program in several configurations.
00801 
00802    Bprint(8,"Trees split and compression modes");
00803 
00804   // First step: make sure the Event shared library exists
00805   // This test dynamic linking when running in interpreted mode
00806    if (!TClassTable::GetDict("Event")) {
00807       Bool_t UNIX = strcmp(gSystem->GetName(), "Unix") == 0;
00808       Int_t st1 = gSystem->Load("$(ROOTSYS)/test/libEvent");
00809       if (st1 == -1) {
00810          st1 = gSystem->Load("test/libEvent");
00811          if (st1 == -1) {
00812             printf("===>stress8 will try to build the libEvent library\n");
00813             if (UNIX) gSystem->Exec("(cd $ROOTSYS/test; make Event)");
00814             else      gSystem->Exec("(cd %ROOTSYS%\\test && nmake libEvent.dll)");
00815             st1 = gSystem->Load("$(ROOTSYS)/test/libEvent");
00816          }
00817       }
00818    }
00819 
00820    // Create the file not compressed, in no-split mode and read it back
00821    gRandom->SetSeed(65539);
00822    Int_t nbw0 = stress8write(100,0,0);
00823    Int_t nbr0 = stress8read(0);
00824    Event::Reset();
00825 
00826    // Create the file compressed, in no-split mode and read it back
00827    gRandom->SetSeed(65539);
00828    Int_t nbw1 = stress8write(100,1,0);
00829    Int_t nbr1 = stress8read(0);
00830    Event::Reset();
00831 
00832    // Create the file compressed, in split mode and read it back
00833    gRandom->SetSeed(65539);
00834    Int_t nbw2 = stress8write(nevent,1,9);
00835    Int_t nbr2 = stress8read(0);
00836    Event::Reset();
00837 
00838    Bool_t OK = kTRUE;
00839    if (nbw0 != nbr0 || nbw1 != nbr1 || nbw2 != nbr2) OK = kFALSE;
00840    if (nbw0 != nbw1) OK = kFALSE;
00841    if (OK) printf("OK\n");
00842    else    {
00843       printf("failed\n");
00844       printf("%-8s nbw0=%d, nbr0=%d, nbw1=%d\n"," ",nbw0,nbr0,nbw1);
00845       printf("%-8s nbr1=%d, nbw2=%d, nbr2=%d\n"," ",nbr1,nbw2,nbr2);
00846    }
00847    if (gPrintSubBench) { printf("Test  8 : "); gBenchmark->Show("stress");gBenchmark->Start("stress"); }
00848 }
00849 
00850 //_______________________________________________________________
00851 Int_t HistCompare(TH1 *h1, TH1 *h2)
00852 {
00853 // Compare histograms h1 and h2
00854 // Check number of entries, mean and rms
00855 // if means differ by more than 1/1000 of the range return -1
00856 // if rms differs in percent by more than 1/1000 return -2
00857 // Otherwise return difference of number of entries
00858 
00859    Int_t n1       = (Int_t)h1->GetEntries();
00860    Double_t mean1 = h1->GetMean();
00861    Double_t rms1  = h1->GetRMS();
00862    Int_t n2       = (Int_t)h2->GetEntries();
00863    Double_t mean2 = h2->GetMean();
00864    Double_t rms2  = h2->GetRMS();
00865    Float_t xrange = h1->GetXaxis()->GetXmax() - h1->GetXaxis()->GetXmin();
00866    if (TMath::Abs((mean1-mean2)/xrange) > 0.001*xrange) return -1;
00867    if (rms1 && TMath::Abs((rms1-rms2)/rms1) > 0.001)    return -2;
00868    return n1-n2;
00869 }
00870 
00871 //_______________________________________________________________
00872 void stress9tree(TTree *tree, Int_t realTestNum)
00873 {
00874 // Test selections via TreeFormula
00875 // tree is a TTree when called by stress9
00876 // tree is a TChain when called from stres11
00877 // This is a quite complex test checking the results of TTree::Draw
00878 // or TChain::Draw with an explicit loop on events.
00879 // Also a good test for the interpreter
00880 
00881    Event *event = 0;
00882    tree->SetBranchAddress("event",&event);
00883    gROOT->cd();
00884    TDirectory *hfile = gDirectory;
00885    Double_t nrsave = TFile::GetFileBytesRead();
00886 
00887    // Each tree->Draw generates an histogram
00888    tree->Draw("fNtrack>>hNtrack",    "","goff");
00889    tree->Draw("fNseg>>hNseg",        "","goff");
00890    tree->Draw("fTemperature>>hTemp", "","goff");
00891    tree->Draw("fH.GetMean()>>hHmean","","goff");
00892    tree->Draw("fTracks.fPx>>hPx","fEvtHdr.fEvtNum%10 == 0","goff");
00893    tree->Draw("fTracks.fPy>>hPy","fEvtHdr.fEvtNum%10 == 0","goff");
00894    tree->Draw("fTracks.fPz>>hPz","fEvtHdr.fEvtNum%10 == 0","goff");
00895    tree->Draw("fRandom>>hRandom","fEvtHdr.fEvtNum%10 == 1","goff");
00896    tree->Draw("fMass2>>hMass2",  "fEvtHdr.fEvtNum%10 == 1","goff");
00897    tree->Draw("fBx>>hBx",        "fEvtHdr.fEvtNum%10 == 1","goff");
00898    tree->Draw("fBy>>hBy",        "fEvtHdr.fEvtNum%10 == 1","goff");
00899    tree->Draw("fXfirst>>hXfirst","fEvtHdr.fEvtNum%10 == 2","goff");
00900    tree->Draw("fYfirst>>hYfirst","fEvtHdr.fEvtNum%10 == 2","goff");
00901    tree->Draw("fZfirst>>hZfirst","fEvtHdr.fEvtNum%10 == 2","goff");
00902    tree->Draw("fXlast>>hXlast",  "fEvtHdr.fEvtNum%10 == 3","goff");
00903    tree->Draw("fYlast>>hYlast",  "fEvtHdr.fEvtNum%10 == 3","goff");
00904    tree->Draw("fZlast>>hZlast",  "fEvtHdr.fEvtNum%10 == 3","goff");
00905    tree->Draw("fCharge>>hCharge","fPx < 0","goff");
00906    tree->Draw("fNpoint>>hNpoint","fPx < 0","goff");
00907    tree->Draw("fValid>>hValid",  "fPx < 0","goff");
00908 
00909    tree->Draw("fMatrix>>hFullMatrix","","goff");
00910    tree->Draw("fMatrix[][0]>>hColMatrix","","goff");
00911    tree->Draw("fMatrix[1][]>>hRowMatrix","","goff");
00912    tree->Draw("fMatrix[2][2]>>hCellMatrix","","goff");
00913 
00914    tree->Draw("fMatrix - fVertex>>hFullOper","","goff");
00915    tree->Draw("fMatrix[2][1] - fVertex[5][1]>>hCellOper","","goff");
00916    tree->Draw("fMatrix[][1]  - fVertex[5][1]>>hColOper","","goff");
00917    tree->Draw("fMatrix[2][]  - fVertex[5][2]>>hRowOper","","goff");
00918    tree->Draw("fMatrix[2][]  - fVertex[5][]>>hMatchRowOper","","goff");
00919    tree->Draw("fMatrix[][2]  - fVertex[][1]>>hMatchColOper","","goff");
00920    tree->Draw("fMatrix[][2]  - fVertex[][]>>hRowMatOper","","goff");
00921    tree->Draw("fMatrix[][2]  - fVertex[5][]>>hMatchDiffOper","","goff");
00922    tree->Draw("fMatrix[][]   - fVertex[][]>>hFullOper2","","goff");
00923 
00924    if (gPrintSubBench) { printf("\n"); printf("Test %2dD: ",realTestNum); gBenchmark->Show("stress");gBenchmark->Start("stress"); }
00925 
00926    ntotin  += TFile::GetFileBytesRead() -nrsave;
00927 
00928    //Get pointers to the histograms generated above
00929    TH1F *hNtrack = (TH1F*)hfile->Get("hNtrack");
00930    TH1F *hNseg   = (TH1F*)hfile->Get("hNseg");
00931    TH1F *hTemp   = (TH1F*)hfile->Get("hTemp");
00932    TH1F *hHmean  = (TH1F*)hfile->Get("hHmean");
00933    TH1F *hPx     = (TH1F*)hfile->Get("hPx");
00934    TH1F *hPy     = (TH1F*)hfile->Get("hPy");
00935    TH1F *hPz     = (TH1F*)hfile->Get("hPz");
00936    TH1F *hRandom = (TH1F*)hfile->Get("hRandom");
00937    TH1F *hMass2  = (TH1F*)hfile->Get("hMass2");
00938    TH1F *hBx     = (TH1F*)hfile->Get("hBx");
00939    TH1F *hBy     = (TH1F*)hfile->Get("hBy");
00940    TH1F *hXfirst = (TH1F*)hfile->Get("hXfirst");
00941    TH1F *hYfirst = (TH1F*)hfile->Get("hYfirst");
00942    TH1F *hZfirst = (TH1F*)hfile->Get("hZfirst");
00943    TH1F *hXlast  = (TH1F*)hfile->Get("hXlast");
00944    TH1F *hYlast  = (TH1F*)hfile->Get("hYlast");
00945    TH1F *hZlast  = (TH1F*)hfile->Get("hZlast");
00946    TH1F *hCharge = (TH1F*)hfile->Get("hCharge");
00947    TH1F *hNpoint = (TH1F*)hfile->Get("hNpoint");
00948    TH1F *hValid  = (TH1F*)hfile->Get("hValid");
00949 
00950    TH1F *hFullMatrix    = (TH1F*)hfile->Get("hFullMatrix");
00951    TH1F *hColMatrix     = (TH1F*)hfile->Get("hColMatrix");
00952    TH1F *hRowMatrix     = (TH1F*)hfile->Get("hRowMatrix");
00953    TH1F *hCellMatrix    = (TH1F*)hfile->Get("hCellMatrix");
00954    TH1F *hFullOper      = (TH1F*)hfile->Get("hFullOper");
00955    TH1F *hCellOper      = (TH1F*)hfile->Get("hCellOper");
00956    TH1F *hColOper       = (TH1F*)hfile->Get("hColOper");
00957    TH1F *hRowOper       = (TH1F*)hfile->Get("hRowOper");
00958    TH1F *hMatchRowOper  = (TH1F*)hfile->Get("hMatchRowOper");
00959    TH1F *hMatchColOper  = (TH1F*)hfile->Get("hMatchColOper");
00960    TH1F *hRowMatOper    = (TH1F*)hfile->Get("hRowMatOper");
00961    TH1F *hMatchDiffOper = (TH1F*)hfile->Get("hMatchDiffOper");
00962    TH1F *hFullOper2     = (TH1F*)hfile->Get("hFullOper2");
00963 
00964    //We make clones of the generated histograms
00965    //We set new names and reset the clones.
00966    //We want to have identical histogram limits
00967    TH1F *bNtrack = (TH1F*)hNtrack->Clone(); bNtrack->SetName("bNtrack"); bNtrack->Reset();
00968    TH1F *bNseg   = (TH1F*)hNseg->Clone();   bNseg->SetName("bNseg");     bNseg->Reset();
00969    TH1F *bTemp   = (TH1F*)hTemp->Clone();   bTemp->SetName("bTemp");     bTemp->Reset();
00970    TH1F *bHmean  = (TH1F*)hHmean->Clone();  bHmean->SetName("bHmean");   bHmean->Reset();
00971    TH1F *bPx     = (TH1F*)hPx->Clone();     bPx->SetName("bPx");         bPx->Reset();
00972    TH1F *bPy     = (TH1F*)hPy->Clone();     bPy->SetName("bPy");         bPy->Reset();
00973    TH1F *bPz     = (TH1F*)hPz->Clone();     bPz->SetName("bPz");         bPz->Reset();
00974    TH1F *bRandom = (TH1F*)hRandom->Clone(); bRandom->SetName("bRandom"); bRandom->Reset();
00975    TH1F *bMass2  = (TH1F*)hMass2->Clone();  bMass2->SetName("bMass2");   bMass2->Reset();
00976    TH1F *bBx     = (TH1F*)hBx->Clone();     bBx->SetName("bBx");         bBx->Reset();
00977    TH1F *bBy     = (TH1F*)hBy->Clone();     bBy->SetName("bBy");         bBy->Reset();
00978    TH1F *bXfirst = (TH1F*)hXfirst->Clone(); bXfirst->SetName("bXfirst"); bXfirst->Reset();
00979    TH1F *bYfirst = (TH1F*)hYfirst->Clone(); bYfirst->SetName("bYfirst"); bYfirst->Reset();
00980    TH1F *bZfirst = (TH1F*)hZfirst->Clone(); bZfirst->SetName("bZfirst"); bZfirst->Reset();
00981    TH1F *bXlast  = (TH1F*)hXlast->Clone();  bXlast->SetName("bXlast");   bXlast->Reset();
00982    TH1F *bYlast  = (TH1F*)hYlast->Clone();  bYlast->SetName("bYlast");   bYlast->Reset();
00983    TH1F *bZlast  = (TH1F*)hZlast->Clone();  bZlast->SetName("bZlast");   bZlast->Reset();
00984    TH1F *bCharge = (TH1F*)hCharge->Clone(); bCharge->SetName("bCharge"); bCharge->Reset();
00985    TH1F *bNpoint = (TH1F*)hNpoint->Clone(); bNpoint->SetName("bNpoint"); bNpoint->Reset();
00986    TH1F *bValid  = (TH1F*)hValid->Clone();  bValid->SetName("bValid");   bValid->Reset();
00987 
00988    TH1F *bFullMatrix    =(TH1F*)hFullMatrix->Clone();    bFullMatrix->SetName("bFullMatrix");       bFullMatrix->Reset();
00989    TH1F *bColMatrix    = (TH1F*)hColMatrix->Clone();     bColMatrix->SetName("bColMatrix");         bColMatrix->Reset();
00990    TH1F *bRowMatrix    = (TH1F*)hRowMatrix->Clone();     bRowMatrix->SetName("bRowMatrix");         bRowMatrix->Reset();
00991    TH1F *bCellMatrix   = (TH1F*)hCellMatrix->Clone();    bCellMatrix->SetName("bCellMatrix");       bCellMatrix->Reset();
00992    TH1F *bFullOper     = (TH1F*)hFullOper->Clone();      bFullOper->SetName("bFullOper");           bFullOper->Reset();
00993    TH1F *bCellOper     = (TH1F*)hCellOper->Clone();      bCellOper->SetName("bCellOper");           bCellOper->Reset();
00994    TH1F *bColOper      = (TH1F*)hColOper->Clone();       bColOper->SetName("bColOper");             bColOper->Reset();
00995    TH1F *bRowOper      = (TH1F*)hRowOper->Clone();       bRowOper->SetName("bRowOper");             bRowOper->Reset();
00996    TH1F *bMatchRowOper = (TH1F*)hMatchRowOper->Clone();  bMatchRowOper->SetName("bMatchRowOper");   bMatchRowOper->Reset();
00997    TH1F *bMatchColOper = (TH1F*)hMatchColOper->Clone();  bMatchColOper->SetName("bMatchColOper");   bMatchColOper->Reset();
00998    TH1F *bRowMatOper   = (TH1F*)hRowMatOper->Clone();    bRowMatOper->SetName("bRowMatOper");       bRowMatOper->Reset();
00999    TH1F *bMatchDiffOper= (TH1F*)hMatchDiffOper->Clone(); bMatchDiffOper->SetName("bMatchDiffOper"); bMatchDiffOper->Reset();
01000    TH1F *bFullOper2    = (TH1F*)hFullOper2->Clone();     bFullOper2->SetName("bFullOper2");         bFullOper2->Reset();
01001 
01002    // Loop with user code on all events and fill the b histograms
01003    // The code below should produce identical results to the tree->Draw above
01004 
01005    TClonesArray *tracks = event->GetTracks();
01006    Int_t nev = (Int_t)tree->GetEntries();
01007    Int_t i, ntracks, evmod,i0,i1;
01008    Track *t;
01009    EventHeader *head;
01010    Int_t nbin = 0;
01011    for (Int_t ev=0;ev<nev;ev++) {
01012       nbin += tree->GetEntry(ev);
01013       head = event->GetHeader();
01014       evmod = head->GetEvtNum()%10;
01015       bNtrack->Fill(event->GetNtrack());
01016       bNseg->Fill(event->GetNseg());
01017       bTemp->Fill(event->GetTemperature());
01018       bHmean->Fill(event->GetHistogram()->GetMean());
01019       ntracks = event->GetNtrack();
01020       for(i0=0;i0<4;i0++) {
01021          for(i1=0;i1<4;i1++) {
01022             bFullMatrix->Fill(event->GetMatrix(i0,i1));
01023          }
01024          bColMatrix->Fill(event->GetMatrix(i0,0));
01025          bRowMatrix->Fill(event->GetMatrix(1,i0)); // done here because the matrix is square!
01026       }
01027       bCellMatrix->Fill(event->GetMatrix(2,2));
01028       if ( 5 < ntracks ) {
01029          t = (Track*)tracks->UncheckedAt(5);
01030          for(i0=0;i0<4;i0++) {
01031             for(i1=0;i1<4;i1++) {
01032             }
01033             bColOper->Fill( event->GetMatrix(i0,1) - t->GetVertex(1) );
01034             bRowOper->Fill( event->GetMatrix(2,i0) - t->GetVertex(2) );
01035          }
01036          for(i0=0;i0<3;i0++) {
01037             bMatchRowOper->Fill( event->GetMatrix(2,i0) - t->GetVertex(i0) );
01038             bMatchDiffOper->Fill( event->GetMatrix(i0,2) - t->GetVertex(i0) );
01039          }
01040          bCellOper->Fill( event->GetMatrix(2,1) - t->GetVertex(1) );
01041       }
01042       for (i=0;i<ntracks;i++) {
01043          t = (Track*)tracks->UncheckedAt(i);
01044          if (evmod == 0) bPx->Fill(t->GetPx());
01045          if (evmod == 0) bPy->Fill(t->GetPy());
01046          if (evmod == 0) bPz->Fill(t->GetPz());
01047          if (evmod == 1) bRandom->Fill(t->GetRandom());
01048          if (evmod == 1) bMass2->Fill(t->GetMass2());
01049          if (evmod == 1) bBx->Fill(t->GetBx());
01050          if (evmod == 1) bBy->Fill(t->GetBy());
01051          if (evmod == 2) bXfirst->Fill(t->GetXfirst());
01052          if (evmod == 2) bYfirst->Fill(t->GetYfirst());
01053          if (evmod == 2) bZfirst->Fill(t->GetZfirst());
01054          if (evmod == 3) bXlast->Fill(t->GetXlast());
01055          if (evmod == 3) bYlast->Fill(t->GetYlast());
01056          if (evmod == 3) bZlast->Fill(t->GetZlast());
01057          if (t->GetPx() < 0) {
01058             bCharge->Fill(t->GetCharge());
01059             bNpoint->Fill(t->GetNpoint());
01060             bValid->Fill(t->GetValid());
01061          }
01062          if (i<4) {
01063             for(i1=0;i1<3;i1++) { // 3 is the min of the 2nd dim of Matrix and Vertex
01064                bFullOper ->Fill( event->GetMatrix(i,i1) - t->GetVertex(i1) );
01065                bFullOper2->Fill( event->GetMatrix(i,i1) - t->GetVertex(i1) );
01066                bRowMatOper->Fill( event->GetMatrix(i,2) - t->GetVertex(i1) );
01067             }
01068             bMatchColOper->Fill( event->GetMatrix(i,2) - t->GetVertex(1) );
01069          }
01070       }
01071    }
01072 
01073    // Compare h and b histograms
01074    Int_t cNtrack = HistCompare(hNtrack,bNtrack);
01075    Int_t cNseg   = HistCompare(hNseg,bNseg);
01076    Int_t cTemp   = HistCompare(hTemp,bTemp);
01077    Int_t cHmean  = HistCompare(hHmean,bHmean);
01078    Int_t cPx     = HistCompare(hPx,bPx);
01079    Int_t cPy     = HistCompare(hPy,bPy);
01080    Int_t cPz     = HistCompare(hPz,bPz);
01081    Int_t cRandom = HistCompare(hRandom,bRandom);
01082    Int_t cMass2  = HistCompare(hMass2,bMass2);
01083    Int_t cBx     = HistCompare(hBx,bBx);
01084    Int_t cBy     = HistCompare(hBy,bBy);
01085    Int_t cXfirst = HistCompare(hXfirst,bXfirst);
01086    Int_t cYfirst = HistCompare(hYfirst,bYfirst);
01087    Int_t cZfirst = HistCompare(hZfirst,bZfirst);
01088    Int_t cXlast  = HistCompare(hXlast,bXlast);
01089    Int_t cYlast  = HistCompare(hYlast,bYlast);
01090    Int_t cZlast  = HistCompare(hZlast,bZlast);
01091    Int_t cCharge = HistCompare(hCharge,bCharge);
01092    Int_t cNpoint = HistCompare(hNpoint,bNpoint);
01093    Int_t cValid  = HistCompare(hValid,bValid);
01094 
01095    Int_t cFullMatrix   = HistCompare(hFullMatrix,bFullMatrix);
01096    Int_t cColMatrix    = HistCompare(hColMatrix,bColMatrix);
01097    Int_t cRowMatrix    = HistCompare(hRowMatrix,bRowMatrix);
01098    Int_t cCellMatrix   = HistCompare(hCellMatrix,bCellMatrix);
01099    Int_t cFullOper     = HistCompare(hFullOper,bFullOper);
01100    Int_t cCellOper     = HistCompare(hCellOper,bCellOper);
01101    Int_t cColOper      = HistCompare(hColOper,bColOper);
01102    Int_t cRowOper      = HistCompare(hRowOper,bRowOper);
01103    Int_t cMatchRowOper = HistCompare(hMatchRowOper,bMatchRowOper);
01104    Int_t cMatchColOper = HistCompare(hMatchColOper,bMatchColOper);
01105    Int_t cRowMatOper   = HistCompare(hRowMatOper,bRowMatOper);
01106    Int_t cMatchDiffOper= HistCompare(hMatchDiffOper,bMatchDiffOper);
01107    Int_t cFullOper2    = HistCompare(hFullOper2,bFullOper2);
01108 
01109    delete event;
01110    Event::Reset();
01111    ntotin += nbin;
01112 
01113    if (gPrintSubBench) { 
01114       printf("Test %2dC: ",realTestNum); 
01115       gBenchmark->Show("stress");gBenchmark->Start("stress");
01116       // Since we disturbed the flow (due to the double benchmark printing),
01117       // let's repeat the header!
01118       printf("Test %2d : ",realTestNum);
01119    }
01120    
01121    Bool_t OK = kTRUE;
01122    if (cNtrack || cNseg   || cTemp  || cHmean || cPx    || cPy     || cPz) OK = kFALSE;
01123    if (cRandom || cMass2  || cBx    || cBy    || cXfirst|| cYfirst || cZfirst) OK = kFALSE;
01124    if (cXlast  || cYlast  || cZlast || cCharge|| cNpoint|| cValid) OK = kFALSE;
01125    if (cFullMatrix || cColMatrix || cRowMatrix || cCellMatrix || cFullOper ) OK = kFALSE;
01126    if (cCellOper || cColOper || cRowOper || cMatchRowOper || cMatchColOper ) OK = kFALSE;
01127    if (cRowMatOper || cMatchDiffOper || cFullOper2 ) OK = kFALSE;
01128    if (OK) printf("OK\n");
01129    else    {
01130       printf("failed\n");
01131       printf("%-8s cNtrak =%d, cNseg  =%d, cTemp  =%d, cHmean =%d\n"," ",cNtrack,cNseg,cTemp,cHmean);
01132       printf("%-8s cPx    =%d, cPy    =%d, cPz    =%d, cRandom=%d\n"," ",cPx,cPy,cPz,cRandom);
01133       printf("%-8s cMass2 =%d, cbx    =%d, cBy    =%d, cXfirst=%d\n"," ",cMass2,cBx,cBy,cXfirst);
01134       printf("%-8s cYfirst=%d, cZfirst=%d, cXlast =%d, cYlast =%d\n"," ",cYfirst,cZfirst,cXlast,cYlast);
01135       printf("%-8s cZlast =%d, cCharge=%d, cNpoint=%d, cValid =%d\n"," ",cZlast,cCharge,cNpoint,cValid);
01136       printf("%-8s cFullMatrix=%d, cColMatrix=%d, cRowMatrix=%d, cCellMatrix=%d\n"," ",cFullMatrix,cColMatrix,cRowMatrix,cCellMatrix);
01137       printf("%-8s cFullOper=%d, cCellOper=%d, cColOper=%d, cRowOper=%d\n"," ",cFullOper,cCellOper,cColOper,cRowOper);
01138       printf("%-8s cMatchRowOper=%d, cMatchColOper=%d, cRowMatOper=%d, cMatchDiffOper=%d\n"," ",cMatchRowOper,cMatchColOper,cRowMatOper,cMatchDiffOper);
01139       printf("%-8s cFullOper2=%d\n"," ",cFullOper2);
01140    }
01141 }
01142 
01143 //_______________________________________________________________
01144 void stress9()
01145 {
01146 // Analyse the file Event.root generated in the last part of test8
01147 
01148    Bprint(9,"Analyze Event.root file of stress 8");
01149 
01150    gROOT->GetList()->Delete();
01151    TFile *hfile = new TFile("Event.root");
01152    TTree *tree; hfile->GetObject("T",tree);
01153 
01154    stress9tree(tree,9);
01155 
01156    // Save test9 histograms
01157    TFile f("stress_test9.root","recreate");
01158    gROOT->GetList()->Write();
01159    gROOT->GetList()->Delete();
01160    ntotout += f.GetBytesWritten();
01161 
01162 
01163    delete hfile;
01164 }
01165 
01166 //_______________________________________________________________
01167 void stress10()
01168 {
01169 // Make 10 Trees starting from the Event.root tree.
01170 // Events for which event_number%10 == 0 go to Event_0.root
01171 // Events for which event_number%10 == 1 go to Event_1.root
01172 //...
01173 // Events for which event_number%10 == 9 go to Event_9.root
01174 
01175    Bprint(10,"Create 10 files starting from Event.root");
01176 
01177    TFile *hfile = new TFile("Event.root");
01178    if (hfile==0 || hfile->IsZombie()) {
01179       delete hfile;
01180       printf("failed\n");
01181       return;
01182    }
01183    TTree *tree; hfile->GetObject("T",tree);
01184 
01185    Event *event = 0;
01186    tree->SetBranchAddress("event",&event);
01187 
01188    // Create 10 clones of this tree
01189    char filename[20];
01190    TTree *chTree[10];
01191    TFile *chfile[10];
01192    Int_t file;
01193    for (file=0;file<10;file++) {
01194       sprintf(filename,"Event_%d.root",file);
01195       chfile[file] = new TFile(filename,"recreate");
01196       chTree[file] = (TTree*)tree->CloneTree(0);
01197    }
01198 
01199    // Fill the small trees
01200    Int_t nev = (Int_t)tree->GetEntries();
01201    Int_t evmod, nbin=0, nbout=0;
01202    EventHeader *head;
01203    for (Int_t ev=0;ev<nev;ev++) {
01204       nbin += tree->GetEntry(ev);
01205       head = event->GetHeader();
01206       evmod = head->GetEvtNum()%10;
01207       nbout += chTree[evmod]->Fill();
01208       event->Clear();
01209    }
01210    // save headers
01211    Int_t ntot = 0;
01212    for (file=0;file<10;file++) {
01213       ntot += (Int_t)chTree[file]->GetEntries();
01214       chfile[file]->Write();
01215       delete chfile[file];
01216    }
01217    delete event;
01218    delete hfile;
01219    Event::Reset();
01220    ntotin  += nbin;
01221    ntotout += nbout;
01222 
01223    //We compare the number of bytes read from the big file
01224    //with the total number of bytes written in the 10 small files
01225    Bool_t OK = kTRUE;
01226    if (nbin != nbout || nev != ntot) OK = kFALSE;
01227    if (OK) printf("OK\n");
01228    else    {
01229       printf("failed\n");
01230       printf("%-8s nbin=%d, nbout=%d, nev=%d, ntot=%d\n"," ",nbin,nbout,nev,ntot);
01231    }
01232    if (gPrintSubBench) { printf("Test 10 : "); gBenchmark->Show("stress");gBenchmark->Start("stress"); }
01233 }
01234 
01235 //_______________________________________________________________
01236 void stress11()
01237 {
01238 // Test chains of Trees
01239 // We make a TChain using the 10 files generated in test10
01240 // We expect the same results when analyzing the chain than
01241 // in the analysis of the original big file Event.root in test9.
01242 // Because TChain derives from TTree, we can use the same
01243 // analysis procedure "stress9tree"
01244 
01245    Bprint(11,"Test chains of Trees using the 10 files");
01246 
01247    gROOT->GetList()->Delete();
01248    TChain *chain = new TChain("T");
01249    char filename[20];
01250    Int_t file;
01251    for (file=0;file<10;file++) {
01252       sprintf(filename,"Event_%d.root",file);
01253       chain->Add(filename);
01254    }
01255 
01256    stress9tree(chain,11);
01257 
01258    // Save test11 histograms
01259    delete chain;
01260    TFile f("stress_test11.root","recreate");
01261    gROOT->GetList()->Write();
01262    gROOT->GetList()->Delete();
01263    ntotout += f.GetBytesWritten();
01264 }
01265 
01266 //_______________________________________________________________
01267 void stress12(Int_t testid)
01268 {
01269 // Compare histograms of stress9 with stress11
01270 
01271    if (testid == 12) Bprint(12,"Compare histograms of test 9 and 11");
01272 
01273    TFile f9("stress_test9.root");
01274    TFile f11("stress_test11.root");
01275    //Let's loop on all keys of second file
01276    //We expect to find the same keys in the original stress9 file
01277    TIter next(f11.GetListOfKeys());
01278    TKey *key;
01279    TH1F *h9, *h11;
01280    Int_t comp, ngood = 0;
01281    while ((key=(TKey*)next())) {
01282       if (strcmp(key->GetClassName(),"TH1F")) continue; //may be a TList of TStreamerInfo
01283       h9  = (TH1F*)f9.Get(key->GetName());
01284       h11 = (TH1F*)f11.Get(key->GetName());
01285       if (h9 == 0 || h11 == 0) continue;
01286       comp = HistCompare(h9,h11);
01287       if (comp == 0) ngood++;
01288    }
01289    ntotin += f9.GetBytesRead();
01290    ntotin += f11.GetBytesRead();
01291    Bool_t OK = kTRUE;
01292    if (ngood < 40) OK = kFALSE;
01293    if (OK) printf("OK\n");
01294    else    {
01295       printf("failed\n");
01296       printf("%-8s ngood=%d\n"," ",ngood);
01297    }
01298    if (gPrintSubBench) { printf("Test 12 : "); gBenchmark->Show("stress");gBenchmark->Start("stress"); }
01299 }
01300 
01301 //_______________________________________________________________
01302 void stress13()
01303 {
01304 // test of TChain::Merge
01305 // The 10 small Tree files generated in stress10 are again merged
01306 // into one single file.
01307 // Should be the same as the file generated in stress8, except
01308 // that events will be in a different order.
01309 // But global analysis histograms should be identical (checked by stress14)
01310 
01311    Bprint(13,"Test merging files of a chain");
01312 
01313    gROOT->GetList()->Delete();
01314    TChain *chain = new TChain("T");
01315    char filename[20];
01316    Int_t file;
01317    for (file=0;file<10;file++) {
01318       sprintf(filename,"Event_%d.root",file);
01319       chain->Add(filename);
01320    }
01321 
01322    chain->Merge("Event.root");
01323 
01324    Double_t chentries = chain->GetEntries();
01325    delete chain;
01326 
01327    Event::Reset();
01328    gROOT->GetList()->Delete();
01329 
01330    TFile f("Event.root");
01331    TTree *tree = (TTree*)f.Get("T");
01332    ntotin  += (Double_t)f.GetEND();
01333    ntotout += (Double_t)f.GetEND();
01334 
01335    Bool_t OK = kTRUE;
01336    if (chentries != tree->GetEntries()) OK = kFALSE;
01337    if (OK) printf("OK\n");
01338    else    {
01339       printf("failed\n");
01340    }
01341    if (gPrintSubBench) { printf("Test 13 : "); gBenchmark->Show("stress");gBenchmark->Start("stress"); }
01342 }
01343 
01344 //_______________________________________________________________
01345 void stress14()
01346 {
01347 // Verify that stress13 has correctly rebuild the original Event.root
01348 
01349    Bprint(14,"Check correct rebuilt of Event.root in test 13");
01350 
01351    stress12(14);
01352 }
01353 
01354 //_______________________________________________________________
01355 void stress15()
01356 {
01357 // Divert some branches to separate files
01358 
01359    Bprint(15,"Divert Tree branches to separate files");
01360 
01361    //Get old file, old tree and set top branch address
01362    //We want to copy only a few branches.
01363    TFile *oldfile = new TFile("Event.root");
01364    if (oldfile->IsZombie()) {
01365       printf("failed\n");
01366       return;
01367    }   
01368    TTree *oldtree; oldfile->GetObject("T",oldtree);
01369    Event *event   = 0;
01370    oldtree->SetBranchAddress("event",&event);
01371    oldtree->SetBranchStatus("*",0);
01372    oldtree->SetBranchStatus("event",1);
01373    oldtree->SetBranchStatus("fNtrack",1);
01374    oldtree->SetBranchStatus("fNseg",1);
01375    oldtree->SetBranchStatus("fH",1);
01376 
01377 
01378    //Create a new file + a clone of old tree header. Do not copy events
01379    TFile *newfile = new TFile("stress_small.root","recreate");
01380    TTree *newtree = oldtree->CloneTree(0);
01381 
01382    //Divert branch fH to a separate file and copy all events
01383    newtree->GetBranch("fH")->SetFile("stress_fH.root");
01384    newtree->CopyEntries(oldtree);
01385 
01386    newfile->Write();
01387    ntotin  += oldfile->GetBytesRead();
01388    ntotout += newfile->GetBytesWritten();
01389    delete event;
01390    delete newfile;
01391    delete oldfile;
01392    Event::Reset();
01393    gROOT->GetList()->Delete();
01394 
01395    // Open small file, histogram fNtrack and fH
01396    newfile = new TFile("stress_small.root");
01397    newfile->GetObject("T", newtree);
01398    newtree->Draw("fNtrack>>hNtrack","","goff");
01399    newtree->Draw("fH.GetMean()>>hHmean","","goff");
01400    TH1F *hNtrack; newfile->GetObject("hNtrack",hNtrack);
01401    TH1F *hHmean; newfile->GetObject("hHmean",hHmean);
01402    ntotin  += newfile->GetBytesRead();
01403 
01404    // Open old reference file of stress9
01405    oldfile = new TFile("stress_test9.root");
01406    if (oldfile->IsZombie()) {
01407       printf("failed\n");
01408       return;
01409    }
01410    TH1F *bNtrack; oldfile->GetObject("bNtrack",bNtrack);
01411    TH1F *bHmean;  oldfile->GetObject("bHmean",bHmean);
01412    Int_t cNtrack = HistCompare(hNtrack,bNtrack);
01413    Int_t cHmean  = HistCompare(hHmean, bHmean);
01414    delete newfile;
01415    delete oldfile;
01416    Event::Reset();
01417    gROOT->GetList()->Delete();
01418 
01419    Bool_t OK = kTRUE;
01420    if (cNtrack || cHmean) OK = kFALSE;
01421    if (OK) printf("OK\n");
01422    else    {
01423       printf("failed\n");
01424       printf("%-8s cNtrack=%d, cHmean=%d\n"," ",cNtrack,cHmean);
01425    }
01426    if (gPrintSubBench) { printf("Test 15 : "); gBenchmark->Show("stress");gBenchmark->Start("stress"); }
01427 }
01428 
01429 void stress16()
01430 {
01431 // Prototype trigger simulation for the LHCb experiment
01432 // This test nested loops with the interpreter.
01433 // Expected to run fast with the compiler, slow with the interpreter.
01434 // This code is extracted from an original macro by Hans Dijkstra (LHCb)
01435 // The program generates histograms and profile histograms.
01436 // A canvas with subpads containing the results is sent to Postscript.
01437 // We check graphics results by counting the number of lines in the ps file.
01438 
01439    Bprint(16,"CINT test (3 nested loops) with LHCb trigger");
01440 
01441    const int nbuf    = 153;    // buffer size
01442    const int nlev    = 4;      // number of trigger levels
01443    const int nstep   = 50000;  // number of steps
01444    const int itt[4]  = { 1000, 4000, 40000, 400000 }; // time needed per trigger
01445    const float a[4]  = { 0.25, 0.04, 0.25, 0 };       // acceptance/trigger (last always 0)
01446 
01447    int   i, il, istep, itim[192], itrig[192], it, im, ipass;
01448    float dead, sum[10];
01449 
01450    // create histogram and array of profile histograms
01451    TCanvas *c = new TCanvas("laten","latency simulation",700,600);
01452    gROOT->LoadClass("TPostScript","Postscript");
01453    TPostScript ps("stress_lhcb.ps",112);
01454    gRandom->SetSeed(65539);
01455    TFile f("stress_lhcb.root", "recreate");
01456    TH1F *pipe = new TH1F("pipe", "free in pipeline", nbuf+1, -0.5, nbuf+0.5);
01457    pipe->SetLineColor(2);
01458    pipe->SetFillColor(2);
01459 
01460    TProfile *hp[nlev+1];
01461    TProfile::Approximate();
01462    for (i = 0; i <= nlev; i++) {
01463       char s[64];
01464       sprintf(s, "buf%d", i);
01465       hp[i] = new TProfile(s, "in buffers", 1000, 0,nstep, -1., 1000.);
01466       hp[i]->SetLineColor(2);
01467    }
01468 
01469    dead   = 0;
01470    sum[0] = nbuf;
01471    for (i = 1; i <= nlev; i++) sum[i] = 0;
01472    for (i = 0; i < nbuf; i++) { itrig[i] = 0; itim[i] = 0; }
01473 
01474    for (istep = 0; istep < nstep; istep++) {
01475       // evaluate status of buffer
01476       pipe->Fill(sum[0]);
01477       if ((istep+1)%10 == 0) {
01478          for (i = 0; i <= nlev; i++)
01479             hp[i]->Fill((float)istep, sum[i], 1.);
01480       }
01481 
01482       ipass = 0;
01483       for (i = 0; i < nbuf; i++) {
01484          it = itrig[i];
01485          if (it >= 1) {
01486             // add 25 ns to all times
01487             itim[i] += 25;
01488             im = itim[i];
01489             // level decisions
01490             for (il = 0; il < nlev; il++) {
01491                if (it == il+1 && im > itt[il]) {
01492                   if (gRandom->Rndm() > a[il]) {
01493                      itrig[i] = -1;
01494                      sum[0]++;
01495                      sum[il+1]--;
01496                   } else {
01497                      itrig[i]++;
01498                      sum[il+1]--;
01499                      sum[il+2]++;
01500                   }
01501                }
01502             }
01503          } else if (ipass == 0) {
01504             itrig[i] = 1;
01505             itim[i]  = 25;
01506             sum[0]--;
01507             sum[1]++;
01508             ipass++;
01509          }
01510       }
01511       if (ipass == 0) dead++;
01512    }
01513 //   Float_t deadTime = 100.*dead/nstep;
01514 
01515    // View results in the canvas and make the Postscript file
01516 
01517    c->Divide(2,3);
01518    c->cd(1); pipe->Draw();
01519    c->cd(2); hp[0]->Draw();
01520    c->cd(3); hp[1]->Draw();
01521    c->cd(4); hp[2]->Draw();
01522    c->cd(5); hp[3]->Draw();
01523    c->cd(6); hp[4]->Draw();
01524    ps.Close();
01525 
01526    f.Write();
01527    ntotout += f.GetBytesWritten();
01528 
01529    // Check length of Postscript file
01530    FILE *fp = fopen("stress_lhcb.ps","r");
01531    char line[260];
01532    Int_t nlines = 0;
01533    Int_t nlinesGood = 2121;
01534    Bool_t counting = kFALSE;
01535    while (fgets(line,255,fp)) {
01536       if (counting) nlines++;
01537       if (strstr(line,"%%EndProlog")) counting = kTRUE;
01538    }
01539    fclose(fp);
01540    delete c;
01541    Bool_t OK = kTRUE;
01542    if (nlines < nlinesGood-100 || nlines > nlinesGood+100) OK = kFALSE;
01543    if (OK) printf("OK\n");
01544    else    {
01545       printf("failed\n");
01546       printf("%-8s nlines in stress_lhcb.ps file = %d\n"," ",nlines);
01547    }
01548    if (gPrintSubBench) { printf("Test 16 : "); gBenchmark->Show("stress");gBenchmark->Start("stress"); }
01549 }
01550 
01551 void cleanup()
01552 {
01553    gSystem->Unlink("Event.root");
01554    gSystem->Unlink("Event_0.root");
01555    gSystem->Unlink("Event_1.root");
01556    gSystem->Unlink("Event_2.root");
01557    gSystem->Unlink("Event_3.root");
01558    gSystem->Unlink("Event_4.root");
01559    gSystem->Unlink("Event_5.root");
01560    gSystem->Unlink("Event_6.root");
01561    gSystem->Unlink("Event_7.root");
01562    gSystem->Unlink("Event_8.root");
01563    gSystem->Unlink("Event_9.root");
01564    gSystem->Unlink("stress.ps");
01565    gSystem->Unlink("stress.root");
01566    gSystem->Unlink("stress_fH.root");
01567    gSystem->Unlink("stress_lhcb.ps");
01568    gSystem->Unlink("stress_lhcb.root");
01569    gSystem->Unlink("stress_small.root");
01570    gSystem->Unlink("stress_test9.root");
01571    gSystem->Unlink("stress_test11.root");
01572 }

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