TSelectorDraw.cxx

Go to the documentation of this file.
00001 // @(#)root/treeplayer:$Id: TSelectorDraw.cxx 37954 2011-02-02 18:30:33Z pcanal $
00002 // Author: Rene Brun   08/01/2003
00003 
00004 /*************************************************************************
00005  * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers.               *
00006  * All rights reserved.                                                  *
00007  *                                                                       *
00008  * For the licensing terms see $ROOTSYS/LICENSE.                         *
00009  * For the list of contributors see $ROOTSYS/README/CREDITS.             *
00010  *************************************************************************/
00011 
00012 //////////////////////////////////////////////////////////////////////////
00013 //                                                                      //
00014 // TSelectorDraw                                                        //
00015 //                                                                      //
00016 // A specialized TSelector for TTree::Draw.                             //
00017 //                                                                      //
00018 //////////////////////////////////////////////////////////////////////////
00019 
00020 #include "TSelectorDraw.h"
00021 #include "TROOT.h"
00022 #include "TH2.h"
00023 #include "TH3.h"
00024 #include "TView.h"
00025 #include "TGraph.h"
00026 #include "TPolyMarker3D.h"
00027 #include "TDirectory.h"
00028 #include "TVirtualPad.h"
00029 #include "TProfile.h"
00030 #include "TProfile2D.h"
00031 #include "TTreeFormulaManager.h"
00032 #include "TEnv.h"
00033 #include "TTree.h"
00034 #include "TCut.h"
00035 #include "TEntryList.h"
00036 #include "TEventList.h"
00037 #include "THLimitsFinder.h"
00038 #include "TStyle.h"
00039 #include "TClass.h"
00040 #include "TColor.h"
00041 
00042 ClassImp(TSelectorDraw)
00043 
00044 const Int_t kCustomHistogram = BIT(17);
00045 
00046 //______________________________________________________________________________
00047 TSelectorDraw::TSelectorDraw()
00048 {
00049    // Default selector constructor.
00050 
00051    fTree           = 0;
00052    fW              = 0;
00053    fValSize        = 4;
00054    fVal            = new Double_t*[fValSize];
00055    fVmin           = new Double_t[fValSize];
00056    fVmax           = new Double_t[fValSize];
00057    fNbins          = new Int_t[fValSize];
00058    fVarMultiple    = new Bool_t[fValSize];
00059    fVar            = new TTreeFormula*[fValSize];
00060    for(Int_t i=0;i<fValSize;++i){
00061       fVal[i] = 0;
00062       fVar[i] = 0;
00063    }
00064    fManager        = 0;
00065    fMultiplicity   = 0;
00066    fSelect         = 0;
00067    fSelectedRows   = 0;
00068    fDraw           = 0;
00069    fObject         = 0;
00070    fOldHistogram   = 0;
00071    fObjEval        = kFALSE;
00072    fSelectMultiple = kFALSE;
00073    fCleanElist     = kFALSE;
00074    fTreeElist      = 0;
00075    fAction         = 0;
00076    fNfill          = 0;
00077    fDimension      = 0;
00078    fOldEstimate    = 0;
00079    fForceRead      = 0;
00080    fWeight         = 1;
00081 }
00082 
00083 //______________________________________________________________________________
00084 TSelectorDraw::~TSelectorDraw()
00085 {
00086    // Selector destructor.
00087 
00088    ClearFormula();
00089    delete [] fVar;
00090    if(fVal){
00091       for(Int_t i=0;i<fValSize;++i)
00092          delete [] fVal[i];
00093       delete [] fVal;
00094    }
00095    if(fVmin) delete [] fVmin;
00096    if(fVmax) delete [] fVmax;
00097    if(fNbins) delete [] fNbins;
00098    if(fVarMultiple) delete [] fVarMultiple;
00099    if (fW)     delete [] fW;
00100 }
00101 
00102 //______________________________________________________________________________
00103 void TSelectorDraw::Begin(TTree *tree)
00104 {
00105    // Called everytime a loop on the tree(s) starts.
00106 
00107    SetStatus(0);
00108    ResetBit(kCustomHistogram);
00109    fSelectedRows   = 0;
00110    fTree = tree;
00111    fDimension = 0;
00112    fAction = 0;
00113 
00114    const char *varexp0   = fInput->FindObject("varexp")->GetTitle();
00115    const char *selection = fInput->FindObject("selection")->GetTitle();
00116    const char *option    = GetOption();
00117 
00118    TString  opt;
00119    char *hdefault = (char *)"htemp";
00120    char *varexp;
00121    Int_t i,j,hkeep;
00122    opt = option;
00123    opt.ToLower();
00124    fOldHistogram = 0;
00125    TEntryList *enlist = 0;
00126    TEventList *evlist = 0;
00127    TString htitle;
00128    Bool_t profile = kFALSE;
00129    Bool_t optSame = kFALSE;
00130    Bool_t optEnlist = kFALSE;
00131    Bool_t optpara = kFALSE;
00132    Bool_t optcandle = kFALSE;
00133    Bool_t opt5d = kFALSE;
00134    if (opt.Contains("same")) {
00135       optSame = kTRUE;
00136       opt.ReplaceAll("same","");
00137    }
00138    if (opt.Contains("entrylist")){
00139       optEnlist = kTRUE;
00140       opt.ReplaceAll("entrylist", "");
00141    }
00142    if (opt.Contains("para")){
00143       optpara = kTRUE;
00144       opt.ReplaceAll("para","");
00145    }
00146    if (opt.Contains("candle")) {
00147       optcandle = kTRUE;
00148       opt.ReplaceAll("candle","");
00149    }
00150    if (opt.Contains("gl5d")) {
00151       opt5d = kTRUE;
00152       opt.ReplaceAll("gl5d","");
00153    }
00154    TCut realSelection(selection);
00155    //input list - only TEntryList
00156    TEntryList *inElist = fTree->GetEntryList();
00157    evlist = fTree->GetEventList();
00158    if (evlist && inElist){
00159       //this is needed because the input entry list was created
00160       //by the fTree from the input TEventList and is owned by the fTree.
00161       //Calling GetEntryList function changes ownership and here
00162       //we want fTree to still delete this entry list
00163 
00164       inElist->SetBit(kCanDelete, kTRUE);
00165    }
00166    fCleanElist = kFALSE;
00167    fTreeElist = inElist;
00168 
00169    if ( inElist && inElist->GetReapplyCut() ) {
00170       realSelection *= inElist->GetTitle();
00171    }
00172 
00173    // what each variable should contain:
00174    //   varexp0   - original expression eg "a:b>>htest"
00175    //   hname     - name of new or old histogram
00176    //   hkeep     - flag if to keep new histogram
00177    //   hnameplus - flag if to add to current histo
00178    //   i         - length of variable expression stipped of everything after ">>"
00179    //   varexp    - variable expression stipped of everything after ">>"
00180    //   fOldHistogram     - pointer to hist hname
00181    //   elist     - pointer to selection list of hname
00182 
00183    Bool_t canRebin = kTRUE;
00184    if (optSame) canRebin = kFALSE;
00185 
00186    Int_t nbinsx=0, nbinsy=0, nbinsz=0;
00187    Double_t xmin=0, xmax=0, ymin=0, ymax=0, zmin=0, zmax=0;
00188 
00189    fObject  = 0;
00190    char *hname = 0;
00191    char *hnamealloc = 0;
00192    i = 0;
00193    if (varexp0 && strlen(varexp0)) {
00194       for(UInt_t k=strlen(varexp0)-1;k>0;k--) {
00195          if (varexp0[k]=='>' && varexp0[k-1]=='>') {
00196             i = (int)( &(varexp0[k-1]) - varexp0 );  //  length of varexp0 before ">>"
00197             hnamealloc = new char[strlen(&(varexp0[k+1]))+1];
00198             hname = hnamealloc;
00199             strcpy(hname,&(varexp0[k+1]));
00200             break;
00201          }
00202       }
00203    }
00204    //   char *hname = (char*)strstr(varexp0,">>");
00205    if (hname) {
00206       hkeep  = 1;
00207       varexp = new char[i+1];
00208       varexp[0] = 0; //necessary if i=0
00209       Bool_t hnameplus = kFALSE;
00210       while (*hname == ' ') hname++;
00211       if (*hname == '+') {
00212          hnameplus = kTRUE;
00213          hname++;
00214          while (*hname == ' ') hname++; //skip ' '
00215       }
00216       j = strlen(hname) - 1;   // skip ' '  at the end
00217       while (j) {
00218          if (hname[j] != ' ') break;
00219          hname[j] = 0;
00220          j--;
00221       }
00222 
00223       if (i) {
00224          strlcpy(varexp,varexp0,i+1); 
00225 
00226          Int_t mustdelete=0;
00227          SetBit(kCustomHistogram);
00228 
00229          // parse things that follow the name of the histo between '(' and ')'.
00230          // At this point hname contains the name of the specified histogram.
00231          //   Now the syntax is exended to handle an hname of the following format
00232          //   hname(nBIN [[,[xlow]][,xhigh]],...)
00233          //   so enclosed in brackets is the binning information, xlow, xhigh, and
00234          //   the same for the other dimensions
00235 
00236          char *pstart;    // pointer to '('
00237          char *pend;      // pointer to ')'
00238          char *cdummy;    // dummy pointer
00239          int ncomma;      // number of commas between '(' and ')', later number of arguments
00240          int ncols;       // number of columns in varexpr
00241          Double_t value;  // parsed value (by sscanf)
00242 
00243          const Int_t maxvalues=9;
00244 
00245          pstart= strchr(hname,'(');
00246          pend =  strchr(hname,')');
00247          if (pstart != 0 ) {  // found the bracket
00248 
00249             mustdelete=1;
00250 
00251             // check that there is only one open and close bracket
00252             if (pstart == strrchr(hname,'(')  &&  pend == strrchr(hname,')')) {
00253 
00254                // count number of ',' between '(' and ')'
00255                ncomma=0;
00256                cdummy = pstart;
00257                cdummy = strchr(&cdummy[1],',');
00258                while (cdummy != 0) {
00259                   cdummy = strchr(&cdummy[1],',');
00260                   ncomma++;
00261                }
00262 
00263                if (ncomma+1 > maxvalues) {
00264                   Error("DrawSelect","ncomma+1>maxvalues, ncomma=%d, maxvalues=%d",ncomma,maxvalues);
00265                   ncomma=maxvalues-1;
00266                }
00267 
00268                ncomma++; // number of arguments
00269                cdummy = pstart;
00270 
00271                //   number of columns
00272                ncols  = 1;
00273                for (j=0;j<i;j++) {
00274                   if (varexp[j] == ':'
00275                       && ! ( (j>0&&varexp[j-1]==':') || varexp[j+1]==':' )
00276                       ) {
00277                      ncols++;
00278                   }
00279                }
00280                if (ncols > 3 ) {  // max 3 columns
00281                   Error("DrawSelect","ncols > 3, ncols=%d",ncols);
00282                   ncols = 0;
00283                }
00284 
00285                // check dimensions before and after ">>"
00286                if (ncols*3 < ncomma) {
00287                   Error("DrawSelect","ncols*3 < ncomma ncols=%d, ncomma=%d",ncols,ncomma);
00288                   ncomma = ncols*3;
00289                }
00290 
00291                // scan the values one after the other
00292                for (j=0;j<ncomma;j++) {
00293                   cdummy++;  // skip '(' or ','
00294                   if (sscanf(cdummy," %lf ",&value) == 1) {
00295                      cdummy=strchr(&cdummy[1],',');
00296 
00297                      switch (j) {  // do certain settings depending on position of argument
00298                         case 0:  // binning x-axis
00299                            nbinsx = (Int_t)value;
00300                            if      (ncols<2) {
00301                               gEnv->SetValue("Hist.Binning.1D.x",nbinsx);
00302                            } else if (ncols<3) {
00303                               gEnv->SetValue("Hist.Binning.2D.x",nbinsx);
00304                               gEnv->SetValue("Hist.Binning.2D.Prof",nbinsx);
00305                            } else {
00306                               gEnv->SetValue("Hist.Binning.3D.x",nbinsx);
00307                               gEnv->SetValue("Hist.Binning.3D.Profx",nbinsx);
00308                            }
00309 
00310                            break;
00311                         case 1:  // lower limit x-axis
00312                            xmin = value;
00313                            break;
00314                         case 2:  // upper limit x-axis
00315                            xmax = value;
00316                            break;
00317                         case 3:  // binning y-axis
00318                            nbinsy = (Int_t)value;
00319                            if (ncols<3) gEnv->SetValue("Hist.Binning.2D.y",nbinsy);
00320                            else {
00321                               gEnv->SetValue("Hist.Binning.3D.y",nbinsy);
00322                               gEnv->SetValue("Hist.Binning.3D.Profy",nbinsy);
00323                            }
00324                            break;
00325                         case 4:  // lower limit y-axis
00326                            ymin = value;
00327                            break;
00328                         case 5:  // upper limit y-axis
00329                            ymax = value;
00330                            break;
00331                         case 6:  // binning z-axis
00332                            nbinsz = (Int_t)value;
00333                            gEnv->SetValue("Hist.Binning.3D.z",nbinsz);
00334                            break;
00335                         case 7:  // lower limit z-axis
00336                            zmin = value;
00337                            break;
00338                         case 8:  // upper limit z-axis
00339                            zmax = value;
00340                            break;
00341                         default:
00342                            Error("DrawSelect","j>8");
00343                            break;
00344                      }
00345                   }  // if sscanf == 1
00346                } // for j=0;j<ncomma;j++
00347             } else {
00348                Error("Begin","Two open or close brackets found, hname=%s",hname);
00349             }
00350 
00351             // fix up hname
00352             pstart[0]='\0';   // removes things after (and including) '('
00353          } // if '(' is found
00354 
00355          j = strlen(hname) - 1; // skip ' '  at the end
00356          while (j) {
00357             if (hname[j] != ' ') break; // skip ' '  at the end
00358             hname[j] = 0;
00359             j--;
00360          }
00361 
00362          TObject *oldObject = gDirectory->Get(hname);  // if hname contains '(...)' the return values is NULL, which is what we want
00363          fOldHistogram = oldObject ? dynamic_cast<TH1*>(oldObject) : 0;
00364 
00365          if (!fOldHistogram && oldObject && !oldObject->InheritsFrom(TH1::Class())) {
00366             Error("Begin","An object of type '%s' has the same name as the requested histo (%s)",oldObject->IsA()->GetName(),hname);
00367             SetStatus(-1);
00368             return;
00369          }
00370          if (fOldHistogram && !hnameplus) fOldHistogram->Reset();  // reset unless adding is wanted
00371 
00372          if (mustdelete) {
00373             if (gDebug) {
00374                Warning("Begin","Deleting old histogram, since (possibly new) limits and binnings have been given");
00375             }
00376             delete fOldHistogram; fOldHistogram=0;
00377          }
00378 
00379       } else {
00380          // make selection list (i.e. varexp0 starts with ">>")
00381          TObject *oldObject = gDirectory->Get(hname);
00382          if (optEnlist){
00383             //write into a TEntryList
00384             enlist = oldObject ? dynamic_cast<TEntryList*>(oldObject) : 0;
00385 
00386             if (!enlist && oldObject) {
00387                Error("Begin","An object of type '%s' has the same name as the requested event list (%s)",
00388                   oldObject->IsA()->GetName(),hname);
00389                SetStatus(-1);
00390                return;
00391             }
00392             if (!enlist) {
00393                enlist = new TEntryList(hname, realSelection.GetTitle());
00394             }
00395             if (enlist) {
00396                if (!hnameplus) {
00397                   if (enlist==inElist) {
00398                      // We have been asked to reset the input list!!
00399                      // Let's set it aside for now ...
00400                      inElist = new TEntryList(*enlist);
00401                      fCleanElist = kTRUE;
00402                      fTree->SetEntryList(inElist);
00403                   }
00404                   enlist->Reset();
00405                   enlist->SetTitle(realSelection.GetTitle());
00406                } else {
00407                   TCut old = enlist->GetTitle();
00408                   TCut upd = old || realSelection.GetTitle();
00409                   enlist->SetTitle(upd.GetTitle());
00410                }
00411             }
00412          }
00413          else {
00414             //write into a TEventList
00415             evlist = oldObject ? dynamic_cast<TEventList*>(oldObject) : 0;
00416 
00417             if (!evlist && oldObject) {
00418                Error("Begin","An object of type '%s' has the same name as the requested event list (%s)",
00419                   oldObject->IsA()->GetName(),hname);
00420                SetStatus(-1);
00421                return;
00422             }
00423             if (!evlist) {
00424                evlist = new TEventList(hname,realSelection.GetTitle(),1000,0);
00425             }
00426             if (evlist) {
00427                if (!hnameplus) {
00428                   if (evlist==fTree->GetEventList()) {
00429                      // We have been asked to reset the input list!!
00430                      // Let's set it aside for now ...
00431                      Error("Begin", "Input and output lists are the same!\n");
00432                      SetStatus(-1);
00433                      delete [] varexp;
00434                      return;
00435                   }
00436                   evlist->Reset();
00437                   evlist->SetTitle(realSelection.GetTitle());
00438                } else {
00439                   TCut old = evlist->GetTitle();
00440                   TCut upd = old || realSelection.GetTitle();
00441                   evlist->SetTitle(upd.GetTitle());
00442                }
00443             }
00444          }
00445 
00446       }  // if (i)
00447    } else { // if (hname)
00448       hname  = hdefault;
00449       hkeep  = 0;
00450       varexp = new char[strlen(varexp0)+1];
00451       strlcpy(varexp,varexp0,strlen(varexp0)+1);
00452       if (gDirectory) {
00453          fOldHistogram = (TH1*)gDirectory->Get(hname);
00454          if (fOldHistogram) { fOldHistogram->Delete(); fOldHistogram = 0;}
00455       }
00456    }
00457 
00458    // Decode varexp and selection
00459    if (!CompileVariables(varexp, realSelection.GetTitle())) {
00460       SetStatus(-1); 
00461       delete [] varexp;
00462       return;
00463    }
00464    if (fDimension > 4 && !(optpara || optcandle || opt5d)) {
00465       Error("Begin","Too many variables. Use the option \"para\", \"gl5d\" or \"candle\" to display more than 4 variables.");
00466       SetStatus(-1);
00467       delete [] varexp;
00468       return;
00469    }
00470 
00471    // In case fOldHistogram exists, check dimensionality
00472    Int_t nsel = strlen(selection);
00473    if (nsel > 1) {
00474       htitle.Form("%s {%s}",varexp,selection);
00475    } else {
00476       htitle = varexp;
00477    }
00478    if (fOldHistogram) {
00479       Int_t olddim = fOldHistogram->GetDimension();
00480       Int_t mustdelete = 0;
00481       if (fOldHistogram->InheritsFrom(TProfile::Class())) {
00482          profile = kTRUE;
00483          olddim = 2;
00484       }
00485       if (fOldHistogram->InheritsFrom(TProfile2D::Class())) {
00486          profile = kTRUE;
00487          olddim = 3;
00488       }
00489       if (opt.Contains("prof") && fDimension>1) {
00490          // ignore "prof" for 1D.
00491          if (!profile || olddim != fDimension) mustdelete = 1;
00492       } else {
00493          if (olddim != fDimension) mustdelete = 1;
00494       }
00495       if (mustdelete) {
00496          Warning("Begin","Deleting old histogram with different dimensions");
00497          delete fOldHistogram; fOldHistogram = 0;
00498       }
00499    }
00500 
00501    // Create a default canvas if none exists
00502    fDraw = 0;
00503    if (!gPad && !opt.Contains("goff") && fDimension > 0) {
00504       gROOT->MakeDefCanvas();
00505       if (!gPad)   {SetStatus(-1); return;}
00506    }
00507 
00508    // 1-D distribution
00509    TH1 *hist;
00510    if (fDimension == 1) {
00511       fAction = 1;
00512       if (!fOldHistogram) {
00513          fNbins[0] = gEnv->GetValue("Hist.Binning.1D.x",100);
00514          if (gPad && optSame) {
00515             TListIter np(gPad->GetListOfPrimitives());
00516             TObject *op;
00517             TH1 *oldhtemp = 0;
00518             while ((op = np()) && !oldhtemp) {
00519                if (op->InheritsFrom(TH1::Class())) oldhtemp = (TH1 *)op;
00520             }
00521             if (oldhtemp) {
00522                fNbins[0] = oldhtemp->GetXaxis()->GetNbins();
00523                fVmin[0]  = oldhtemp->GetXaxis()->GetXmin();
00524                fVmax[0]  = oldhtemp->GetXaxis()->GetXmax();
00525             } else {
00526                fVmin[0]  = gPad->GetUxmin();
00527                fVmax[0]  = gPad->GetUxmax();
00528             }
00529          } else {
00530             fAction   = -1;
00531             fVmin[0] = xmin;
00532             fVmax[0] = xmax;
00533             if (xmin < xmax) canRebin = kFALSE;
00534          }
00535       }
00536       if (fOldHistogram) {
00537          hist    = fOldHistogram;
00538          fNbins[0] = hist->GetXaxis()->GetNbins();
00539       } else {
00540          hist = new TH1F(hname,htitle.Data(),fNbins[0],fVmin[0],fVmax[0]);
00541          hist->SetLineColor(fTree->GetLineColor());
00542          hist->SetLineWidth(fTree->GetLineWidth());
00543          hist->SetLineStyle(fTree->GetLineStyle());
00544          hist->SetFillColor(fTree->GetFillColor());
00545          hist->SetFillStyle(fTree->GetFillStyle());
00546          hist->SetMarkerStyle(fTree->GetMarkerStyle());
00547          hist->SetMarkerColor(fTree->GetMarkerColor());
00548          hist->SetMarkerSize(fTree->GetMarkerSize());
00549          if (canRebin)hist->SetBit(TH1::kCanRebin);
00550          if (!hkeep) {
00551             hist->GetXaxis()->SetTitle(fVar[0]->GetTitle());
00552             hist->SetBit(kCanDelete);
00553             if (!opt.Contains("goff")) hist->SetDirectory(0);
00554          }
00555          if (opt.Length() && opt.Contains("e")) hist->Sumw2();
00556       }
00557       fVar[0]->SetAxis(hist->GetXaxis());
00558       fObject = hist;
00559 
00560       // 2-D distribution
00561    } else if (fDimension == 2 && !(optpara || optcandle)) {
00562       fAction = 2;
00563       if (!fOldHistogram || !optSame) {
00564          fNbins[0] = gEnv->GetValue("Hist.Binning.2D.y",40);
00565          fNbins[1] = gEnv->GetValue("Hist.Binning.2D.x",40);
00566          if (opt.Contains("prof")) fNbins[1] = gEnv->GetValue("Hist.Binning.2D.Prof",100);
00567          if (optSame) {
00568             TH1 *oldhtemp = (TH1*)gPad->FindObject(hdefault);
00569             if (oldhtemp) {
00570                fNbins[1] = oldhtemp->GetXaxis()->GetNbins();
00571                fVmin[1]  = oldhtemp->GetXaxis()->GetXmin();
00572                fVmax[1]  = oldhtemp->GetXaxis()->GetXmax();
00573                fNbins[0] = oldhtemp->GetYaxis()->GetNbins();
00574                fVmin[0]  = oldhtemp->GetYaxis()->GetXmin();
00575                fVmax[0]  = oldhtemp->GetYaxis()->GetXmax();
00576             } else {
00577                fNbins[1] = gEnv->GetValue("Hist.Binning.2D.x",40);
00578                fVmin[1]  = gPad->GetUxmin();
00579                fVmax[1]  = gPad->GetUxmax();
00580                fNbins[0] = gEnv->GetValue("Hist.Binning.2D.y",40);
00581                fVmin[0]  = gPad->GetUymin();
00582                fVmax[0]  = gPad->GetUymax();
00583             }
00584          } else {
00585             if (!fOldHistogram) fAction = -2;
00586             fVmin[1] = xmin;
00587             fVmax[1] = xmax;
00588             fVmin[0] = ymin;
00589             fVmax[0] = ymax;
00590             if (xmin < xmax && ymin < ymax) canRebin = kFALSE;
00591          }
00592       }
00593       if (profile || opt.Contains("prof")) {
00594          TProfile *hp;
00595          if (fOldHistogram) {
00596             fAction = 4;
00597             hp = (TProfile*)fOldHistogram;
00598          } else {
00599             if (fAction < 0) {
00600                fAction = -4;
00601                fVmin[1] = xmin;
00602                fVmax[1] = xmax;
00603                if (xmin < xmax) canRebin = kFALSE;
00604             }
00605             if (fAction == 2) {
00606                //we come here when option = "same prof"
00607                fAction = -4;
00608                TH1 *oldhtemp = (TH1*)gPad->FindObject(hdefault);
00609                if (oldhtemp) {
00610                   fNbins[1] = oldhtemp->GetXaxis()->GetNbins();
00611                   fVmin[1]  = oldhtemp->GetXaxis()->GetXmin();
00612                   fVmax[1]  = oldhtemp->GetXaxis()->GetXmax();
00613                }
00614             }
00615             if (opt.Contains("profs")) {
00616                hp = new TProfile(hname,htitle.Data(),fNbins[1],fVmin[1], fVmax[1],"s");
00617             } else if (opt.Contains("profi")) {
00618                hp = new TProfile(hname,htitle.Data(),fNbins[1],fVmin[1], fVmax[1],"i");
00619             } else if (opt.Contains("profg")) {
00620                hp = new TProfile(hname,htitle.Data(),fNbins[1],fVmin[1], fVmax[1],"g");
00621             } else {
00622                hp = new TProfile(hname,htitle.Data(),fNbins[1],fVmin[1], fVmax[1],"");
00623             }
00624             if (!hkeep) {
00625                hp->SetBit(kCanDelete);
00626                if (!opt.Contains("goff")) hp->SetDirectory(0);
00627             }
00628             hp->SetLineColor(fTree->GetLineColor());
00629             hp->SetLineWidth(fTree->GetLineWidth());
00630             hp->SetLineStyle(fTree->GetLineStyle());
00631             hp->SetFillColor(fTree->GetFillColor());
00632             hp->SetFillStyle(fTree->GetFillStyle());
00633             hp->SetMarkerStyle(fTree->GetMarkerStyle());
00634             hp->SetMarkerColor(fTree->GetMarkerColor());
00635             hp->SetMarkerSize(fTree->GetMarkerSize());
00636             if (canRebin)hp->SetBit(TH1::kCanRebin);
00637          }
00638          fVar[1]->SetAxis(hp->GetXaxis());
00639          fObject = hp;
00640 
00641       } else {
00642          TH2F *h2;
00643          if (fOldHistogram) {
00644             h2 = (TH2F*)fOldHistogram;
00645          } else {
00646             h2 = new TH2F(hname,htitle.Data(),fNbins[1],fVmin[1], fVmax[1], fNbins[0], fVmin[0], fVmax[0]);
00647             h2->SetLineColor(fTree->GetLineColor());
00648             h2->SetFillColor(fTree->GetFillColor());
00649             h2->SetFillStyle(fTree->GetFillStyle());
00650             h2->SetMarkerStyle(fTree->GetMarkerStyle());
00651             h2->SetMarkerColor(fTree->GetMarkerColor());
00652             h2->SetMarkerSize(fTree->GetMarkerSize());
00653             if (canRebin)h2->SetBit(TH1::kCanRebin);
00654             if (!hkeep) {
00655                h2->GetXaxis()->SetTitle(fVar[1]->GetTitle());
00656                h2->GetYaxis()->SetTitle(fVar[0]->GetTitle());
00657                h2->SetBit(TH1::kNoStats);
00658                h2->SetBit(kCanDelete);
00659                if (!opt.Contains("goff")) h2->SetDirectory(0);
00660             }
00661          }
00662          fVar[0]->SetAxis(h2->GetYaxis());
00663          fVar[1]->SetAxis(h2->GetXaxis());
00664          Bool_t graph = kFALSE;
00665          Int_t l = opt.Length();
00666          if (l == 0 || optSame) graph = kTRUE;
00667          if (opt.Contains("p")     || opt.Contains("*")    || opt.Contains("l"))    graph = kTRUE;
00668          if (opt.Contains("surf")  || opt.Contains("lego") || opt.Contains("cont")) graph = kFALSE;
00669          if (opt.Contains("col")   || opt.Contains("hist") || opt.Contains("scat")) graph = kFALSE;
00670          if (opt.Contains("box"))                                                   graph = kFALSE;
00671          fObject = h2;
00672          if (graph) {
00673             fAction = 12;
00674             if (!fOldHistogram && !optSame) fAction = -12;
00675          }
00676       }
00677 
00678       // 3-D distribution
00679    } else if ((fDimension == 3 || fDimension == 4) && !(optpara || optcandle)) {
00680       fAction = 3;
00681       if (fDimension == 4) fAction = 40;
00682       if (!fOldHistogram || !optSame) {
00683          fNbins[0] = gEnv->GetValue("Hist.Binning.3D.z",20);
00684          fNbins[1] = gEnv->GetValue("Hist.Binning.3D.y",20);
00685          fNbins[2] = gEnv->GetValue("Hist.Binning.3D.x",20);
00686          if (fDimension == 3 && opt.Contains("prof")) {
00687             fNbins[1] = gEnv->GetValue("Hist.Binning.3D.Profy",20);
00688             fNbins[2] = gEnv->GetValue("Hist.Binning.3D.Profx",20);
00689          }
00690          if (optSame) {
00691             TH1 *oldhtemp = (TH1*)gPad->FindObject(hdefault);
00692             if (oldhtemp) {
00693                fNbins[2] = oldhtemp->GetXaxis()->GetNbins();
00694                fVmin[2]  = oldhtemp->GetXaxis()->GetXmin();
00695                fVmax[2]  = oldhtemp->GetXaxis()->GetXmax();
00696                fNbins[1] = oldhtemp->GetYaxis()->GetNbins();
00697                fVmin[1]  = oldhtemp->GetYaxis()->GetXmin();
00698                fVmax[1]  = oldhtemp->GetYaxis()->GetXmax();
00699                fNbins[0] = oldhtemp->GetZaxis()->GetNbins();
00700                fVmin[0]  = oldhtemp->GetZaxis()->GetXmin();
00701                fVmax[0]  = oldhtemp->GetZaxis()->GetXmax();
00702             } else {
00703                TView *view = gPad->GetView();
00704                if (!view) {
00705                   Error("Begin","You cannot use option same when no 3D view exists");
00706                   fVmin[0]=fVmin[1]=fVmin[2]=-1;
00707                   fVmax[0]=fVmax[1]=fVmax[2]= 1;
00708                   view = TView::CreateView(1,fVmin,fVmax);
00709                }
00710                Double_t *rmin = view->GetRmin();
00711                Double_t *rmax = view->GetRmax();
00712                fNbins[2] = gEnv->GetValue("Hist.Binning.3D.z",20);
00713                fVmin[2]  = rmin[0];
00714                fVmax[2]  = rmax[0];
00715                fNbins[1] = gEnv->GetValue("Hist.Binning.3D.y",20);
00716                fVmin[1]  = rmin[1];
00717                fVmax[1]  = rmax[1];
00718                fNbins[0] = gEnv->GetValue("Hist.Binning.3D.x",20);
00719                fVmin[0]  = rmin[2];
00720                fVmax[0]  = rmax[2];
00721             }
00722          } else {
00723             if (!fOldHistogram && fDimension ==3) fAction = -3;
00724             fVmin[2] = xmin;
00725             fVmax[2] = xmax;
00726             fVmin[1] = ymin;
00727             fVmax[1] = ymax;
00728             fVmin[0] = zmin;
00729             fVmax[0] = zmax;
00730             if (xmin < xmax && ymin < ymax && zmin < zmax) canRebin = kFALSE;
00731          }
00732       }
00733       if ((fDimension == 3) && (profile || opt.Contains("prof"))) {
00734          TProfile2D *hp;
00735          if (fOldHistogram) {
00736             fAction = 23;
00737             hp = (TProfile2D*)fOldHistogram;
00738          } else {
00739             if (fAction < 0) {
00740                fAction = -23;
00741                fVmin[2] = xmin;
00742                fVmax[2] = xmax;
00743                fVmin[1] = ymin;
00744                fVmax[1] = ymax;
00745                if (xmin < xmax && ymin < ymax) canRebin = kFALSE;
00746             }
00747             if (opt.Contains("profs")) {
00748                hp = new TProfile2D(hname,htitle.Data(),fNbins[2],fVmin[2], fVmax[2],fNbins[1],fVmin[1], fVmax[1],"s");
00749             } else if (opt.Contains("profi")) {
00750                hp = new TProfile2D(hname,htitle.Data(),fNbins[2],fVmin[2], fVmax[2],fNbins[1],fVmin[1], fVmax[1],"i");
00751             } else if (opt.Contains("profg")) {
00752                hp = new TProfile2D(hname,htitle.Data(),fNbins[2],fVmin[2], fVmax[2],fNbins[1],fVmin[1], fVmax[1],"g");
00753             } else {
00754                hp = new TProfile2D(hname,htitle.Data(),fNbins[2],fVmin[2], fVmax[2],fNbins[1],fVmin[1], fVmax[1],"");
00755             }
00756             if (!hkeep) {
00757                hp->SetBit(kCanDelete);
00758                if (!opt.Contains("goff")) hp->SetDirectory(0);
00759             }
00760             hp->SetLineColor(fTree->GetLineColor());
00761             hp->SetLineWidth(fTree->GetLineWidth());
00762             hp->SetLineStyle(fTree->GetLineStyle());
00763             hp->SetFillColor(fTree->GetFillColor());
00764             hp->SetFillStyle(fTree->GetFillStyle());
00765             hp->SetMarkerStyle(fTree->GetMarkerStyle());
00766             hp->SetMarkerColor(fTree->GetMarkerColor());
00767             hp->SetMarkerSize(fTree->GetMarkerSize());
00768             if (canRebin)hp->SetBit(TH1::kCanRebin);
00769          }
00770          fVar[1]->SetAxis(hp->GetYaxis());
00771          fVar[2]->SetAxis(hp->GetXaxis());
00772          fObject = hp;
00773       } else if (fDimension == 3 && opt.Contains("col")) {
00774          TH2F *h2;
00775          if (fOldHistogram) {
00776             h2 = (TH2F*)fOldHistogram;
00777          } else {
00778             h2 = new TH2F(hname,htitle.Data(),fNbins[1],fVmin[1], fVmax[1], fNbins[0], fVmin[0], fVmax[0]);
00779             h2->SetLineColor(fTree->GetLineColor());
00780             h2->SetFillColor(fTree->GetFillColor());
00781             h2->SetFillStyle(fTree->GetFillStyle());
00782             h2->SetMarkerStyle(fTree->GetMarkerStyle());
00783             h2->SetMarkerColor(fTree->GetMarkerColor());
00784             h2->SetMarkerSize(fTree->GetMarkerSize());
00785             if (canRebin)h2->SetBit(TH1::kCanRebin);
00786             if (!hkeep) {
00787                h2->GetXaxis()->SetTitle(fVar[1]->GetTitle());
00788                h2->GetZaxis()->SetTitle(fVar[0]->GetTitle());
00789                h2->SetBit(TH1::kNoStats);
00790                h2->SetBit(kCanDelete);
00791                if (!opt.Contains("goff")) h2->SetDirectory(0);
00792             }
00793          }
00794          fVar[0]->SetAxis(h2->GetYaxis());
00795          fVar[1]->SetAxis(h2->GetXaxis());
00796          fObject = h2;
00797          fAction = 33;
00798       } else {
00799          TH3F *h3;
00800          if (fOldHistogram) {
00801             h3 = (TH3F*)fOldHistogram;
00802          } else {
00803             h3 = new TH3F(hname,htitle.Data(),fNbins[2],fVmin[2], fVmax[2],fNbins[1],fVmin[1], fVmax[1], fNbins[0], fVmin[0], fVmax[0]);
00804             h3->SetLineColor(fTree->GetLineColor());
00805             h3->SetFillColor(fTree->GetFillColor());
00806             h3->SetFillStyle(fTree->GetFillStyle());
00807             h3->SetMarkerStyle(fTree->GetMarkerStyle());
00808             h3->SetMarkerColor(fTree->GetMarkerColor());
00809             h3->SetMarkerSize(fTree->GetMarkerSize());
00810             if (canRebin)h3->SetBit(TH1::kCanRebin);
00811             if (!hkeep) {
00812                //small correction for the title offsets in x,y to take into account the angles
00813                Double_t xoffset = h3->GetXaxis()->GetTitleOffset();
00814                Double_t yoffset = h3->GetYaxis()->GetTitleOffset();
00815                h3->GetXaxis()->SetTitleOffset(1.2*xoffset);
00816                h3->GetYaxis()->SetTitleOffset(1.2*yoffset);
00817                h3->GetXaxis()->SetTitle(fVar[2]->GetTitle());
00818                h3->GetYaxis()->SetTitle(fVar[1]->GetTitle());
00819                h3->GetZaxis()->SetTitle(fVar[0]->GetTitle());
00820                h3->SetBit(kCanDelete);
00821                h3->SetBit(TH1::kNoStats);
00822                if (!opt.Contains("goff")) h3->SetDirectory(0);
00823             }
00824          }
00825          fVar[0]->SetAxis(h3->GetZaxis());
00826          fVar[1]->SetAxis(h3->GetYaxis());
00827          fVar[2]->SetAxis(h3->GetXaxis());
00828          fObject = h3;
00829          Int_t noscat = strlen(option);
00830          if (optSame) noscat -= 4;
00831          if (!noscat && fDimension ==3) {
00832             fAction = 13;
00833             if (!fOldHistogram && !optSame) fAction = -13;
00834          }
00835       }
00836       // An Event List
00837    } else if (enlist) {
00838       fAction = 5;
00839       fOldEstimate = fTree->GetEstimate();
00840       fTree->SetEstimate(1);
00841       fObject = enlist;
00842    } else if (evlist) {
00843       fAction = 5;
00844       fOldEstimate = fTree->GetEstimate();
00845       fTree->SetEstimate(1);
00846       fObject = evlist;
00847    } else if (optcandle || optpara || opt5d) {
00848       if (optcandle)  fAction = 7;
00849       else if (opt5d) fAction = 8;
00850       else            fAction = 6;
00851    }
00852    if (hkeep) delete [] varexp;
00853    if (hnamealloc) delete [] hnamealloc;
00854    for(i=0;i<fValSize;++i)
00855       fVarMultiple[i] = kFALSE;
00856    fSelectMultiple = kFALSE;
00857    for(i=0;i<fDimension;++i){
00858       if(fVar[i] && fVar[i]->GetMultiplicity()) fVarMultiple[i] = kTRUE;
00859    }
00860 
00861    if (fSelect && fSelect->GetMultiplicity()) fSelectMultiple = kTRUE;
00862 
00863    fForceRead = fTree->TestBit(TTree::kForceRead);
00864    fWeight  = fTree->GetWeight();
00865    fNfill   = 0;
00866 
00867    for(i=0;i<fDimension;++i){
00868       if(!fVal[i] && fVar[i]){
00869          fVal[i] = new Double_t[(Int_t)fTree->GetEstimate()];
00870       }
00871    }
00872 
00873    if (!fW)             fW  = new Double_t[(Int_t)fTree->GetEstimate()];
00874 
00875    for(i=0;i<fValSize;++i){
00876       fVmin[i] = FLT_MAX;
00877       fVmax[i] = -FLT_MAX;
00878    }
00879 }
00880 
00881 //______________________________________________________________________________
00882 void TSelectorDraw::ClearFormula()
00883 {
00884    // Delete internal buffers.
00885 
00886    ResetBit(kWarn);
00887    for (Int_t i=0;i<fValSize;++i){
00888       delete fVar[i];
00889       fVar[i] = 0;
00890    }
00891    delete fSelect; fSelect = 0;
00892    fManager = 0;
00893    fMultiplicity = 0;
00894 }
00895 
00896 //______________________________________________________________________________
00897 Bool_t TSelectorDraw::CompileVariables(const char *varexp, const char *selection)
00898 {
00899    // Compile input variables and selection expression.
00900    //
00901    //  varexp is an expression of the general form e1:e2:e3
00902    //    where e1,etc is a formula referencing a combination of the columns
00903    //  Example:
00904    //     varexp = x  simplest case: draw a 1-Dim distribution of column named x
00905    //            = sqrt(x)         : draw distribution of sqrt(x)
00906    //            = x*y/z
00907    //            = y:sqrt(x) 2-Dim dsitribution of y versus sqrt(x)
00908    //
00909    //  selection is an expression with a combination of the columns
00910    //  Example:
00911    //      selection = "x<y && sqrt(z)>3.2"
00912    //       in a selection all the C++ operators are authorized
00913    //
00914    //  Return kFALSE if any of the variable is not compilable.
00915 
00916    Int_t i,nch,ncols;
00917 
00918    // Compile selection expression if there is one
00919    fDimension = 0;
00920    ClearFormula();
00921    fMultiplicity = 0;
00922    fObjEval = kFALSE;
00923 
00924    if (strlen(selection)) {
00925       fSelect = new TTreeFormula("Selection",selection,fTree);
00926       fSelect->SetQuickLoad(kTRUE);
00927       if (!fSelect->GetNdim()) {delete fSelect; fSelect = 0; return kFALSE; }
00928    }
00929 
00930    // if varexp is empty, take first column by default
00931    nch = strlen(varexp);
00932    if (nch == 0) {
00933       fDimension = 0;
00934       fManager = new TTreeFormulaManager();
00935       if (fSelect) fManager->Add(fSelect);
00936       fTree->ResetBit(TTree::kForceRead);
00937 
00938       fManager->Sync();
00939 
00940       if (fManager->GetMultiplicity()==-1) fTree->SetBit(TTree::kForceRead);
00941       if (fManager->GetMultiplicity()>=1) fMultiplicity = fManager->GetMultiplicity();
00942 
00943       return kTRUE;
00944    }
00945 
00946    // otherwise select only the specified columns
00947    std::vector<TString> varnames;
00948    ncols  = SplitNames(varexp,varnames);
00949 
00950    InitArrays(ncols);
00951 
00952    fManager = new TTreeFormulaManager();
00953    if (fSelect) fManager->Add(fSelect);
00954    fTree->ResetBit(TTree::kForceRead);
00955    for(i=0; i<ncols;++i){
00956       fVar[i] = new TTreeFormula(Form("Var%i",i+1),varnames[i].Data(),fTree);
00957       fVar[i]->SetQuickLoad(kTRUE);
00958       if(!fVar[i]->GetNdim()) { ClearFormula(); return kFALSE; }
00959       fManager->Add(fVar[i]);
00960    }
00961    fManager->Sync();
00962 
00963    if (fManager->GetMultiplicity()==-1) fTree->SetBit(TTree::kForceRead);
00964    if (fManager->GetMultiplicity()>=1) fMultiplicity = fManager->GetMultiplicity();
00965 
00966    fDimension    = ncols;
00967 
00968    if (ncols==1) {
00969       TClass *cl = fVar[0]->EvalClass();
00970       if (cl) {
00971          fObjEval = kTRUE;
00972       }
00973    }
00974    return kTRUE;
00975 }
00976 
00977 //______________________________________________________________________________
00978 Double_t* TSelectorDraw::GetVal(Int_t i) const
00979 {
00980    //Get variable buffer.
00981 
00982    if(i<0 || i >= fDimension)
00983       return 0;
00984    else
00985       return fVal[i];
00986 }
00987 
00988 //______________________________________________________________________________
00989 TTreeFormula* TSelectorDraw::GetVar(Int_t i) const
00990 {
00991    //Get variable formula.
00992 
00993    if(i<0 || i>=fDimension)
00994       return 0;
00995    else
00996       return fVar[i];
00997 }
00998 
00999 //______________________________________________________________________________
01000 void TSelectorDraw::InitArrays(Int_t newsize)
01001 {
01002    // Initialization of the primitive type arrays if the new size is bigger than the available space.
01003 
01004    if(newsize>fValSize){
01005       Int_t oldsize = fValSize;
01006       while(fValSize<newsize)
01007          fValSize*=2;            // Double the available space until it matches the new size.
01008       if(fNbins) delete [] fNbins;
01009       if(fVmin) delete [] fVmin;
01010       if(fVmax) delete [] fVmax;
01011       if(fVarMultiple) delete [] fVarMultiple;
01012 
01013       fNbins = new Int_t[fValSize];
01014       fVmin = new Double_t[fValSize];
01015       fVmax = new Double_t[fValSize];
01016       fVarMultiple = new Bool_t[fValSize];
01017 
01018       for(Int_t i=0;i<oldsize;++i)
01019          delete [] fVal[i];
01020       delete [] fVal;
01021       delete [] fVar;
01022       fVal = new Double_t*[fValSize];
01023       fVar = new TTreeFormula*[fValSize];
01024       for(Int_t i=0;i<fValSize;++i){
01025          fVal[i] = 0;
01026          fVar[i] = 0;
01027       }
01028    }
01029 }
01030 
01031 //______________________________________________________________________________
01032 UInt_t TSelectorDraw::SplitNames(const TString &varexp, std::vector<TString> &names)
01033 {
01034    // Build Index array for names in varexp.
01035    // This will allocated a C style array of TString and Ints
01036 
01037    names.clear();
01038 
01039    Bool_t ternary = kFALSE;
01040    Int_t prev = 0;
01041    for (Int_t i=0;i<varexp.Length();i++) {
01042       if (varexp[i] == ':'
01043           && ! ( (i>0&&varexp[i-1]==':') || varexp[i+1]==':' )
01044           ) {
01045          if (ternary) {
01046             ternary = kFALSE;
01047          } else {
01048             names.push_back( varexp(prev,i-prev) );
01049             prev = i+1;
01050          }
01051       }
01052       if (varexp[i] == '?') {
01053          ternary = kTRUE;
01054       }
01055    }
01056    names.push_back( varexp(prev, varexp.Length()-prev) );
01057    return names.size();
01058 }
01059 
01060 
01061 //______________________________________________________________________________
01062 Bool_t TSelectorDraw::Notify()
01063 {
01064    // This function is called at the first entry of a new tree in a chain.
01065 
01066    if (fTree) fWeight  = fTree->GetWeight();
01067    if(fVar){
01068       for(Int_t i=0;i<fDimension;++i){
01069          if(fVar[i]) fVar[i]->UpdateFormulaLeaves();
01070       }
01071    }
01072    if (fSelect) fSelect->UpdateFormulaLeaves();
01073    return kTRUE;
01074 }
01075 
01076 //______________________________________________________________________________
01077 void TSelectorDraw::ProcessFill(Long64_t entry)
01078 {
01079    // Called in the entry loop for all entries accepted by Select.
01080 
01081    if (fObjEval) {
01082       ProcessFillObject(entry);
01083       return;
01084    }
01085 
01086    if (fMultiplicity) {
01087       ProcessFillMultiple(entry);
01088       return;
01089    }
01090 
01091    // simple case with no multiplicity
01092    if ( fForceRead && fManager->GetNdata() <= 0) return;
01093 
01094    if (fSelect) {
01095       fW[fNfill] = fWeight*fSelect->EvalInstance(0);
01096       if (!fW[fNfill]) return;
01097    } else fW[fNfill] = fWeight;
01098    if(fVal){
01099       for(Int_t i=0;i<fDimension;++i){
01100          if(fVar[i]) fVal[i][fNfill] = fVar[i]->EvalInstance(0);
01101       }
01102    }
01103    fNfill++;
01104    if (fNfill >= fTree->GetEstimate()) {
01105       TakeAction();
01106       fNfill = 0;
01107    }
01108 }
01109 
01110 //______________________________________________________________________________
01111 void TSelectorDraw::ProcessFillMultiple(Long64_t /*entry*/)
01112 {
01113    // Called in the entry loop for all entries accepted by Select.
01114    // Complex case with multiplicity.
01115 
01116    // Grab the array size of the formulas for this entry
01117    Int_t ndata = fManager->GetNdata();
01118 
01119    // No data at all, let's move on to the next entry.
01120    if (!ndata) return;
01121 
01122    Int_t nfill0 = fNfill;
01123 
01124    // Calculate the first values
01125    if (fSelect) {
01126       fW[fNfill] = fWeight*fSelect->EvalInstance(0);
01127       if (!fW[fNfill] && !fSelectMultiple) return;
01128    } else fW[fNfill] = fWeight;
01129 
01130    // Always call EvalInstance(0) to insure the loading
01131    // of the branches.
01132    if (fW[fNfill]) {
01133       for(Int_t i=0;i<fDimension;++i){
01134          if(fVar[i]) fVal[i][fNfill] = fVar[i]->EvalInstance(0);
01135       }
01136       fNfill++;
01137       if (fNfill >= fTree->GetEstimate()) {
01138          TakeAction();
01139          fNfill = 0;
01140       }
01141    } else {
01142       for(Int_t i=0;i<fDimension;++i){
01143          if(fVar[i]) fVar[i]->ResetLoading();
01144       }
01145    }
01146    Double_t ww = fW[nfill0];
01147 
01148    for (Int_t i=1;i<ndata;i++) {
01149       if (fSelectMultiple) {
01150          ww = fWeight*fSelect->EvalInstance(i);
01151          if (ww == 0) continue;
01152          if (fNfill == nfill0) {
01153             for(Int_t k=0;k<fDimension;++k){
01154                if(!fVarMultiple[k]) fVal[k][fNfill] = fVar[k]->EvalInstance(0);
01155             }
01156          }
01157       }
01158       for(Int_t k=0;k<fDimension;++k){
01159          if(fVarMultiple[k]) fVal[k][fNfill] = fVar[k]->EvalInstance(i);
01160          else                fVal[k][fNfill] = fVal[k][nfill0];
01161       }
01162       fW[fNfill] = ww;
01163 
01164       fNfill++;
01165       if (fNfill >= fTree->GetEstimate()) {
01166          TakeAction();
01167          fNfill = 0;
01168       }
01169    }
01170 }
01171 
01172 //______________________________________________________________________________
01173 void TSelectorDraw::ProcessFillObject(Long64_t /*entry*/)
01174 {
01175    // Called in the entry loop for all entries accepted by Select.
01176    // Case where the only variable returns an object (or pointer to).
01177 
01178    // Complex case with multiplicity.
01179 
01180    // Grab the array size of the formulas for this entry
01181    Int_t ndata = fManager->GetNdata();
01182 
01183    // No data at all, let's move on to the next entry.
01184    if (!ndata) return;
01185 
01186    Int_t nfill0 = fNfill;
01187    Double_t ww = 0;
01188 
01189    for (Int_t i=0;i<ndata;i++) {
01190       if (i==0) {
01191          if (fSelect) {
01192             fW[fNfill] = fWeight*fSelect->EvalInstance(0);
01193             if (!fW[fNfill] && !fSelectMultiple) return;
01194          } else fW[fNfill] = fWeight;
01195          ww = fW[nfill0];
01196       } else if (fSelectMultiple) {
01197          ww = fWeight*fSelect->EvalInstance(i);
01198          if (ww == 0) continue;
01199       }
01200       if (fDimension >= 1 && fVar[0]) {
01201          TClass *cl = fVar[0]->EvalClass();
01202          if (cl==TBits::Class()) {
01203 
01204             void *obj = fVar[0]->EvalObject(i);
01205 
01206             TBits *bits = (TBits*)obj;
01207             Int_t nbits = bits->GetNbits();
01208 
01209             Int_t nextbit = -1;
01210             while(1) {
01211                nextbit = bits->FirstSetBit(nextbit+1);
01212                if (nextbit >= nbits) break;
01213                fVal[0][fNfill] = nextbit;
01214                fW[fNfill] =  ww;
01215                fNfill++;
01216             }
01217 
01218          } else {
01219 
01220             if (!TestBit(kWarn)) {
01221                Warning("ProcessFillObject",
01222                        "Not implemented for %s",
01223                        cl?cl->GetName():"unknown class");
01224                SetBit(kWarn);
01225             }
01226 
01227          }
01228       }
01229       if (fNfill >= fTree->GetEstimate()) {
01230          TakeAction();
01231          fNfill = 0;
01232       }
01233    }
01234 
01235 }
01236 
01237 //_______________________________________________________________________
01238 void TSelectorDraw::SetEstimate(Long64_t)
01239 {
01240    // Set number of entries to estimate variable limits.
01241 
01242    if(fVal){
01243       for(Int_t i=0;i<fValSize;++i){
01244          delete [] fVal[i];
01245          fVal[i] = 0;
01246       }
01247    }
01248    delete [] fW;   fW  = 0;
01249 }
01250 
01251 //______________________________________________________________________________
01252 void TSelectorDraw::TakeAction()
01253 {
01254    // Execute action for object obj fNfill times.
01255 
01256    Int_t i;
01257    //__________________________1-D histogram_______________________
01258    if      (fAction ==  1) ((TH1*)fObject)->FillN(fNfill,fVal[0],fW);
01259    //__________________________2-D histogram_______________________
01260    else if (fAction ==  2) {
01261       TH2 *h2 = (TH2*)fObject;
01262       for(i=0;i<fNfill;i++) h2->Fill(fVal[1][i],fVal[0][i],fW[i]);
01263    }
01264    //__________________________Profile histogram_______________________
01265    else if (fAction ==  4) ((TProfile*)fObject)->FillN(fNfill,fVal[1],fVal[0],fW);
01266    //__________________________Event List______________________________
01267    else if (fAction ==  5) {
01268       if (fObject->InheritsFrom(TEntryList::Class())){
01269          TEntryList *enlist = (TEntryList*)fObject;
01270          Long64_t enumb = fTree->GetTree()->GetReadEntry();
01271          enlist->Enter(enumb);
01272       }
01273       else {
01274          TEventList *evlist = (TEventList*)fObject;
01275          Long64_t enumb = fTree->GetChainOffset() + fTree->GetTree()->GetReadEntry();
01276          if (evlist->GetIndex(enumb) < 0) evlist->Enter(enumb);
01277       }
01278    }
01279    //__________________________2D scatter plot_______________________
01280    else if (fAction == 12) {
01281       TH2 *h2 = (TH2*)fObject;
01282       if (h2->TestBit(TH1::kCanRebin) && h2->TestBit(kCanDelete)) {
01283          for (i=0;i<fNfill;i++) {
01284             if (fVmin[0] > fVal[0][i]) fVmin[0] = fVal[0][i];
01285             if (fVmax[0] < fVal[0][i]) fVmax[0] = fVal[0][i];
01286             if (fVmin[1] > fVal[1][i]) fVmin[1] = fVal[1][i];
01287             if (fVmax[1] < fVal[1][i]) fVmax[1] = fVal[1][i];
01288          }
01289          THLimitsFinder::GetLimitsFinder()->FindGoodLimits(h2,fVmin[1],fVmax[1],fVmin[0],fVmax[0]);
01290       }
01291       TGraph *pm = new TGraph(fNfill,fVal[1],fVal[0]);
01292       pm->SetEditable(kFALSE);
01293       pm->SetBit(kCanDelete);
01294       pm->SetMarkerStyle(fTree->GetMarkerStyle());
01295       pm->SetMarkerColor(fTree->GetMarkerColor());
01296       pm->SetMarkerSize(fTree->GetMarkerSize());
01297       pm->SetLineColor(fTree->GetLineColor());
01298       pm->SetLineStyle(fTree->GetLineStyle());
01299       pm->SetFillColor(fTree->GetFillColor());
01300       pm->SetFillStyle(fTree->GetFillStyle());
01301 
01302       if (!fDraw && !strstr(fOption.Data(),"goff")) {
01303          if (fOption.Length() == 0 || strcasecmp(fOption.Data(),"same")==0)  pm->Draw("p");
01304          else                                                                pm->Draw(fOption.Data());
01305       }
01306       if (!h2->TestBit(kCanDelete)) {
01307          for (i=0;i<fNfill;i++) h2->Fill(fVal[1][i],fVal[0][i],fW[i]);
01308       }
01309    }
01310    //__________________________3D scatter plot_______________________
01311    else if (fAction ==  3) {
01312       TH3 *h3 =(TH3*)fObject;
01313       for(i=0;i<fNfill;i++) h3->Fill(fVal[2][i],fVal[1][i],fVal[0][i],fW[i]);
01314    }
01315    else if (fAction == 13) {
01316       TPolyMarker3D *pm3d = new TPolyMarker3D(fNfill);
01317       pm3d->SetMarkerStyle(fTree->GetMarkerStyle());
01318       pm3d->SetMarkerColor(fTree->GetMarkerColor());
01319       pm3d->SetMarkerSize(fTree->GetMarkerSize());
01320       for (i=0;i<fNfill;i++) { pm3d->SetPoint(i,fVal[2][i],fVal[1][i],fVal[0][i]);}
01321       pm3d->Draw();
01322       TH3 *h3 =(TH3*)fObject;
01323       for(i=0;i<fNfill;i++) h3->Fill(fVal[2][i],fVal[1][i],fVal[0][i],fW[i]);
01324    }
01325    //__________________________3D scatter plot (3rd variable = col)__
01326    else if (fAction == 33) {
01327       TH2 *h2 = (TH2*)fObject;
01328       TakeEstimate();
01329       Int_t ncolors  = gStyle->GetNumberOfColors();
01330       TObjArray *grs = (TObjArray*)h2->GetListOfFunctions()->FindObject("graphs");
01331       Int_t col;
01332       TGraph *gr;
01333       if (!grs) {
01334          grs = new TObjArray(ncolors);
01335          grs->SetOwner();
01336          grs->SetName("graphs");
01337          h2->GetListOfFunctions()->Add(grs, "P");
01338          for (col=0;col<ncolors;col++) {
01339             gr = new TGraph();
01340             gr->SetMarkerColor(gStyle->GetColorPalette(col));
01341             gr->SetMarkerStyle(fTree->GetMarkerStyle());
01342             gr->SetMarkerSize(fTree->GetMarkerSize());
01343             grs->AddAt(gr,col);
01344          }
01345       }
01346       h2->SetEntries(fNfill);
01347       h2->SetMinimum(fVmin[2]);
01348       h2->SetMaximum(fVmax[2]);
01349       // Fill the graphs acording to the color
01350       for (i=0;i<fNfill;i++) {
01351          col = Int_t(ncolors*((fVal[2][i]-fVmin[2])/(fVmax[2]-fVmin[2])));
01352          if (col < 0) col = 0;
01353          if (col > ncolors-1) col = ncolors-1;
01354          gr = (TGraph*)grs->UncheckedAt(col);
01355          if (gr) gr->SetPoint(gr->GetN(),fVal[1][i],fVal[0][i]);
01356       }
01357       // Remove potential empty graphs
01358       for (col=0;col<ncolors;col++) {
01359          gr = (TGraph*)grs->At(col);
01360          if (gr && gr->GetN() <= 0) grs->Remove(gr);
01361       }
01362    }
01363    //__________________________2D Profile Histogram__________________
01364    else if (fAction == 23) {
01365       TProfile2D *hp2 =(TProfile2D*)fObject;
01366       for(i=0;i<fNfill;i++) hp2->Fill(fVal[2][i],fVal[1][i],fVal[0][i],fW[i]);
01367    }
01368    //__________________________4D scatter plot_______________________
01369    else if (fAction ==  40) {
01370       TakeEstimate();
01371       TH3 *h3 =(TH3*)fObject;
01372       Int_t ncolors  = gStyle->GetNumberOfColors();
01373       if (ncolors == 0) {
01374          TColor::InitializeColors();
01375          ncolors  = gStyle->GetNumberOfColors();
01376       }
01377       TObjArray *pms = (TObjArray*)h3->GetListOfFunctions()->FindObject("polymarkers");
01378       Int_t col;
01379       TPolyMarker3D *pm3d;
01380       if (!pms) {
01381          pms = new TObjArray(ncolors);
01382          pms->SetOwner();
01383          pms->SetName("polymarkers");
01384          h3->GetListOfFunctions()->Add(pms);
01385          for (col=0;col<ncolors;col++) {
01386             pm3d = new TPolyMarker3D();
01387             pm3d->SetMarkerColor(gStyle->GetColorPalette(col));
01388             pm3d->SetMarkerStyle(fTree->GetMarkerStyle());
01389             pm3d->SetMarkerSize(fTree->GetMarkerSize());
01390             pms->AddAt(pm3d,col);
01391          }
01392       }
01393       h3->SetEntries(fNfill);
01394       h3->SetMinimum(fVmin[3]);
01395       h3->SetMaximum(fVmax[3]);
01396       for (i=0;i<fNfill;i++) {
01397          col = Int_t(ncolors*((fVal[3][i]-fVmin[3])/(fVmax[3]-fVmin[3])));
01398          if (col > ncolors-1) col = ncolors-1;
01399          if (col < 0) col = 0;
01400          pm3d = (TPolyMarker3D*)pms->UncheckedAt(col);
01401          pm3d->SetPoint(pm3d->GetLastPoint()+1,fVal[2][i],fVal[1][i],fVal[0][i]);
01402       }
01403    }
01404    //__________________________Parallel coordinates / candle chart_______________________
01405    else if (fAction == 6 || fAction == 7) {
01406       TakeEstimate();
01407       Bool_t candle = (fAction==7);
01408       // Using CINT to avoid a dependency in TParallelCoord
01409       if (!fOption.Contains("goff"))
01410          gROOT->ProcessLineFast(Form("TParallelCoord::BuildParallelCoord((TSelectorDraw*)0x%lx,0x%lx",
01411                                 (ULong_t)this, (ULong_t)candle));
01412    } else if (fAction == 8) {
01413       //gROOT->ProcessLineFast(Form("(new TGL5DDataSet((TTree *)0x%1x))->Draw(\"%s\");", fTree, fOption.Data()));
01414    }
01415    //__________________________something else_______________________
01416    else if (fAction < 0) {
01417       fAction = -fAction;
01418       TakeEstimate();
01419    }
01420 
01421    // Do we need to update screen?
01422    fSelectedRows += fNfill;
01423    if (!fTree->GetUpdate()) return;
01424    if (fSelectedRows > fDraw+fTree->GetUpdate()) {
01425       if (fDraw) gPad->Modified();
01426       else       fObject->Draw(fOption.Data());
01427       gPad->Update();
01428       fDraw = fSelectedRows;
01429    }
01430 }
01431 
01432 //______________________________________________________________________________
01433 void TSelectorDraw::TakeEstimate()
01434 {
01435    // Estimate limits for 1-D, 2-D or 3-D objects.
01436 
01437    Int_t i;
01438    Double_t rmin[3],rmax[3];
01439    Double_t vminOld[4], vmaxOld[4];
01440    for (i = 0; i < fValSize && i < 4; i++) {
01441       vminOld[i] = fVmin[i];
01442       vmaxOld[i] = fVmax[i];
01443    }
01444    for(i=0;i<fValSize;++i){
01445       fVmin[i] = FLT_MAX;
01446       fVmax[i] = - FLT_MAX;
01447    }
01448    //__________________________1-D histogram_______________________
01449    if      (fAction ==  1) {
01450       TH1 *h1 = (TH1*)fObject;
01451       if (fObject->TestBit(TH1::kCanRebin)) {
01452          for (i=0;i<fNfill;i++) {
01453             if (fVmin[0] > fVal[0][i]) fVmin[0] = fVal[0][i];
01454             if (fVmax[0] < fVal[0][i]) fVmax[0] = fVal[0][i];
01455          }
01456          THLimitsFinder::GetLimitsFinder()->FindGoodLimits(h1,fVmin[0],fVmax[0]);
01457       }
01458       h1->FillN(fNfill, fVal[0], fW);
01459    //__________________________2-D histogram_______________________
01460    } else if (fAction ==  2) {
01461       TH2 *h2 = (TH2*)fObject;
01462       if (fObject->TestBit(TH1::kCanRebin)) {
01463          for (i=0;i<fNfill;i++) {
01464             if (fVmin[0] > fVal[0][i]) fVmin[0] = fVal[0][i];
01465             if (fVmax[0] < fVal[0][i]) fVmax[0] = fVal[0][i];
01466             if (fVmin[1] > fVal[1][i]) fVmin[1] = fVal[1][i];
01467             if (fVmax[1] < fVal[1][i]) fVmax[1] = fVal[1][i];
01468          }
01469          THLimitsFinder::GetLimitsFinder()->FindGoodLimits(h2,fVmin[1],fVmax[1],fVmin[0],fVmax[0]);
01470       }
01471       for(i=0;i<fNfill;i++) h2->Fill(fVal[1][i],fVal[0][i],fW[i]);
01472    //__________________________Profile histogram_______________________
01473    } else if (fAction ==  4) {
01474       TProfile *hp = (TProfile*)fObject;
01475       if (fObject->TestBit(TH1::kCanRebin)) {
01476          for (i=0;i<fNfill;i++) {
01477             if (fVmin[0] > fVal[0][i]) fVmin[0] = fVal[0][i];
01478             if (fVmax[0] < fVal[0][i]) fVmax[0] = fVal[0][i];
01479             if (fVmin[1] > fVal[1][i]) fVmin[1] = fVal[1][i];
01480             if (fVmax[1] < fVal[1][i]) fVmax[1] = fVal[1][i];
01481          }
01482          THLimitsFinder::GetLimitsFinder()->FindGoodLimits(hp,fVmin[1],fVmax[1]);
01483       }
01484       hp->FillN(fNfill, fVal[1], fVal[0], fW);
01485    //__________________________2D scatter plot_______________________
01486    } else if (fAction == 12) {
01487       TH2 *h2 = (TH2*)fObject;
01488       if (h2->TestBit(TH1::kCanRebin)) {
01489          for (i=0;i<fNfill;i++) {
01490             if (fVmin[0] > fVal[0][i]) fVmin[0] = fVal[0][i];
01491             if (fVmax[0] < fVal[0][i]) fVmax[0] = fVal[0][i];
01492             if (fVmin[1] > fVal[1][i]) fVmin[1] = fVal[1][i];
01493             if (fVmax[1] < fVal[1][i]) fVmax[1] = fVal[1][i];
01494          }
01495          THLimitsFinder::GetLimitsFinder()->FindGoodLimits(h2,fVmin[1],fVmax[1],fVmin[0],fVmax[0]);
01496       }
01497 
01498       if (!strstr(fOption.Data(),"same") && !strstr(fOption.Data(),"goff")) {
01499          if (!h2->TestBit(kCanDelete)) {
01500             // case like: T.Draw("y:x>>myhist")
01501             // we must draw a copy before filling the histogram h2=myhist
01502             // because h2 will be filled below and we do not want to show
01503             // the binned scatter-plot, the TGraph being better.
01504             TH1 *h2c = h2->DrawCopy(fOption.Data());
01505             h2c->SetStats(kFALSE);
01506          } else {
01507             // case like: T.Draw("y:x")
01508             // h2 is a temporary histogram (htemp). This histogram
01509             // will be automatically deleted by TPad::Clear
01510             h2->Draw();
01511          }
01512          gPad->Update();
01513       }
01514       TGraph *pm = new TGraph(fNfill,fVal[1],fVal[0]);
01515       pm->SetEditable(kFALSE);
01516       pm->SetBit(kCanDelete);
01517       pm->SetMarkerStyle(fTree->GetMarkerStyle());
01518       pm->SetMarkerColor(fTree->GetMarkerColor());
01519       pm->SetMarkerSize(fTree->GetMarkerSize());
01520       pm->SetLineColor(fTree->GetLineColor());
01521       pm->SetLineStyle(fTree->GetLineStyle());
01522       pm->SetFillColor(fTree->GetFillColor());
01523       pm->SetFillStyle(fTree->GetFillStyle());
01524       if (!fDraw && !strstr(fOption.Data(),"goff")) {
01525          if (fOption.Length() == 0 || strcasecmp(fOption.Data(),"same")==0) {
01526             pm->Draw("p");
01527          } 
01528          else {
01529             if (fOption.Contains("a")) {
01530                TString temp(fOption);
01531                temp.ReplaceAll("same","");
01532                if (temp.Contains("a")) {
01533                   if (h2->TestBit(kCanDelete)) {
01534                      // h2 will be deleted, the axis setting is delegated to only
01535                      // the TGraph.
01536                      h2 = 0;
01537                   }
01538                }
01539             }
01540             pm->Draw(fOption.Data());
01541          }
01542       }
01543       if (h2 && !h2->TestBit(kCanDelete)) {
01544          for (i=0;i<fNfill;i++) h2->Fill(fVal[1][i],fVal[0][i],fW[i]);
01545       }
01546    //__________________________3D scatter plot with option col_______________________
01547    } else if (fAction == 33) {
01548       TH2 *h2 = (TH2*)fObject;
01549       Bool_t process2 = kFALSE;
01550       if (h2->TestBit(TH1::kCanRebin)) {
01551          if (vminOld[2] == FLT_MAX)
01552             process2 = kTRUE;
01553          for (i = 0; i < fValSize && i < 4; i++) {
01554             fVmin[i] = vminOld[i];
01555             fVmax[i] = vmaxOld[i];
01556          }
01557          for (i=0;i<fNfill;i++) {
01558             if (fVmin[0] > fVal[0][i]) fVmin[0] = fVal[0][i];
01559             if (fVmax[0] < fVal[0][i]) fVmax[0] = fVal[0][i];
01560             if (fVmin[1] > fVal[1][i]) fVmin[1] = fVal[1][i];
01561             if (fVmax[1] < fVal[1][i]) fVmax[1] = fVal[1][i];
01562             if (process2) {
01563                if (fVmin[2] > fVal[2][i]) fVmin[2] = fVal[2][i];
01564                if (fVmax[2] < fVal[2][i]) fVmax[2] = fVal[2][i];
01565             }
01566          }
01567          THLimitsFinder::GetLimitsFinder()->FindGoodLimits(h2,fVmin[1],fVmax[1],fVmin[0],fVmax[0]);
01568       }
01569    //__________________________3D scatter plot_______________________
01570    } else if (fAction == 3 || fAction == 13) {
01571       TH3 *h3 = (TH3*)fObject;
01572       if (fObject->TestBit(TH1::kCanRebin)) {
01573          for (i=0;i<fNfill;i++) {
01574             if (fVmin[0] > fVal[0][i]) fVmin[0] = fVal[0][i];
01575             if (fVmax[0] < fVal[0][i]) fVmax[0] = fVal[0][i];
01576             if (fVmin[1] > fVal[1][i]) fVmin[1] = fVal[1][i];
01577             if (fVmax[1] < fVal[1][i]) fVmax[1] = fVal[1][i];
01578             if (fVmin[2] > fVal[2][i]) fVmin[2] = fVal[2][i];
01579             if (fVmax[2] < fVal[2][i]) fVmax[2] = fVal[2][i];
01580          }
01581          THLimitsFinder::GetLimitsFinder()->FindGoodLimits(h3,fVmin[2],fVmax[2],fVmin[1],fVmax[1],fVmin[0],fVmax[0]);
01582       }
01583       if (fAction == 3) {
01584          for (i=0;i<fNfill;i++) h3->Fill(fVal[2][i],fVal[1][i],fVal[0][i],fW[i]);
01585          return;
01586       }
01587       if (!strstr(fOption.Data(),"same") && !strstr(fOption.Data(),"goff")) {
01588          if (!h3->TestBit(kCanDelete)) {
01589             // case like: T.Draw("y:x>>myhist")
01590             // we must draw a copy before filling the histogram h3=myhist
01591             // because h3 will be filled below and we do not want to show
01592             // the binned scatter-plot, the TGraph being better.
01593             TH1 *h3c = h3->DrawCopy(fOption.Data());
01594             h3c->SetStats(kFALSE);
01595          } else {
01596             // case like: T.Draw("y:x")
01597             // h3 is a temporary histogram (htemp). This histogram
01598             // will be automatically deleted by TPad::Clear
01599             h3->Draw(fOption.Data());
01600          }
01601          gPad->Update();
01602       } else {
01603          rmin[0] = fVmin[2]; rmin[1] = fVmin[1]; rmin[2] = fVmin[0];
01604          rmax[0] = fVmax[2]; rmax[1] = fVmax[1]; rmax[2] = fVmax[0];
01605          gPad->Clear();
01606          gPad->Range(-1,-1,1,1);
01607          TView::CreateView(1, rmin,rmax);
01608       }
01609       TPolyMarker3D *pm3d = new TPolyMarker3D(fNfill);
01610       pm3d->SetMarkerStyle(fTree->GetMarkerStyle());
01611       pm3d->SetMarkerColor(fTree->GetMarkerColor());
01612       pm3d->SetMarkerSize(fTree->GetMarkerSize());
01613       for (i=0;i<fNfill;i++) { pm3d->SetPoint(i,fVal[2][i],fVal[1][i],fVal[0][i]);}
01614       if (!fDraw && !strstr(fOption.Data(),"goff")) pm3d->Draw();
01615       if (!h3->TestBit(kCanDelete)) {
01616          for (i=0;i<fNfill;i++) h3->Fill(fVal[2][i],fVal[1][i],fVal[0][i],fW[i]);
01617       }
01618 
01619    //__________________________2D Profile Histogram__________________
01620    } else if (fAction == 23) {
01621       TProfile2D *hp = (TProfile2D*)fObject;
01622       if (hp->TestBit(TH1::kCanRebin)) {
01623          for (i=0;i<fNfill;i++) {
01624             if (fVmin[0] > fVal[0][i]) fVmin[0] = fVal[0][i];
01625             if (fVmax[0] < fVal[0][i]) fVmax[0] = fVal[0][i];
01626             if (fVmin[1] > fVal[1][i]) fVmin[1] = fVal[1][i];
01627             if (fVmax[1] < fVal[1][i]) fVmax[1] = fVal[1][i];
01628             if (fVmin[2] > fVal[2][i]) fVmin[2] = fVal[2][i];
01629             if (fVmax[2] < fVal[2][i]) fVmax[2] = fVal[2][i];
01630          }
01631          THLimitsFinder::GetLimitsFinder()->FindGoodLimits(hp,fVmin[2],fVmax[2],fVmin[1],fVmax[1]);
01632       }
01633       for (i=0;i<fNfill;i++) hp->Fill(fVal[2][i],fVal[1][i],fVal[0][i],fW[i]);
01634    //__________________________4D scatter plot_______________________
01635    } else if (fAction == 40) {
01636       TH3 *h3 = (TH3*)fObject;
01637       if (fObject->TestBit(TH1::kCanRebin)) {
01638          for (i = 0; i < fValSize && i < 4; i++) {
01639             fVmin[i] = vminOld[i];
01640             fVmax[i] = vmaxOld[i];
01641          }
01642          for (i=0;i<fNfill;i++) {
01643             if (fVmin[0] > fVal[0][i]) fVmin[0] = fVal[0][i];
01644             if (fVmax[0] < fVal[0][i]) fVmax[0] = fVal[0][i];
01645             if (fVmin[1] > fVal[1][i]) fVmin[1] = fVal[1][i];
01646             if (fVmax[1] < fVal[1][i]) fVmax[1] = fVal[1][i];
01647             if (fVmin[2] > fVal[2][i]) fVmin[2] = fVal[2][i];
01648             if (fVmax[2] < fVal[2][i]) fVmax[2] = fVal[2][i];
01649             if (fVmin[3] > fVal[3][i]) fVmin[3] = fVal[3][i];
01650             if (fVmax[3] < fVal[3][i]) fVmax[3] = fVal[3][i];
01651          }
01652          THLimitsFinder::GetLimitsFinder()->FindGoodLimits(h3,fVmin[2],fVmax[2],fVmin[1],fVmax[1],fVmin[0],fVmax[0]);
01653       }
01654    }
01655    //__________________________Parallel coordinates plot / candle chart_______________________
01656    else if (fAction == 6 || fAction == 7){
01657       for(i=0;i<fDimension;++i){
01658          for (Long64_t entry=0;entry<fNfill;entry++){
01659             if (fVmin[i] > fVal[i][entry]) fVmin[i] = fVal[i][entry];
01660             if (fVmax[i] < fVal[i][entry]) fVmax[i] = fVal[i][entry];
01661          }
01662       }
01663    }
01664 }
01665 
01666 //______________________________________________________________________________
01667 void TSelectorDraw::Terminate()
01668 {
01669    // Called at the end of a loop on a TTree.
01670 
01671    if (fNfill) TakeAction();
01672 
01673    if ((fSelectedRows == 0) && (TestBit(kCustomHistogram) == 0)) fDraw = 1; // do not draw
01674 
01675    SetStatus(fSelectedRows);
01676 }

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