TChain.cxx

Go to the documentation of this file.
00001 // @(#)root/tree:$Id: TChain.cxx 37837 2011-01-21 16:52:15Z pcanal $
00002 // Author: Rene Brun   03/02/97
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 // TChain                                                               //
00015 //                                                                      //
00016 // A chain is a collection of files containg TTree objects.             //
00017 // When the chain is created, the first parameter is the default name   //
00018 // for the Tree to be processed later on.                               //
00019 //                                                                      //
00020 // Enter a new element in the chain via the TChain::Add function.       //
00021 // Once a chain is defined, one can use the normal TTree functions      //
00022 // to Draw,Scan,etc.                                                    //
00023 //                                                                      //
00024 // Use TChain::SetBranchStatus to activate one or more branches for all //
00025 // the trees in the chain.                                              //
00026 //                                                                      //
00027 //////////////////////////////////////////////////////////////////////////
00028 
00029 #include "TChain.h"
00030 
00031 #include "TBranch.h"
00032 #include "TBrowser.h"
00033 #include "TChainElement.h"
00034 #include "TClass.h"
00035 #include "TCut.h"
00036 #include "TError.h"
00037 #include "TMath.h"
00038 #include "TFile.h"
00039 #include "TFileInfo.h"
00040 #include "TFriendElement.h"
00041 #include "TLeaf.h"
00042 #include "TList.h"
00043 #include "TObjString.h"
00044 #include "TPluginManager.h"
00045 #include "TROOT.h"
00046 #include "TRegexp.h"
00047 #include "TSelector.h"
00048 #include "TSystem.h"
00049 #include "TTree.h"
00050 #include "TTreeCloner.h"
00051 #include "TTreeCache.h"
00052 #include "TUrl.h"
00053 #include "TVirtualIndex.h"
00054 #include "TEventList.h"
00055 #include "TEntryList.h"
00056 #include "TEntryListFromFile.h"
00057 #include "TFileStager.h"
00058 
00059 const Long64_t theBigNumber = Long64_t(1234567890)<<28;
00060 
00061 ClassImp(TChain)
00062 
00063 //______________________________________________________________________________
00064 TChain::TChain()
00065 : TTree()
00066 , fTreeOffsetLen(100)
00067 , fNtrees(0)
00068 , fTreeNumber(-1)
00069 , fTreeOffset(0)
00070 , fCanDeleteRefs(kFALSE)
00071 , fTree(0)
00072 , fFile(0)
00073 , fFiles(0)
00074 , fStatus(0)
00075 , fProofChain(0)
00076 {
00077    // -- Default constructor.
00078 
00079    fTreeOffset = new Long64_t[fTreeOffsetLen];
00080    fFiles = new TObjArray(fTreeOffsetLen);
00081    fStatus = new TList();
00082    fTreeOffset[0]  = 0;
00083    gDirectory->Remove(this);
00084    gROOT->GetListOfSpecials()->Add(this);
00085    fFile = 0;
00086    fDirectory = 0;
00087 
00088    // Reset PROOF-related bits
00089    ResetBit(kProofUptodate);
00090    ResetBit(kProofLite);
00091 
00092    // Add to the global list
00093    gROOT->GetListOfDataSets()->Add(this);
00094 
00095    // Make sure we are informed if the TFile is deleted.
00096    gROOT->GetListOfCleanups()->Add(this);
00097 }
00098 
00099 //______________________________________________________________________________
00100 TChain::TChain(const char* name, const char* title)
00101 :TTree(name, title)
00102 , fTreeOffsetLen(100)
00103 , fNtrees(0)
00104 , fTreeNumber(-1)
00105 , fTreeOffset(0)
00106 , fCanDeleteRefs(kFALSE)
00107 , fTree(0)
00108 , fFile(0)
00109 , fFiles(0)
00110 , fStatus(0)
00111 , fProofChain(0)
00112 {
00113    // -- Create a chain.
00114    //
00115    //   A TChain is a collection of TFile objects.
00116    //    the first parameter "name" is the name of the TTree object
00117    //    in the files added with Add.
00118    //   Use TChain::Add to add a new element to this chain.
00119    //
00120    //   In case the Tree is in a subdirectory, do, eg:
00121    //     TChain ch("subdir/treename");
00122    //
00123    //    Example:
00124    //  Suppose we have 3 files f1.root, f2.root and f3.root. Each file
00125    //  contains a TTree object named "T".
00126    //     TChain ch("T");  creates a chain to process a Tree called "T"
00127    //     ch.Add("f1.root");
00128    //     ch.Add("f2.root");
00129    //     ch.Add("f3.root");
00130    //     ch.Draw("x");
00131    //       The Draw function above will process the variable "x" in Tree "T"
00132    //       reading sequentially the 3 files in the chain ch.
00133    //
00134    //   The TChain data structure
00135    //       Each TChainElement has a name equal to the tree name of this TChain
00136    //       and a title equal to the file name. So, to loop over the
00137    //       TFiles that have been added to this chain:
00138    //
00139    //         TObjArray *fileElements=chain->GetListOfFiles();
00140    //         TIter next(fileElements);
00141    //         TChainElement *chEl=0;
00142    //         while (( chEl=(TChainElement*)next() )) {
00143    //            TFile f(chEl->GetTitle());
00144    //            ... do something with f ...
00145    //         }
00146 
00147    //
00148    //*-*
00149 
00150    fTreeOffset = new Long64_t[fTreeOffsetLen];
00151    fFiles = new TObjArray(fTreeOffsetLen);
00152    fStatus = new TList();
00153    fTreeOffset[0]  = 0;
00154    gDirectory->Remove(this);
00155    gROOT->GetListOfSpecials()->Add(this);
00156    fFile = 0;
00157    fDirectory = 0;
00158 
00159    // Reset PROOF-related bits
00160    ResetBit(kProofUptodate);
00161    ResetBit(kProofLite);
00162 
00163    // Add to the global list
00164    gROOT->GetListOfDataSets()->Add(this);
00165 
00166    // Make sure we are informed if the TFile is deleted.
00167    gROOT->GetListOfCleanups()->Add(this);
00168 }
00169 
00170 //______________________________________________________________________________
00171 TChain::~TChain()
00172 {
00173    // -- Destructor.
00174    gROOT->GetListOfCleanups()->Remove(this);
00175 
00176    SafeDelete(fProofChain);
00177    fStatus->Delete();
00178    delete fStatus;
00179    fStatus = 0;
00180    fFiles->Delete();
00181    delete fFiles;
00182    fFiles = 0;
00183    delete fFile;
00184    fFile = 0;
00185    // Note: We do *not* own the tree.
00186    fTree = 0;
00187    delete[] fTreeOffset;
00188    fTreeOffset = 0;
00189 
00190    gROOT->GetListOfSpecials()->Remove(this);
00191 
00192    // Remove from the global list
00193    gROOT->GetListOfDataSets()->Remove(this);
00194 
00195    // This is the same as fFile, don't delete it a second time.
00196    fDirectory = 0;
00197 }
00198 
00199 //______________________________________________________________________________
00200 Int_t TChain::Add(TChain* chain)
00201 {
00202    // -- Add all files referenced by the passed chain to this chain.
00203    // The function returns the total number of files connected.
00204 
00205    if (!chain) return 0;
00206 
00207    // Check for enough space in fTreeOffset.
00208    if ((fNtrees + chain->GetNtrees()) >= fTreeOffsetLen) {
00209       fTreeOffsetLen += 2 * chain->GetNtrees();
00210       Long64_t* trees = new Long64_t[fTreeOffsetLen];
00211       for (Int_t i = 0; i <= fNtrees; i++) {
00212          trees[i] = fTreeOffset[i];
00213       }
00214       delete[] fTreeOffset;
00215       fTreeOffset = trees;
00216    }
00217    chain->GetEntries(); //to force the computation of nentries
00218    TIter next(chain->GetListOfFiles());
00219    Int_t nf = 0;
00220    TChainElement* element = 0;
00221    while ((element = (TChainElement*) next())) {
00222       Long64_t nentries = element->GetEntries();
00223       if (fTreeOffset[fNtrees] == theBigNumber) {
00224          fTreeOffset[fNtrees+1] = theBigNumber;
00225       } else {
00226          fTreeOffset[fNtrees+1] = fTreeOffset[fNtrees] + nentries;
00227       }
00228       fNtrees++;
00229       fEntries += nentries;
00230       TChainElement* newelement = new TChainElement(element->GetName(), element->GetTitle());
00231       newelement->SetPacketSize(element->GetPacketSize());
00232       newelement->SetNumberEntries(nentries);
00233       fFiles->Add(newelement);
00234       nf++;
00235    }
00236    if (fProofChain)
00237       // This updates the proxy chain when we will really use PROOF
00238       ResetBit(kProofUptodate);
00239 
00240    return nf;
00241 }
00242 
00243 //______________________________________________________________________________
00244 Int_t TChain::Add(const char* name, Long64_t nentries /* = kBigNumber */)
00245 {
00246    // -- Add a new file to this chain.
00247    //
00248    // Argument name may have the following format:
00249    //   //machine/file_name.root/subdir/tree_name
00250    // machine, subdir and tree_name are optional. If tree_name is missing,
00251    // the chain name will be assumed.
00252    // In the file name part (but not in preceding directories) wildcarding
00253    // notation may be used, eg. specifying "xxx*.root" adds all files starting
00254    // with xxx in the current file system directory.
00255    // NB. To add all the files of a TChain to a chain, use Add(TChain *chain).
00256    //
00257    //    A- if nentries <= 0, the file is connected and the tree header read
00258    //       in memory to get the number of entries.
00259    //
00260    //    B- if (nentries > 0, the file is not connected, nentries is assumed to be
00261    //       the number of entries in the file. In this case, no check is made that
00262    //       the file exists and the Tree existing in the file. This second mode
00263    //       is interesting in case the number of entries in the file is already stored
00264    //       in a run data base for example.
00265    //
00266    //    C- if (nentries == kBigNumber) (default), the file is not connected.
00267    //       the number of entries in each file will be read only when the file
00268    //       will need to be connected to read an entry.
00269    //       This option is the default and very efficient if one process
00270    //       the chain sequentially. Note that in case TChain::GetEntry(entry)
00271    //       is called and entry refers to an entry in the 3rd file, for example,
00272    //       this forces the Tree headers in the first and second file
00273    //       to be read to find the number of entries in these files.
00274    //       Note that if one calls TChain::GetEntriesFast() after having created
00275    //       a chain with this default, GetEntriesFast will return kBigNumber!
00276    //       TChain::GetEntries will force of the Tree headers in the chain to be
00277    //       read to read the number of entries in each Tree.
00278    //
00279    //
00280    //    D- The TChain data structure
00281    //       Each TChainElement has a name equal to the tree name of this TChain
00282    //       and a title equal to the file name. So, to loop over the
00283    //       TFiles that have been added to this chain:
00284    //
00285    //         TObjArray *fileElements=chain->GetListOfFiles();
00286    //         TIter next(fileElements);
00287    //         TChainElement *chEl=0;
00288    //         while (( chEl=(TChainElement*)next() )) {
00289    //            TFile f(chEl->GetTitle());
00290    //            ... do something with f ...
00291    //         }
00292    //
00293    // The function returns the total number of files connected.
00294 
00295    // case with one single file
00296    if (!TString(name).MaybeWildcard()) {
00297       return AddFile(name, nentries);
00298    }
00299 
00300    // wildcarding used in name
00301    Int_t nf = 0;
00302    TString basename(name);
00303 
00304    Int_t dotslashpos = -1;
00305    {
00306       Int_t next_dot = basename.Index(".root");
00307       while(next_dot>=0) {
00308          dotslashpos = next_dot;
00309          next_dot = basename.Index(".root",dotslashpos+1);
00310       }
00311       if (basename[dotslashpos+5]!='/') {
00312          // We found the 'last' .root in the name and it is not followed by
00313          // a '/', so the tree name is _not_ specified in the name.
00314          dotslashpos = -1;
00315       }
00316    }
00317    //Int_t dotslashpos = basename.Index(".root/");
00318    TString behind_dot_root;
00319    if (dotslashpos>=0) {
00320       // Copy the tree name specification
00321       behind_dot_root = basename(dotslashpos+6,basename.Length()-dotslashpos+6);
00322       // and remove it from basename
00323       basename.Remove(dotslashpos+5);
00324    }
00325 
00326    Int_t slashpos = basename.Last('/');
00327    TString directory;
00328    if (slashpos>=0) {
00329       directory = basename(0,slashpos); // Copy the directory name
00330       basename.Remove(0,slashpos+1);      // and remove it from basename
00331    } else {
00332       directory = gSystem->UnixPathName(gSystem->WorkingDirectory());
00333    }
00334 
00335    const char *file;
00336    const char *epath = gSystem->ExpandPathName(directory.Data());
00337    void *dir = gSystem->OpenDirectory(epath);
00338    delete [] epath;
00339    if (dir) {
00340       //create a TList to store the file names (not yet sorted)
00341       TList l;
00342       TRegexp re(basename,kTRUE);
00343       while ((file = gSystem->GetDirEntry(dir))) {
00344          if (!strcmp(file,".") || !strcmp(file,"..")) continue;
00345          TString s = file;
00346          if ( (basename!=file) && s.Index(re) == kNPOS) continue;
00347          l.Add(new TObjString(file));
00348       }
00349       gSystem->FreeDirectory(dir);
00350       //sort the files in alphanumeric order
00351       l.Sort();
00352       TIter next(&l);
00353       TObjString *obj;
00354       while ((obj = (TObjString*)next())) {
00355          file = obj->GetName();
00356          if (behind_dot_root.Length() != 0)
00357             nf += AddFile(Form("%s/%s/%s",directory.Data(),file,behind_dot_root.Data()),nentries);
00358          else
00359             nf += AddFile(Form("%s/%s",directory.Data(),file),nentries);
00360       }
00361       l.Delete();
00362    }
00363    if (fProofChain)
00364       // This updates the proxy chain when we will really use PROOF
00365       ResetBit(kProofUptodate);
00366 
00367    return nf;
00368 }
00369 
00370 //______________________________________________________________________________
00371 Int_t TChain::AddFile(const char* name, Long64_t nentries /* = kBigNumber */, const char* tname /* = "" */)
00372 {
00373    // -- Add a new file to this chain.
00374    //
00375    //    If tname is specified, the chain will load the tree named tname
00376    //    from the file, otherwise the original treename specified in the
00377    //    TChain constructor will be used.
00378    //
00379    // A. If nentries <= 0, the file is opened and the tree header read
00380    //    into memory to get the number of entries.
00381    //
00382    // B. If nentries > 0, the file is not opened, and nentries is assumed
00383    //    to be the number of entries in the file. In this case, no check
00384    //    is made that the file exists nor that the tree exists in the file.
00385    //    This second mode is interesting in case the number of entries in
00386    //    the file is already stored in a run database for example.
00387    //
00388    // C. If nentries == kBigNumber (default), the file is not opened.
00389    //    The number of entries in each file will be read only when the file
00390    //    is opened to read an entry.  This option is the default and very
00391    //    efficient if one processes the chain sequentially.  Note that in
00392    //    case GetEntry(entry) is called and entry refers to an entry in the
00393    //    third file, for example, this forces the tree headers in the first
00394    //    and second file to be read to find the number of entries in those
00395    //    files.  Note that if one calls GetEntriesFast() after having created
00396    //    a chain with this default, GetEntriesFast() will return kBigNumber!
00397    //    Using the GetEntries() function instead will force all of the tree
00398    //    headers in the chain to be read to read the number of entries in
00399    //    each tree.
00400    //
00401    // D. The TChain data structure
00402    //    Each TChainElement has a name equal to the tree name of this TChain
00403    //    and a title equal to the file name. So, to loop over the
00404    //    TFiles that have been added to this chain:
00405    //
00406    //      TObjArray *fileElements=chain->GetListOfFiles();
00407    //      TIter next(fileElements);
00408    //      TChainElement *chEl=0;
00409    //      while (( chEl=(TChainElement*)next() )) {
00410    //         TFile f(chEl->GetTitle());
00411    //         ... do something with f ...
00412    //      }
00413    //
00414    // The function returns 1 if the file is successfully connected, 0 otherwise.
00415 
00416 
00417    const char *treename = GetName();
00418    if (tname && strlen(tname) > 0) treename = tname;
00419    char *dot = 0;
00420    {
00421       char *nextdot =  (char*)strstr(name,".root");
00422       while (nextdot) {
00423          dot = nextdot;
00424          nextdot = (char*)strstr(dot+1,".root");
00425       }
00426    }
00427    //the ".root" is mandatory only if one wants to specify a treename
00428    //if (!dot) {
00429    //   Error("AddFile","a chain element name must contain the string .root");
00430    //   return 0;
00431    //}
00432 
00433    //Check enough space in fTreeOffset
00434    if (fNtrees+1 >= fTreeOffsetLen) {
00435       fTreeOffsetLen *= 2;
00436       Long64_t *trees = new Long64_t[fTreeOffsetLen];
00437       for (Int_t i=0;i<=fNtrees;i++) trees[i] = fTreeOffset[i];
00438       delete [] fTreeOffset;
00439       fTreeOffset = trees;
00440    }
00441 
00442    //Search for a a slash between the .root and the end
00443    Int_t nch = strlen(name) + strlen(treename);
00444    char *filename = new char[nch+1];
00445    strlcpy(filename,name,nch+1); 
00446    if (dot) {
00447       char *pos = filename + (dot-name) + 5;
00448       while (*pos) {
00449          if (*pos == '/') {
00450             treename = pos+1;
00451             *pos = 0;
00452             break;
00453          }
00454          pos++;
00455       }
00456    }
00457 
00458    // Open the file to get the number of entries.
00459    Int_t pksize = 0;
00460    if (nentries <= 0) {
00461       TFile* file;
00462       {
00463          TDirectory::TContext ctxt(0);
00464          file = TFile::Open(filename);
00465       }
00466       if (!file || file->IsZombie()) {
00467          delete file;
00468          file = 0;
00469          delete[] filename;
00470          filename = 0;
00471          return 0;
00472       }
00473 
00474       // Check that tree with the right name exists in the file.
00475       // Note: We are not the owner of obj, the file is!
00476       TObject* obj = file->Get(treename);
00477       if (!obj || !obj->InheritsFrom(TTree::Class())) {
00478          Error("AddFile", "cannot find tree with name %s in file %s", treename, filename);
00479          delete file;
00480          file = 0;
00481          delete[] filename;
00482          filename = 0;
00483          return 0;
00484       }
00485       TTree* tree = (TTree*) obj;
00486       nentries = tree->GetEntries();
00487       pksize = tree->GetPacketSize();
00488       // Note: This deletes the tree we fetched.
00489       delete file;
00490       file = 0;
00491    }
00492 
00493    if (nentries > 0) {
00494       if (nentries != kBigNumber) {
00495          fTreeOffset[fNtrees+1] = fTreeOffset[fNtrees] + nentries;
00496          fEntries += nentries;
00497       } else {
00498          fTreeOffset[fNtrees+1] = theBigNumber;
00499          fEntries = nentries;
00500       }
00501       fNtrees++;
00502 
00503       TChainElement* element = new TChainElement(treename, filename);
00504       element->SetPacketSize(pksize);
00505       element->SetNumberEntries(nentries);
00506       fFiles->Add(element);
00507    } else {
00508       Warning("AddFile", "Adding tree with no entries from file: %s", filename);
00509    }
00510 
00511    delete [] filename;
00512    if (fProofChain)
00513       // This updates the proxy chain when we will really use PROOF
00514       ResetBit(kProofUptodate);
00515 
00516    return 1;
00517 }
00518 
00519 //______________________________________________________________________________
00520 Int_t TChain::AddFileInfoList(TCollection* filelist, Long64_t nfiles /* = kBigNumber */)
00521 {
00522    // Add all files referenced in the list to the chain. The object type in the
00523    // list must be either TFileInfo or TObjString or TUrl .
00524    // The function return 1 if successful, 0 otherwise.
00525    if (!filelist)
00526       return 0;
00527    TIter next(filelist);
00528 
00529    TObject *o = 0;
00530    Long64_t cnt=0;
00531    while ((o = next())) {
00532       // Get the url
00533       TString cn = o->ClassName();
00534       const char *url = 0;
00535       if (cn == "TFileInfo") {
00536          TFileInfo *fi = (TFileInfo *)o;
00537          url = (fi->GetCurrentUrl()) ? fi->GetCurrentUrl()->GetUrl() : 0;
00538          if (!url) {
00539             Warning("AddFileInfoList", "found TFileInfo with empty Url - ignoring");
00540             continue;
00541          }
00542       } else if (cn == "TUrl") {
00543          url = ((TUrl*)o)->GetUrl();
00544       } else if (cn == "TObjString") {
00545          url = ((TObjString*)o)->GetName();
00546       }
00547       if (!url) {
00548          Warning("AddFileInfoList", "object is of type %s : expecting TFileInfo, TUrl"
00549                                  " or TObjString - ignoring", o->ClassName());
00550          continue;
00551       }
00552       // Good entry
00553       cnt++;
00554       AddFile(url);
00555       if (cnt >= nfiles)
00556          break;
00557    }
00558    if (fProofChain) {
00559       // This updates the proxy chain when we will really use PROOF
00560       ResetBit(kProofUptodate);
00561    }
00562 
00563    return 1;
00564 }
00565 
00566 //______________________________________________________________________________
00567 TFriendElement* TChain::AddFriend(const char* chain, const char* dummy /* = "" */)
00568 {
00569    // -- Add a TFriendElement to the list of friends of this chain.
00570    //
00571    // A TChain has a list of friends similar to a tree (see TTree::AddFriend).
00572    // You can add a friend to a chain with the TChain::AddFriend method, and you
00573    // can retrieve the list of friends with TChain::GetListOfFriends.
00574    // This example has four chains each has 20 ROOT trees from 20 ROOT files.
00575    //
00576    // TChain ch("t"); // a chain with 20 trees from 20 files
00577    // TChain ch1("t1");
00578    // TChain ch2("t2");
00579    // TChain ch3("t3");
00580    // Now we can add the friends to the first chain.
00581    //
00582    // ch.AddFriend("t1")
00583    // ch.AddFriend("t2")
00584    // ch.AddFriend("t3")
00585    //
00586    //Begin_Html
00587    /*
00588    <img src="gif/chain_friend.gif">
00589    */
00590    //End_Html
00591    //
00592    // The parameter is the name of friend chain (the name of a chain is always
00593    // the name of the tree from which it was created).
00594    // The original chain has access to all variable in its friends.
00595    // We can use the TChain::Draw method as if the values in the friends were
00596    // in the original chain.
00597    // To specify the chain to use in the Draw method, use the syntax:
00598    //
00599    // <chainname>.<branchname>.<varname>
00600    // If the variable name is enough to uniquely identify the variable, you can
00601    // leave out the chain and/or branch name.
00602    // For example, this generates a 3-d scatter plot of variable "var" in the
00603    // TChain ch versus variable v1 in TChain t1 versus variable v2 in TChain t2.
00604    //
00605    // ch.Draw("var:t1.v1:t2.v2");
00606    // When a TChain::Draw is executed, an automatic call to TTree::AddFriend
00607    // connects the trees in the chain. When a chain is deleted, its friend
00608    // elements are also deleted.
00609    //
00610    // The number of entries in the friend must be equal or greater to the number
00611    // of entries of the original chain. If the friend has fewer entries a warning
00612    // is given and the resulting histogram will have missing entries.
00613    // For additional information see TTree::AddFriend.
00614 
00615    if (!fFriends) {
00616       fFriends = new TList();
00617    }
00618    TFriendElement* fe = new TFriendElement(this, chain, dummy);
00619 
00620    R__ASSERT(fe); // There used to be a "if (fe)" test ... Keep this assert until we are sure that fe is never null
00621 
00622    fFriends->Add(fe);
00623 
00624    if (fProofChain)
00625       // This updates the proxy chain when we will really use PROOF
00626       ResetBit(kProofUptodate);
00627 
00628    // We need to invalidate the loading of the current tree because its list
00629    // of real friends is now obsolete.  It is repairable only from LoadTree.
00630    fTreeNumber = -1;
00631 
00632    TTree* tree = fe->GetTree();
00633    if (!tree) {
00634       Warning("AddFriend", "Unknown TChain %s", chain);
00635    }
00636    return fe;
00637 }
00638 
00639 //______________________________________________________________________________
00640 TFriendElement* TChain::AddFriend(const char* chain, TFile* dummy)
00641 {
00642    // -- Add the whole chain or tree as a friend of this chain.
00643 
00644    if (!fFriends) fFriends = new TList();
00645    TFriendElement *fe = new TFriendElement(this,chain,dummy);
00646 
00647    R__ASSERT(fe); // There used to be a "if (fe)" test ... Keep this assert until we are sure that fe is never null
00648 
00649    fFriends->Add(fe);
00650 
00651    if (fProofChain)
00652       // This updates the proxy chain when we will really use PROOF
00653       ResetBit(kProofUptodate);
00654 
00655    // We need to invalidate the loading of the current tree because its list
00656    // of real friend is now obsolete.  It is repairable only from LoadTree
00657    fTreeNumber = -1;
00658 
00659    TTree *t = fe->GetTree();
00660    if (!t) {
00661       Warning("AddFriend","Unknown TChain %s",chain);
00662    }
00663    return fe;
00664 }
00665 
00666 //______________________________________________________________________________
00667 TFriendElement* TChain::AddFriend(TTree* chain, const char* alias, Bool_t /* warn = kFALSE */)
00668 {
00669    // -- Add the whole chain or tree as a friend of this chain.
00670 
00671    if (!fFriends) fFriends = new TList();
00672    TFriendElement *fe = new TFriendElement(this,chain,alias);
00673    R__ASSERT(fe);
00674 
00675    fFriends->Add(fe);
00676 
00677    if (fProofChain)
00678       // This updates the proxy chain when we will really use PROOF
00679       ResetBit(kProofUptodate);
00680 
00681    // We need to invalidate the loading of the current tree because its list
00682    // of real friend is now obsolete.  It is repairable only from LoadTree
00683    fTreeNumber = -1;
00684 
00685    TTree *t = fe->GetTree();
00686    if (!t) {
00687       Warning("AddFriend","Unknown TChain %s",chain->GetName());
00688    }
00689    return fe;
00690 }
00691 
00692 //______________________________________________________________________________
00693 void TChain::Browse(TBrowser* b)
00694 {
00695    // -- Browse the contents of the chain.
00696 
00697    TTree::Browse(b);
00698 }
00699 
00700 //_______________________________________________________________________
00701 void TChain::CanDeleteRefs(Bool_t flag /* = kTRUE */)
00702 {
00703    // When closing a file during the chain processing, the file
00704    // may be closed with option "R" if flag is set to kTRUE.
00705    // by default flag is kTRUE.
00706    // When closing a file with option "R", all TProcessIDs referenced by this
00707    // file are deleted.
00708    // Calling TFile::Close("R") might be necessary in case one reads a long list
00709    // of files having TRef, writing some of the referenced objects or TRef
00710    // to a new file. If the TRef or referenced objects of the file being closed
00711    // will not be referenced again, it is possible to minimize the size
00712    // of the TProcessID data structures in memory by forcing a delete of
00713    // the unused TProcessID.
00714 
00715    fCanDeleteRefs = flag;
00716 }
00717 
00718 //_______________________________________________________________________
00719 void TChain::CreatePackets()
00720 {
00721    // -- Initialize the packet descriptor string.
00722 
00723    TIter next(fFiles);
00724    TChainElement* element = 0;
00725    while ((element = (TChainElement*) next())) {
00726       element->CreatePackets();
00727    }
00728 }
00729 
00730 //______________________________________________________________________________
00731 void TChain::DirectoryAutoAdd(TDirectory * /* dir */)
00732 {
00733    // Override the TTree::DirectoryAutoAdd behavior:
00734    // we never auto add.
00735 
00736 }
00737 
00738 //______________________________________________________________________________
00739 Long64_t TChain::Draw(const char* varexp, const TCut& selection,
00740                       Option_t* option, Long64_t nentries, Long64_t firstentry)
00741 {
00742    // Draw expression varexp for selected entries.
00743    // Returns -1 in case of error or number of selected events in case of success.
00744    //
00745    // This function accepts TCut objects as arguments.
00746    // Useful to use the string operator +, example:
00747    //    ntuple.Draw("x",cut1+cut2+cut3);
00748    //
00749    if (fProofChain) {
00750       // Make sure the element list is uptodate
00751       if (!TestBit(kProofUptodate))
00752          SetProof(kTRUE, kTRUE);
00753       fProofChain->SetEventList(fEventList);
00754       fProofChain->SetEntryList(fEntryList);
00755       return fProofChain->Draw(varexp, selection, option, nentries, firstentry);
00756    }
00757 
00758    return TChain::Draw(varexp, selection.GetTitle(), option, nentries, firstentry);
00759 }
00760 
00761 //______________________________________________________________________________
00762 Long64_t TChain::Draw(const char* varexp, const char* selection,
00763                       Option_t* option,Long64_t nentries, Long64_t firstentry)
00764 {
00765    // Process all entries in this chain and draw histogram corresponding to
00766    // expression varexp.
00767    // Returns -1 in case of error or number of selected events in case of success.
00768 
00769    if (fProofChain) {
00770       // Make sure the element list is uptodate
00771       if (!TestBit(kProofUptodate))
00772          SetProof(kTRUE, kTRUE);
00773       fProofChain->SetEventList(fEventList);
00774       fProofChain->SetEntryList(fEntryList);
00775       return fProofChain->Draw(varexp, selection, option, nentries, firstentry);
00776    }
00777    GetPlayer();
00778    if (LoadTree(firstentry) < 0) return 0;
00779    return TTree::Draw(varexp,selection,option,nentries,firstentry);
00780 }
00781 
00782 //______________________________________________________________________________
00783 TBranch* TChain::FindBranch(const char* branchname)
00784 {
00785    // -- See TTree::GetReadEntry().
00786 
00787    if (fProofChain && !(fProofChain->TestBit(kProofLite))) {
00788       // Make sure the element list is uptodate
00789       if (!TestBit(kProofUptodate))
00790          SetProof(kTRUE, kTRUE);
00791       return fProofChain->FindBranch(branchname);
00792    }
00793    if (fTree) {
00794       return fTree->FindBranch(branchname);
00795    }
00796    LoadTree(0);
00797    if (fTree) {
00798       return fTree->FindBranch(branchname);
00799    }
00800    return 0;
00801 }
00802 
00803 //______________________________________________________________________________
00804 TLeaf* TChain::FindLeaf(const char* searchname)
00805 {
00806    // -- See TTree::GetReadEntry().
00807 
00808    if (fProofChain && !(fProofChain->TestBit(kProofLite))) {
00809       // Make sure the element list is uptodate
00810       if (!TestBit(kProofUptodate))
00811          SetProof(kTRUE, kTRUE);
00812       return fProofChain->FindLeaf(searchname);
00813    }
00814    if (fTree) {
00815       return fTree->FindLeaf(searchname);
00816    }
00817    LoadTree(0);
00818    if (fTree) {
00819       return fTree->FindLeaf(searchname);
00820    }
00821    return 0;
00822 }
00823 
00824 //______________________________________________________________________________
00825 const char* TChain::GetAlias(const char* aliasName) const
00826 {
00827    // -- Returns the expanded value of the alias.  Search in the friends if any.
00828 
00829    const char* alias = TTree::GetAlias(aliasName);
00830    if (alias) {
00831       return alias;
00832    }
00833    if (fTree) {
00834       return fTree->GetAlias(aliasName);
00835    }
00836    const_cast<TChain*>(this)->LoadTree(0);
00837    if (fTree) {
00838       return fTree->GetAlias(aliasName);
00839    }
00840    return 0;
00841 }
00842 
00843 //______________________________________________________________________________
00844 TBranch* TChain::GetBranch(const char* name)
00845 {
00846    // -- Return pointer to the branch name in the current tree.
00847 
00848    if (fProofChain && !(fProofChain->TestBit(kProofLite))) {
00849       // Make sure the element list is uptodate
00850       if (!TestBit(kProofUptodate))
00851          SetProof(kTRUE, kTRUE);
00852       return fProofChain->GetBranch(name);
00853    }
00854    if (fTree) {
00855       return fTree->GetBranch(name);
00856    }
00857    LoadTree(0);
00858    if (fTree) {
00859       return fTree->GetBranch(name);
00860    }
00861    return 0;
00862 }
00863 
00864 //______________________________________________________________________________
00865 Bool_t TChain::GetBranchStatus(const char* branchname) const
00866 {
00867    // -- See TTree::GetReadEntry().
00868    if (fProofChain && !(fProofChain->TestBit(kProofLite))) {
00869       // Make sure the element list is uptodate
00870       if (!TestBit(kProofUptodate))
00871          Warning("GetBranchStatus", "PROOF proxy not up-to-date:"
00872                                     " run TChain::SetProof(kTRUE, kTRUE) first");
00873       return fProofChain->GetBranchStatus(branchname);
00874    }
00875    return TTree::GetBranchStatus(branchname);
00876 }
00877 
00878 //______________________________________________________________________________
00879 Long64_t TChain::GetChainEntryNumber(Long64_t entry) const
00880 {
00881    // -- Return absolute entry number in the chain.
00882    // The input parameter entry is the entry number in
00883    // the current tree of this chain.
00884 
00885    return entry + fTreeOffset[fTreeNumber];
00886 }
00887 
00888 //______________________________________________________________________________
00889 Long64_t TChain::GetEntries() const
00890 {
00891    // -- Return the total number of entries in the chain.
00892    // In case the number of entries in each tree is not yet known,
00893    // the offset table is computed.
00894 
00895    if (fProofChain && !(fProofChain->TestBit(kProofLite))) {
00896       // Make sure the element list is uptodate
00897       if (!TestBit(kProofUptodate))
00898          Warning("GetEntries", "PROOF proxy not up-to-date:"
00899                                " run TChain::SetProof(kTRUE, kTRUE) first");
00900       return fProofChain->GetEntries();
00901    }
00902    if (fEntries >= theBigNumber || fEntries==kBigNumber) {
00903       const_cast<TChain*>(this)->LoadTree(theBigNumber-1);
00904    }
00905    return fEntries;
00906 }
00907 
00908 //______________________________________________________________________________
00909 Int_t TChain::GetEntry(Long64_t entry, Int_t getall)
00910 {
00911    // -- Get entry from the file to memory.
00912    //
00913    //     getall = 0 : get only active branches
00914    //     getall = 1 : get all branches
00915    //
00916    // Return the total number of bytes read,
00917    // 0 bytes read indicates a failure.
00918 
00919    Long64_t treeReadEntry = LoadTree(entry);
00920    if (treeReadEntry < 0) {
00921       return 0;
00922    }
00923    if (!fTree) {
00924       return 0;
00925    }
00926    return fTree->GetEntry(treeReadEntry, getall);
00927 }
00928 
00929 //______________________________________________________________________________
00930 Long64_t TChain::GetEntryNumber(Long64_t entry) const
00931 {
00932    // -- Return entry number corresponding to entry.
00933    //
00934    // if no TEntryList set returns entry
00935    // else returns entry #entry from this entry list and
00936    // also computes the global entry number (loads all tree headers)
00937 
00938 
00939    if (fEntryList){
00940       Int_t treenum = 0;
00941       Long64_t localentry = fEntryList->GetEntryAndTree(entry, treenum);
00942       //find the global entry number
00943       //same const_cast as in the GetEntries() function
00944       if (localentry<0) return -1;
00945       if (treenum != fTreeNumber){
00946          if (fTreeOffset[treenum]==theBigNumber){
00947             for (Int_t i=0; i<=treenum; i++){
00948                if (fTreeOffset[i]==theBigNumber)
00949                   (const_cast<TChain*>(this))->LoadTree(fTreeOffset[i-1]);
00950             }
00951          }
00952          //(const_cast<TChain*>(this))->LoadTree(fTreeOffset[treenum]);
00953       }
00954       Long64_t globalentry = fTreeOffset[treenum] + localentry;
00955       return globalentry;
00956    }
00957    return entry;
00958 }
00959 
00960 //______________________________________________________________________________
00961 Int_t TChain::GetEntryWithIndex(Int_t major, Int_t minor)
00962 {
00963    // -- Return entry corresponding to major and minor number.
00964    //
00965    //  The function returns the total number of bytes read.
00966    //  If the Tree has friend trees, the corresponding entry with
00967    //  the index values (major,minor) is read. Note that the master Tree
00968    //  and its friend may have different entry serial numbers corresponding
00969    //  to (major,minor).
00970 
00971    Long64_t serial = GetEntryNumberWithIndex(major, minor);
00972    if (serial < 0) return -1;
00973    return GetEntry(serial);
00974 }
00975 
00976 //______________________________________________________________________________
00977 TFile* TChain::GetFile() const
00978 {
00979    // -- Return a pointer to the current file.
00980    // If no file is connected, the first file is automatically loaded.
00981 
00982    if (fFile) {
00983       return fFile;
00984    }
00985    // Force opening the first file in the chain.
00986    const_cast<TChain*>(this)->LoadTree(0);
00987    return fFile;
00988 }
00989 
00990 //______________________________________________________________________________
00991 TLeaf* TChain::GetLeaf(const char* name)
00992 {
00993    // -- Return a pointer to the leaf name in the current tree.
00994 
00995    if (fProofChain && !(fProofChain->TestBit(kProofLite))) {
00996       // Make sure the element list is uptodate
00997       if (!TestBit(kProofUptodate))
00998          SetProof(kTRUE, kTRUE);
00999       return fProofChain->GetLeaf(name);
01000    }
01001    if (fTree) {
01002       return fTree->GetLeaf(name);
01003    }
01004    LoadTree(0);
01005    if (fTree) {
01006       return fTree->GetLeaf(name);
01007    }
01008    return 0;
01009 }
01010 
01011 //______________________________________________________________________________
01012 TObjArray* TChain::GetListOfBranches()
01013 {
01014    // -- Return a pointer to the list of branches of the current tree.
01015    //
01016    // Warning: If there is no current TTree yet, this routine will open the
01017    //     first in the chain.
01018    //
01019    // Returns 0 on failure.
01020 
01021    if (fProofChain && !(fProofChain->TestBit(kProofLite))) {
01022       // Make sure the element list is uptodate
01023       if (!TestBit(kProofUptodate))
01024          SetProof(kTRUE, kTRUE);
01025       return fProofChain->GetListOfBranches();
01026    }
01027    if (fTree) {
01028       return fTree->GetListOfBranches();
01029    }
01030    LoadTree(0);
01031    if (fTree) {
01032       return fTree->GetListOfBranches();
01033    }
01034    return 0;
01035 }
01036 
01037 //______________________________________________________________________________
01038 TObjArray* TChain::GetListOfLeaves()
01039 {
01040    // -- Return a pointer to the list of leaves of the current tree.
01041    //
01042    // Warning: May set the current tree!
01043    //
01044 
01045    if (fProofChain && !(fProofChain->TestBit(kProofLite))) {
01046       // Make sure the element list is uptodate
01047       if (!TestBit(kProofUptodate))
01048          SetProof(kTRUE, kTRUE);
01049       return fProofChain->GetListOfLeaves();
01050    }
01051    if (fTree) {
01052       return fTree->GetListOfLeaves();
01053    }
01054    LoadTree(0);
01055    if (fTree) {
01056       return fTree->GetListOfLeaves();
01057    }
01058    return 0;
01059 }
01060 
01061 //______________________________________________________________________________
01062 Double_t TChain::GetMaximum(const char* columname)
01063 {
01064    // -- Return maximum of column with name columname.
01065 
01066    Double_t theMax = -FLT_MAX;
01067    for (Int_t file = 0; file < fNtrees; file++) {
01068       Long64_t first = fTreeOffset[file];
01069       LoadTree(first);
01070       Double_t curmax = fTree->GetMaximum(columname);
01071       if (curmax > theMax) {
01072          theMax = curmax;
01073       }
01074    }
01075    return theMax;
01076 }
01077 
01078 //______________________________________________________________________________
01079 Double_t TChain::GetMinimum(const char* columname)
01080 {
01081    // -- Return minimum of column with name columname.
01082 
01083    Double_t theMin = FLT_MAX;
01084    for (Int_t file = 0; file < fNtrees; file++) {
01085       Long64_t first = fTreeOffset[file];
01086       LoadTree(first);
01087       Double_t curmin = fTree->GetMinimum(columname);
01088       if (curmin < theMin) {
01089          theMin = curmin;
01090       }
01091    }
01092    return theMin;
01093 }
01094 
01095 //______________________________________________________________________________
01096 Int_t TChain::GetNbranches()
01097 {
01098    // -- Return the number of branches of the current tree.
01099    //
01100    // Warning: May set the current tree!
01101    //
01102 
01103    if (fTree) {
01104       return fTree->GetNbranches();
01105    }
01106    LoadTree(0);
01107    if (fTree) {
01108       return fTree->GetNbranches();
01109    }
01110    return 0;
01111 }
01112 
01113 //______________________________________________________________________________
01114 Long64_t TChain::GetReadEntry() const
01115 {
01116    // -- See TTree::GetReadEntry().
01117 
01118    if (fProofChain && !(fProofChain->TestBit(kProofLite))) {
01119       // Make sure the element list is uptodate
01120       if (!TestBit(kProofUptodate))
01121          Warning("GetBranchStatus", "PROOF proxy not up-to-date:"
01122                                     " run TChain::SetProof(kTRUE, kTRUE) first");
01123       return fProofChain->GetReadEntry();
01124    }
01125    return TTree::GetReadEntry();
01126 }
01127 
01128 //______________________________________________________________________________
01129 Double_t TChain::GetWeight() const
01130 {
01131    // -- Return the chain weight.
01132    //
01133    // By default the weight is the weight of the current tree.
01134    // However, if the weight has been set in TChain::SetWeight()
01135    // with the option "global", then that weight will be returned.
01136    //
01137    // Warning: May set the current tree!
01138    //
01139 
01140    if (TestBit(kGlobalWeight)) {
01141       return fWeight;
01142    } else {
01143       if (fTree) {
01144          return fTree->GetWeight();
01145       }
01146       const_cast<TChain*>(this)->LoadTree(0);
01147       if (fTree) {
01148          return fTree->GetWeight();
01149       }
01150       return 0;
01151    }
01152 }
01153 
01154 //______________________________________________________________________________
01155 Int_t TChain::LoadBaskets(Long64_t /*maxmemory*/)
01156 {
01157    // -- Dummy function.
01158    // It could be implemented and load all baskets of all trees in the chain.
01159    // For the time being use TChain::Merge and TTree::LoadBasket
01160    // on the resulting tree.
01161 
01162    Error("LoadBaskets", "Function not yet implemented for TChain.");
01163    return 0;
01164 }
01165 
01166 //______________________________________________________________________________
01167 Long64_t TChain::LoadTree(Long64_t entry)
01168 {
01169    // -- Find the tree which contains entry, and set it as the current tree.
01170    //
01171    // Returns the entry number in that tree.
01172    //
01173    // The input argument entry is the entry serial number in the whole chain.
01174    //
01175    // Note: This is the only routine which sets the value of fTree to
01176    //       a non-zero pointer.
01177    //
01178 
01179    // We already have been visited while recursively looking
01180    // through the friends tree, let's return.
01181    if (kLoadTree & fFriendLockStatus) {
01182       return 0;
01183    }
01184 
01185    if (!fNtrees) {
01186       // -- The chain is empty.
01187       return 1;
01188    }
01189 
01190    if ((entry < 0) || ((entry > 0) && (entry >= fEntries && entry!=(theBigNumber-1) ))) {
01191       // -- Invalid entry number.
01192       if (fTree) fTree->LoadTree(-1);
01193       fReadEntry = -1;
01194       return -2;
01195    }
01196 
01197    // Find out which tree in the chain contains the passed entry.
01198    Int_t treenum = fTreeNumber;
01199    if ((fTreeNumber == -1) || (entry < fTreeOffset[fTreeNumber]) || (entry >= fTreeOffset[fTreeNumber+1]) || (entry==theBigNumber-1)) {
01200       // -- Entry is *not* in the chain's current tree.
01201       // Do a linear search of the tree offset array.
01202       // FIXME: We could be smarter by starting at the
01203       //        current tree number and going forwards,
01204       //        then wrapping around at the end.
01205       for (treenum = 0; treenum < fNtrees; treenum++) {
01206          if (entry < fTreeOffset[treenum+1]) {
01207             break;
01208          }
01209       }
01210    }
01211 
01212    // Calculate the entry number relative to the found tree.
01213    Long64_t treeReadEntry = entry - fTreeOffset[treenum];
01214    fReadEntry = entry;
01215 
01216    // If entry belongs to the current tree return entry.
01217    if (fTree && treenum == fTreeNumber) {
01218       // First set the entry the tree on its owns friends
01219       // (the friends of the chain will be updated in the
01220       // next loop).
01221       fTree->LoadTree(treeReadEntry);
01222       if (fFriends) {
01223          // The current tree has not changed but some of its friends might.
01224          //
01225          // An alternative would move this code to each of
01226          // the functions calling LoadTree (and to overload a few more).
01227          TIter next(fFriends);
01228          TFriendLock lock(this, kLoadTree);
01229          TFriendElement* fe = 0;
01230          TFriendElement* fetree = 0;
01231          Bool_t needUpdate = kFALSE;
01232          while ((fe = (TFriendElement*) next())) {
01233             TObjLink* lnk = 0;
01234             if (fTree->GetListOfFriends()) {
01235                lnk = fTree->GetListOfFriends()->FirstLink();
01236             }
01237             fetree = 0;
01238             while (lnk) {
01239                TObject* obj = lnk->GetObject();
01240                if (obj->TestBit(TFriendElement::kFromChain) && obj->GetName() && !strcmp(fe->GetName(), obj->GetName())) {
01241                   fetree = (TFriendElement*) obj;
01242                   break;
01243                }
01244                lnk = lnk->Next();
01245             }
01246             TTree* at = fe->GetTree();
01247             if (at->InheritsFrom(TChain::Class())) {
01248                Int_t oldNumber = ((TChain*) at)->GetTreeNumber();
01249                TTree* old = at->GetTree();
01250                TTree* oldintree = fetree ? fetree->GetTree() : 0;
01251                at->LoadTreeFriend(entry, this);
01252                Int_t newNumber = ((TChain*) at)->GetTreeNumber();
01253                if ((oldNumber != newNumber) || (old != at->GetTree()) || (oldintree && (oldintree != at->GetTree()))) {
01254                   // We can not compare just the tree pointers because
01255                   // they could be reused. So we compare the tree
01256                   // number instead.
01257                   needUpdate = kTRUE;
01258                   fTree->RemoveFriend(oldintree);
01259                   fTree->AddFriend(at->GetTree(), fe->GetName())->SetBit(TFriendElement::kFromChain);
01260                }
01261             } else {
01262                // else we assume it is a simple tree If the tree is a
01263                // direct friend of the chain, it should be scanned
01264                // used the chain entry number and NOT the tree entry
01265                // number (treeReadEntry) hence we redo:
01266                at->LoadTreeFriend(entry, this);
01267             }
01268          }
01269          if (needUpdate) {
01270             // Update the branch/leaf addresses and
01271             // thelist of leaves in all TTreeFormula of the TTreePlayer (if any).
01272 
01273             // Set the branch statuses for the newly opened file.
01274             TChainElement *frelement;
01275             TIter fnext(fStatus);
01276             while ((frelement = (TChainElement*) fnext())) {
01277                Int_t status = frelement->GetStatus();
01278                fTree->SetBranchStatus(frelement->GetName(), status);
01279             }
01280 
01281             // Set the branch addresses for the newly opened file.
01282             fnext.Reset();
01283             while ((frelement = (TChainElement*) fnext())) {
01284                void* addr = frelement->GetBaddress();
01285                if (addr) {
01286                   TBranch* br = fTree->GetBranch(frelement->GetName());
01287                   TBranch** pp = frelement->GetBranchPtr();
01288                   if (pp) {
01289                      // FIXME: What if br is zero here?
01290                      *pp = br;
01291                   }
01292                   if (br) {
01293                      // FIXME: We may have to tell the branch it should
01294                      //        not be an owner of the object pointed at.
01295                      br->SetAddress(addr);
01296                      if (TestBit(kAutoDelete)) {
01297                         br->SetAutoDelete(kTRUE);
01298                      }
01299                   }
01300                }
01301             }
01302             if (fPlayer) {
01303                fPlayer->UpdateFormulaLeaves();
01304             }
01305             // Notify user if requested.
01306             if (fNotify) {
01307                fNotify->Notify();
01308             }
01309          }
01310       }
01311       return treeReadEntry;
01312    }
01313 
01314    // If the tree has clones, copy them into the chain
01315    // clone list so we can change their branch addresses
01316    // when necessary.
01317    //
01318    // This is to support the syntax:
01319    //
01320    //      TTree* clone = chain->GetTree()->CloneTree(0);
01321    //
01322    if (fTree && fTree->GetListOfClones()) {
01323       for (TObjLink* lnk = fTree->GetListOfClones()->FirstLink(); lnk; lnk = lnk->Next()) {
01324          TTree* clone = (TTree*) lnk->GetObject();
01325          AddClone(clone);
01326       }
01327    }
01328 
01329    // Delete the current tree and open the new tree.
01330    TTreeCache* tpf = 0;
01331 
01332    // Delete file unless the file owns this chain!
01333    // FIXME: The "unless" case here causes us to leak memory.
01334    if (fFile) {
01335       if (!fDirectory->GetList()->FindObject(this)) {
01336          tpf = (TTreeCache*) fFile->GetCacheRead();
01337          if (tpf) tpf->ResetCache();
01338          fFile->SetCacheRead(0);
01339          if (fCanDeleteRefs) {
01340             fFile->Close("R");
01341          }
01342          delete fFile;
01343          fFile = 0;
01344          // Note: We do *not* own fTree.
01345          fTree = 0;
01346       }
01347    }
01348 
01349    TChainElement* element = (TChainElement*) fFiles->At(treenum);
01350    if (!element) {
01351       if (treeReadEntry) {
01352          return -4;
01353       }
01354       // Last attempt, just in case all trees in the chain have 0 entries.
01355       element = (TChainElement*) fFiles->At(0);
01356       if (!element) {
01357          return -4;
01358       }
01359    }
01360 
01361    // FIXME: We leak memory here, we've just lost the open file
01362    //        if we did not delete it above.
01363    {
01364       TDirectory::TContext ctxt(0);
01365       fFile = TFile::Open(element->GetTitle());
01366       if (fFile) fFile->SetBit(kMustCleanup);
01367    }
01368 
01369    // ----- Begin of modifications by MvL
01370    Int_t returnCode = 0;
01371    if (!fFile || fFile->IsZombie()) {
01372       if (fFile) {
01373          delete fFile;
01374          fFile = 0;
01375       }
01376       // Note: We do *not* own fTree.
01377       fTree = 0;
01378       returnCode = -3;
01379    } else {
01380       // Note: We do *not* own fTree after this, the file does!
01381       fTree = (TTree*) fFile->Get(element->GetName());
01382       if (!fTree) {
01383          // Now that we do not check during the addition, we need to check here!
01384          Error("LoadTree", "Cannot find tree with name %s in file %s", element->GetName(), element->GetTitle());
01385          delete fFile;
01386          fFile = 0;
01387          // We do not return yet so that 'fEntries' can be updated with the
01388          // sum of the entries of all the other trees.
01389          returnCode = -4;
01390       }
01391    }
01392 
01393    fTreeNumber = treenum;
01394    // FIXME: We own fFile, we must be careful giving away a pointer to it!
01395    // FIXME: We may set fDirectory to zero here!
01396    fDirectory = fFile;
01397 
01398    // Reuse cache from previous file (if any).
01399    if (tpf) {
01400       if (fFile) {
01401          tpf->ResetCache();
01402          fFile->SetCacheRead(tpf);
01403          tpf->SetFile(fFile);
01404          // FIXME: fTree may be zero here.
01405          tpf->UpdateBranches(fTree);
01406       } else {
01407          // FIXME: One of the file in the chain is missing
01408          // we have no place to hold the pointer to the
01409          // TTreeCache.
01410          delete tpf;
01411          tpf = 0;
01412          this->SetCacheSize(fCacheSize);
01413       }
01414    } else {
01415       this->SetCacheSize(fCacheSize);
01416    }
01417 
01418    // Check if fTreeOffset has really been set.
01419    Long64_t nentries = 0;
01420    if (fTree) {
01421       nentries = fTree->GetEntries();
01422    }
01423 
01424    if (fTreeOffset[fTreeNumber+1] != (fTreeOffset[fTreeNumber] + nentries)) {
01425       fTreeOffset[fTreeNumber+1] = fTreeOffset[fTreeNumber] + nentries;
01426       fEntries = fTreeOffset[fNtrees];
01427       element->SetNumberEntries(nentries);
01428       // Below we must test >= in case the tree has no entries.
01429       if (entry >= fTreeOffset[fTreeNumber+1]) {
01430          if ((fTreeNumber < (fNtrees - 1)) && (entry < fTreeOffset[fTreeNumber+2])) {
01431             return LoadTree(entry);
01432          } else {
01433             treeReadEntry = fReadEntry = -2;
01434          }
01435       }
01436    }
01437 
01438    if (!fTree) {
01439       // The Error message already issued.  However if we reach here
01440       // we need to make sure that we do not use fTree.
01441       //
01442       // Force a reload of the tree next time.
01443       fTreeNumber = -1;
01444       return returnCode;
01445    }
01446    // ----- End of modifications by MvL
01447 
01448    // Copy the chain's clone list into the new tree's
01449    // clone list so that branch addresses stay synchronized.
01450    if (fClones) {
01451       for (TObjLink* lnk = fClones->FirstLink(); lnk; lnk = lnk->Next()) {
01452          TTree* clone = (TTree*) lnk->GetObject();
01453          ((TChain*) fTree)->TTree::AddClone(clone);
01454       }
01455    }
01456 
01457    // Since some of the friends of this chain might simple trees
01458    // (i.e., not really chains at all), we need to execute this
01459    // before calling LoadTree(entry) on the friends (so that
01460    // they use the correct read entry number).
01461 
01462    // Change the new current tree to the new entry.
01463    fTree->LoadTree(treeReadEntry);
01464 
01465    // Change the chain friends to the new entry.
01466    if (fFriends) {
01467       // An alternative would move this code to each of the function
01468       // calling LoadTree (and to overload a few more).
01469       TIter next(fFriends);
01470       TFriendLock lock(this, kLoadTree);
01471       TFriendElement* fe = 0;
01472       while ((fe = (TFriendElement*) next())) {
01473          TTree* t = fe->GetTree();
01474          if (!t) continue;
01475          if (t->GetTreeIndex()) {
01476             t->GetTreeIndex()->UpdateFormulaLeaves(0);
01477          }
01478          if (t->GetTree() && t->GetTree()->GetTreeIndex()) {
01479             t->GetTree()->GetTreeIndex()->UpdateFormulaLeaves(GetTree());
01480          }
01481          t->LoadTreeFriend(entry, this);
01482          TTree* friend_t = t->GetTree();
01483          if (friend_t) {
01484             fTree->AddFriend(friend_t, fe->GetName())->SetBit(TFriendElement::kFromChain);
01485          }
01486       }
01487    }
01488 
01489    fTree->SetMakeClass(fMakeClass);
01490    fTree->SetMaxVirtualSize(fMaxVirtualSize);
01491 
01492    SetChainOffset(fTreeOffset[fTreeNumber]);
01493 
01494    // Set the branch statuses for the newly opened file.
01495    TIter next(fStatus);
01496    while ((element = (TChainElement*) next())) {
01497       Int_t status = element->GetStatus();
01498       fTree->SetBranchStatus(element->GetName(), status);
01499    }
01500 
01501    // Set the branch addresses for the newly opened file.
01502    next.Reset();
01503    while ((element = (TChainElement*) next())) {
01504       void* addr = element->GetBaddress();
01505       if (addr) {
01506          TBranch* br = fTree->GetBranch(element->GetName());
01507          TBranch** pp = element->GetBranchPtr();
01508          if (pp) {
01509             // FIXME: What if br is zero here?
01510             *pp = br;
01511          }
01512          if (br) {
01513             // FIXME: We may have to tell the branch it should
01514             //        not be an owner of the object pointed at.
01515             br->SetAddress(addr);
01516             if (TestBit(kAutoDelete)) {
01517                br->SetAutoDelete(kTRUE);
01518             }
01519          }
01520       }
01521    }
01522 
01523    // Update the addresses of the chain's cloned trees, if any.
01524    if (fClones) {
01525       for (TObjLink* lnk = fClones->FirstLink(); lnk; lnk = lnk->Next()) {
01526          TTree* clone = (TTree*) lnk->GetObject();
01527          CopyAddresses(clone);
01528       }
01529    }
01530 
01531    // Update list of leaves in all TTreeFormula's of the TTreePlayer (if any).
01532    if (fPlayer) {
01533       fPlayer->UpdateFormulaLeaves();
01534    }
01535 
01536    // Notify user we have switched trees if requested.
01537    if (fNotify) {
01538       fNotify->Notify();
01539    }
01540 
01541    // Return the new local entry number.
01542    return treeReadEntry;
01543 }
01544 
01545 //______________________________________________________________________________
01546 void TChain::Lookup(Bool_t force)
01547 {
01548    // Check / locate the files in the chain.
01549    // By default only the files not yet looked up are checked.
01550    // Use force = kTRUE to check / re-check every file.
01551 
01552    TIter next(fFiles);
01553    TChainElement* element = 0;
01554    Int_t nelements = fFiles->GetEntries();
01555    printf("\n");
01556    printf("TChain::Lookup - Looking up %d files .... \n", nelements);
01557    Int_t nlook = 0;
01558    TFileStager *stg = 0;
01559    while ((element = (TChainElement*) next())) {
01560       // Do not do it more than needed
01561       if (element->HasBeenLookedUp() && !force) continue;
01562       // Count
01563       nlook++;
01564       // Get the Url
01565       TUrl elemurl(element->GetTitle(), kTRUE);
01566       // Save current options and anchor
01567       TString anchor = elemurl.GetAnchor();
01568       TString options = elemurl.GetOptions();
01569       // Reset options and anchor
01570       elemurl.SetOptions("");
01571       elemurl.SetAnchor("");
01572       // Locate the file
01573       TString eurl(elemurl.GetUrl());
01574       if (!stg || !stg->Matches(eurl)) {
01575          SafeDelete(stg);
01576          {
01577             TDirectory::TContext ctxt(0);
01578             stg = TFileStager::Open(eurl);
01579          }
01580          if (!stg) {
01581             Error("Lookup", "TFileStager instance cannot be instantiated");
01582             break;
01583          }
01584       }
01585       Int_t n1 = (nelements > 100) ? (Int_t) nelements / 100 : 1;
01586       if (stg->Locate(eurl.Data(), eurl) == 0) {
01587          if (nlook > 0 && !(nlook % n1)) {
01588             printf("Lookup | %3d %% finished\r", 100 * nlook / nelements);
01589             fflush(stdout);
01590          }
01591          // Get the effective end-point Url
01592          elemurl.SetUrl(eurl);
01593          // Restore original options and anchor, if any
01594          elemurl.SetOptions(options);
01595          elemurl.SetAnchor(anchor);
01596          // Save it into the element
01597          element->SetTitle(elemurl.GetUrl());
01598          // Remember
01599          element->SetLookedUp();
01600       } else {
01601          // Failure: remove
01602          fFiles->Remove(element);
01603          if (gSystem->AccessPathName(eurl))
01604             Error("Lookup", "file %s does not exist\n", eurl.Data());
01605          else
01606             Error("Lookup", "file %s cannot be read\n", eurl.Data());
01607       }
01608    }
01609    if (nelements > 0)
01610       printf("Lookup | %3d %% finished\n", 100 * nlook / nelements);
01611    else
01612       printf("\n");
01613    fflush(stdout);
01614    SafeDelete(stg);
01615 }
01616 
01617 //______________________________________________________________________________
01618 void TChain::Loop(Option_t* option, Long64_t nentries, Long64_t firstentry)
01619 {
01620    // -- Loop on nentries of this chain starting at firstentry.  (NOT IMPLEMENTED)
01621 
01622    Error("Loop", "Function not yet implemented");
01623 
01624    if (option || nentries || firstentry) { }  // keep warnings away
01625 
01626 #if 0
01627    if (LoadTree(firstentry) < 0) return;
01628 
01629    if (firstentry < 0) firstentry = 0;
01630    Long64_t lastentry = firstentry + nentries -1;
01631    if (lastentry > fEntries-1) {
01632       lastentry = fEntries -1;
01633    }
01634 
01635    GetPlayer();
01636    GetSelector();
01637    fSelector->Start(option);
01638 
01639    Long64_t entry = firstentry;
01640    Int_t tree,e0,en;
01641    for (tree=0;tree<fNtrees;tree++) {
01642       e0 = fTreeOffset[tree];
01643       en = fTreeOffset[tree+1] - 1;
01644       if (en > lastentry) en = lastentry;
01645       if (entry > en) continue;
01646 
01647       LoadTree(entry);
01648       fSelector->BeginFile();
01649 
01650       while (entry <= en) {
01651          fSelector->Execute(fTree, entry - e0);
01652          entry++;
01653       }
01654       fSelector->EndFile();
01655    }
01656 
01657    fSelector->Finish(option);
01658 #endif
01659 }
01660 
01661 //______________________________________________________________________________
01662 void TChain::ls(Option_t* option) const
01663 {
01664    // -- List the chain.
01665    TIter next(fFiles);
01666    TChainElement* file = 0;
01667    while ((file = (TChainElement*)next())) {
01668       file->ls(option);
01669    }
01670 }
01671 
01672 //______________________________________________________________________________
01673 Long64_t TChain::Merge(const char* name, Option_t* option)
01674 {
01675    // Merge all the entries in the chain into a new tree in a new file.
01676    //
01677    // See important note in the following function Merge().
01678    //
01679    // If the chain is expecting the input tree inside a directory,
01680    // this directory is NOT created by this routine.
01681    //
01682    // So in a case where we have:
01683    //
01684    //      TChain ch("mydir/mytree");
01685    //      ch.Merge("newfile.root");
01686    //
01687    // The resulting file will have not subdirectory. To recreate
01688    // the directory structure do:
01689    //
01690    //      TFile* file = TFile::Open("newfile.root", "RECREATE");
01691    //      file->mkdir("mydir")->cd();
01692    //      ch.Merge(file);
01693    //
01694 
01695    TFile *file = TFile::Open(name, "recreate", "chain files", 1);
01696    return Merge(file, 0, option);
01697 }
01698 
01699 //______________________________________________________________________________
01700 Long64_t TChain::Merge(TCollection* /* list */, Option_t* /* option */ )
01701 {
01702    // Merge all chains in the collection.  (NOT IMPLEMENTED)
01703 
01704    Error("Merge", "not implemented");
01705    return -1;
01706 }
01707 
01708 //______________________________________________________________________________
01709 Long64_t TChain::Merge(TFile* file, Int_t basketsize, Option_t* option)
01710 {
01711    // Merge all the entries in the chain into a new tree in the current file.
01712    //
01713    // Note: The "file" parameter is *not* the file where the new
01714    //       tree will be inserted.  The new tree is inserted into
01715    //       gDirectory, which is usually the most recently opened
01716    //       file, or the directory most recently cd()'d to.
01717    //
01718    // If option = "C" is given, the compression level for all branches
01719    // in the new Tree is set to the file compression level.  By default,
01720    // the compression level of all branches is the original compression
01721    // level in the old trees.
01722    //
01723    // If basketsize > 1000, the basket size for all branches of the
01724    // new tree will be set to basketsize.
01725    //
01726    // Example using the file generated in $ROOTSYS/test/Event
01727    // merge two copies of Event.root
01728    //
01729    //        gSystem.Load("libEvent");
01730    //        TChain ch("T");
01731    //        ch.Add("Event1.root");
01732    //        ch.Add("Event2.root");
01733    //        ch.Merge("all.root");
01734    //
01735    // If the chain is expecting the input tree inside a directory,
01736    // this directory is NOT created by this routine.
01737    //
01738    // So if you do:
01739    //
01740    //      TChain ch("mydir/mytree");
01741    //      ch.Merge("newfile.root");
01742    //
01743    // The resulting file will not have subdirectories.  In order to
01744    // preserve the directory structure do the following instead:
01745    //
01746    //      TFile* file = TFile::Open("newfile.root", "RECREATE");
01747    //      file->mkdir("mydir")->cd();
01748    //      ch.Merge(file);
01749    //
01750    // If 'option' contains the word 'fast' the merge will be done without
01751    // unzipping or unstreaming the baskets (i.e., a direct copy of the raw
01752    // bytes on disk).
01753    //
01754    // When 'fast' is specified, 'option' can also contains a
01755    // sorting order for the baskets in the output file.
01756    //
01757    // There is currently 3 supported sorting order:
01758    //    SortBasketsByOffset (the default)
01759    //    SortBasketsByBranch
01760    //    SortBasketsByEntry
01761    //
01762    // When using SortBasketsByOffset the baskets are written in
01763    // the output file in the same order as in the original file
01764    // (i.e. the basket are sorted on their offset in the original
01765    // file; Usually this also means that the baskets are sorted
01766    // on the index/number of the _last_ entry they contain)
01767    //
01768    // When using SortBasketsByBranch all the baskets of each
01769    // individual branches are stored contiguously.  This tends to
01770    // optimize reading speed when reading a small number (1->5) of
01771    // branches, since all their baskets will be clustered together
01772    // instead of being spread across the file.  However it might
01773    // decrease the performance when reading more branches (or the full
01774    // entry).
01775    //
01776    // When using SortBasketsByEntry the baskets with the lowest
01777    // starting entry are written first.  (i.e. the baskets are
01778    // sorted on the index/number of the first entry they contain).
01779    // This means that on the file the baskets will be in the order
01780    // in which they will be needed when reading the whole tree
01781    // sequentially.
01782    //
01783    // IMPORTANT Note 1: AUTOMATIC FILE OVERFLOW
01784    // -----------------------------------------
01785    // When merging many files, it may happen that the resulting file
01786    // reaches a size > TTree::fgMaxTreeSize (default = 1.9 GBytes).
01787    // In this case the current file is automatically closed and a new
01788    // file started.  If the name of the merged file was "merged.root",
01789    // the subsequent files will be named "merged_1.root", "merged_2.root",
01790    // etc.  fgMaxTreeSize may be modified via the static function
01791    // TTree::SetMaxTreeSize.
01792    // When in fast mode, the check and switch is only done in between each
01793    // input file.
01794    //
01795    // IMPORTANT Note 2: The output file is automatically closed and deleted.
01796    // ----------------------------------------------------------------------
01797    // This is required because in general the automatic file overflow described
01798    // above may happen during the merge.
01799    // If only the current file is produced (the file passed as first argument),
01800    // one can instruct Merge to not close and delete the file by specifying
01801    // the option "keep".
01802    //
01803    // The function returns the total number of files produced.
01804    // To check that all files have been merged use something like:
01805    //    if (newchain->GetEntries()!=oldchain->GetEntries()) {
01806    //      ... not all the file have been copied ...
01807    //    }
01808 
01809    // We must have been passed a file, we will use it
01810    // later to reset the compression level of the branches.
01811    if (!file) {
01812       // FIXME: We need an error message here.
01813       return 0;
01814    }
01815 
01816    // Options
01817    Bool_t fastClone = kFALSE;
01818    TString opt = option;
01819    opt.ToLower();
01820    if (opt.Contains("fast")) {
01821       fastClone = kTRUE;
01822    }
01823 
01824    // The chain tree must have a list of branches
01825    // because we may try to change their basket
01826    // size later.
01827    TObjArray* lbranches = GetListOfBranches();
01828    if (!lbranches) {
01829       // FIXME: We need an error message here.
01830       return 0;
01831    }
01832 
01833    // The chain must have a current tree because
01834    // that is the one we will clone.
01835    if (!fTree) {
01836       // -- LoadTree() has not yet been called, no current tree.
01837       // FIXME: We need an error message here.
01838       return 0;
01839    }
01840 
01841    // Copy the chain's current tree without
01842    // copying any entries, we will do that later.
01843    TTree* newTree = CloneTree(0);
01844    if (!newTree) {
01845       // FIXME: We need an error message here.
01846       return 0;
01847    }
01848 
01849    // Strip out the (potential) directory name.
01850    // FIXME: The merged chain may or may not have the
01851    //        same name as the original chain.  This is
01852    //        bad because the chain name determines the
01853    //        names of the trees in the chain by default.
01854    newTree->SetName(gSystem->BaseName(GetName()));
01855 
01856    // FIXME: Why do we do this?
01857    newTree->SetAutoSave(2000000000);
01858 
01859    // Circularity is incompatible with merging, it may
01860    // force us to throw away entries, which is not what
01861    // we are supposed to do.
01862    newTree->SetCircular(0);
01863 
01864    // Reset the compression level of the branches.
01865    if (opt.Contains("c")) {
01866       TBranch* branch = 0;
01867       TIter nextb(newTree->GetListOfBranches());
01868       while ((branch = (TBranch*) nextb())) {
01869          branch->SetCompressionLevel(file->GetCompressionLevel());
01870       }
01871    }
01872 
01873    // Reset the basket size of the branches.
01874    if (basketsize > 1000) {
01875       TBranch* branch = 0;
01876       TIter nextb(newTree->GetListOfBranches());
01877       while ((branch = (TBranch*) nextb())) {
01878          branch->SetBasketSize(basketsize);
01879       }
01880    }
01881 
01882    // Copy the entries.
01883    if (fastClone) {
01884       if ( newTree->CopyEntries( this, -1, option ) < 0 ) {
01885          // There was a problem!
01886          Error("Merge", "TTree has not been cloned\n");
01887       }
01888    } else {
01889       newTree->CopyEntries( this, -1, option );
01890    }
01891 
01892    // Write the new tree header.
01893    newTree->Write();
01894 
01895    // Get our return value.
01896    Int_t nfiles = newTree->GetFileNumber() + 1;
01897 
01898    // Close and delete the current file of the new tree.
01899    if (!opt.Contains("keep")) {
01900       // Delete the currentFile and the TTree object.
01901       delete newTree->GetCurrentFile();
01902    }
01903    return nfiles;
01904 }
01905 
01906 //______________________________________________________________________________
01907 void TChain::Print(Option_t *option) const
01908 {
01909    // -- Print the header information of each tree in the chain.
01910    // See TTree::Print for a list of options.
01911 
01912    TIter next(fFiles);
01913    TChainElement *element;
01914    while ((element = (TChainElement*)next())) {
01915       TFile *file = TFile::Open(element->GetTitle());
01916       if (file && !file->IsZombie()) {
01917          TTree *tree = (TTree*)file->Get(element->GetName());
01918          if (tree) tree->Print(option);
01919       }
01920       delete file;
01921    }
01922 }
01923 
01924 //______________________________________________________________________________
01925 Long64_t TChain::Process(const char *filename, Option_t *option, Long64_t nentries, Long64_t firstentry)
01926 {
01927    // Process all entries in this chain, calling functions in filename.
01928    // The return value is -1 in case of error and TSelector::GetStatus() in
01929    // in case of success.
01930    // See TTree::Process.
01931 
01932    if (fProofChain) {
01933       // Make sure the element list is uptodate
01934       if (!TestBit(kProofUptodate))
01935          SetProof(kTRUE, kTRUE);
01936       fProofChain->SetEventList(fEventList);
01937       fProofChain->SetEntryList(fEntryList);
01938       return fProofChain->Process(filename, option, nentries, firstentry);
01939    }
01940 
01941    if (LoadTree(firstentry) < 0) {
01942       return 0;
01943    }
01944    return TTree::Process(filename, option, nentries, firstentry);
01945 }
01946 
01947 //______________________________________________________________________________
01948 Long64_t TChain::Process(TSelector* selector, Option_t* option, Long64_t nentries, Long64_t firstentry)
01949 {
01950    // Process this chain executing the code in selector.
01951    // The return value is -1 in case of error and TSelector::GetStatus() in
01952    // in case of success.
01953 
01954    if (fProofChain) {
01955       // Make sure the element list is uptodate
01956       if (!TestBit(kProofUptodate))
01957          SetProof(kTRUE, kTRUE);
01958       fProofChain->SetEventList(fEventList);
01959       fProofChain->SetEntryList(fEntryList);
01960       return fProofChain->Process(selector, option, nentries, firstentry);
01961    }
01962 
01963    return TTree::Process(selector, option, nentries, firstentry);
01964 }
01965 
01966 //______________________________________________________________________________
01967 void TChain::RecursiveRemove(TObject *obj)
01968 {
01969    // Make sure that obj (which is being deleted or will soon be) is no
01970    // longer referenced by this TTree.
01971 
01972    if (fFile == obj) {
01973       fFile = 0;
01974       fDirectory = 0;
01975       fTree = 0;
01976    }
01977    if (fDirectory == obj) {
01978       fDirectory = 0;
01979       fTree = 0;
01980    }
01981    if (fTree == obj) {
01982       fTree = 0;
01983    }
01984 }
01985 
01986 //______________________________________________________________________________
01987 void TChain::Reset(Option_t*)
01988 {
01989    // -- Resets the state of this chain.
01990 
01991    delete fFile;
01992    fFile = 0;
01993    fNtrees         = 0;
01994    fTreeNumber     = -1;
01995    fTree           = 0;
01996    fFile           = 0;
01997    fFiles->Delete();
01998    fStatus->Delete();
01999    fTreeOffset[0]  = 0;
02000    TChainElement* element = new TChainElement("*", "");
02001    fStatus->Add(element);
02002    fDirectory = 0;
02003 
02004    TTree::Reset();
02005 }
02006 
02007 //_______________________________________________________________________
02008 Long64_t TChain::Scan(const char* varexp, const char* selection, Option_t* option, Long64_t nentries, Long64_t firstentry)
02009 {
02010    // -- Loop on tree and print entries passing selection.
02011    // If varexp is 0 (or "") then print only first 8 columns.
02012    // If varexp = "*" print all columns.
02013    // Otherwise a columns selection can be made using "var1:var2:var3".
02014    // See TTreePlayer::Scan for more information.
02015 
02016    if (LoadTree(firstentry) < 0) {
02017       return 0;
02018    }
02019    return TTree::Scan(varexp, selection, option, nentries, firstentry);
02020 }
02021 
02022 //_______________________________________________________________________
02023 void TChain::SetAutoDelete(Bool_t autodelete)
02024 {
02025    // -- Set the global branch kAutoDelete bit.
02026    //
02027    //  When LoadTree loads a new Tree, the branches for which
02028    //  the address is set will have the option AutoDelete set
02029    //  For more details on AutoDelete, see TBranch::SetAutoDelete.
02030 
02031    if (autodelete) {
02032       SetBit(kAutoDelete, 1);
02033    } else {
02034       SetBit(kAutoDelete, 0);
02035    }
02036 }
02037 
02038 //______________________________________________________________________________
02039 void TChain::ResetBranchAddress(TBranch *branch)
02040 {
02041    // -- Reset the addresses of the branch.
02042 
02043    TChainElement* element = (TChainElement*) fStatus->FindObject(branch->GetName());
02044    if (element) {
02045       element->SetBaddress(0);
02046    }   
02047    if (fTree) {
02048       fTree->ResetBranchAddress(branch);
02049    }
02050 }
02051 
02052 //______________________________________________________________________________
02053 void TChain::ResetBranchAddresses()
02054 {
02055    // Reset the addresses of the branches.
02056 
02057    TIter next(fStatus);
02058    TChainElement* element = 0;
02059    while ((element = (TChainElement*) next())) {
02060       element->SetBaddress(0);
02061    }
02062    if (fTree) {
02063       fTree->ResetBranchAddresses();
02064    }
02065 }
02066 
02067 //_______________________________________________________________________
02068 Int_t TChain::SetBranchAddress(const char *bname, void* add, TBranch** ptr)
02069 {
02070    // Set branch address.
02071    //
02072    //      bname is the name of a branch.
02073    //      add is the address of the branch.
02074    //
02075    //    Note: See the comments in TBranchElement::SetAddress() for a more
02076    //          detailed discussion of the meaning of the add parameter.
02077    //
02078    // IMPORTANT REMARK:
02079    // In case TChain::SetBranchStatus is called, it must be called
02080    // BEFORE calling this function.
02081    //
02082    // See TTree::CheckBranchAddressType for the semantic of the return value.
02083 
02084    Int_t res = kNoCheck;
02085 
02086    // Check if bname is already in the status list.
02087    // If not, create a TChainElement object and set its address.
02088    TChainElement* element = (TChainElement*) fStatus->FindObject(bname);
02089    if (!element) {
02090       element = new TChainElement(bname, "");
02091       fStatus->Add(element);
02092    }
02093    element->SetBaddress(add);
02094    element->SetBranchPtr(ptr);
02095    // Also set address in current tree.
02096    // FIXME: What about the chain clones?
02097    if (fTreeNumber >= 0) {
02098       TBranch* branch = fTree->GetBranch(bname);
02099       if (ptr) {
02100          *ptr = branch;
02101       }
02102       if (branch) {
02103          res = CheckBranchAddressType(branch, TClass::GetClass(element->GetBaddressClassName()), (EDataType) element->GetBaddressType(), element->GetBaddressIsPtr());
02104          if (fClones) {
02105             void* oldAdd = branch->GetAddress();
02106             for (TObjLink* lnk = fClones->FirstLink(); lnk; lnk = lnk->Next()) {
02107                TTree* clone = (TTree*) lnk->GetObject();
02108                TBranch* cloneBr = clone->GetBranch(bname);
02109                if (cloneBr && (cloneBr->GetAddress() == oldAdd)) {
02110                   // the clone's branch is still pointing to us
02111                   cloneBr->SetAddress(add);
02112                }
02113             }
02114          }
02115          branch->SetAddress(add);
02116       }
02117    } else {
02118       if (ptr) {
02119          *ptr = 0;
02120       }
02121    }
02122    return res;
02123 }
02124 
02125 //_______________________________________________________________________
02126 Int_t TChain::SetBranchAddress(const char* bname, void* add, TClass* realClass, EDataType datatype, Bool_t isptr)
02127 {
02128    // Check if bname is already in the status list, and if not, create a TChainElement object and set its address.
02129    // See TTree::CheckBranchAddressType for the semantic of the return value.
02130    //
02131    //    Note: See the comments in TBranchElement::SetAddress() for a more
02132    //          detailed discussion of the meaning of the add parameter.
02133    //
02134    return SetBranchAddress(bname, add, 0, realClass, datatype, isptr);
02135 }
02136 
02137 //_______________________________________________________________________
02138 Int_t TChain::SetBranchAddress(const char* bname, void* add, TBranch** ptr, TClass* realClass, EDataType datatype, Bool_t isptr)
02139 {
02140    // Check if bname is already in the status list, and if not, create a TChainElement object and set its address.
02141    // See TTree::CheckBranchAddressType for the semantic of the return value.
02142    //
02143    //    Note: See the comments in TBranchElement::SetAddress() for a more
02144    //          detailed discussion of the meaning of the add parameter.
02145    //
02146 
02147    TChainElement* element = (TChainElement*) fStatus->FindObject(bname);
02148    if (!element) {
02149       element = new TChainElement(bname, "");
02150       fStatus->Add(element);
02151    }
02152    if (realClass) {
02153       element->SetBaddressClassName(realClass->GetName());
02154    }
02155    element->SetBaddressType((UInt_t) datatype);
02156    element->SetBaddressIsPtr(isptr);
02157    element->SetBranchPtr(ptr);
02158    return SetBranchAddress(bname, add, ptr);
02159 }
02160 
02161 //_______________________________________________________________________
02162 void TChain::SetBranchStatus(const char* bname, Bool_t status, UInt_t* found)
02163 {
02164    // -- Set branch status to Process or DoNotProcess
02165    //
02166    //      bname is the name of a branch. if bname="*", apply to all branches.
02167    //      status = 1  branch will be processed
02168    //             = 0  branch will not be processed
02169    //  See IMPORTANT REMARKS in TTree::SetBranchStatus and TChain::SetBranchAddress
02170    //
02171    //  If found is not 0, the number of branch(es) found matching the regular
02172    //  expression is returned in *found AND the error message 'unknown branch'
02173    //  is suppressed.
02174 
02175    // FIXME: We never explicitly set found to zero!
02176 
02177    // Check if bname is already in the status list,
02178    // if not create a TChainElement object and set its status.
02179    TChainElement* element = (TChainElement*) fStatus->FindObject(bname);
02180    if (element) {
02181       fStatus->Remove(element);
02182    } else {
02183       element = new TChainElement(bname, "");
02184    }
02185    fStatus->Add(element);
02186    element->SetStatus(status);
02187    // Also set status in current tree.
02188    if (fTreeNumber >= 0) {
02189       fTree->SetBranchStatus(bname, status, found);
02190    } else if (found) {
02191       *found = 1;
02192    }
02193 }
02194 
02195 //______________________________________________________________________________
02196 void TChain::SetDirectory(TDirectory* dir)
02197 {
02198    // Remove reference to this chain from current directory and add
02199    // reference to new directory dir. dir can be 0 in which case the chain
02200    // does not belong to any directory.
02201 
02202    if (fDirectory == dir) return;
02203    if (fDirectory) fDirectory->Remove(this);
02204    fDirectory = dir;
02205    if (fDirectory) {
02206       fDirectory->Append(this);
02207       fFile = fDirectory->GetFile();
02208    } else {
02209       fFile = 0;
02210    }
02211 }
02212 
02213 //_______________________________________________________________________
02214 void TChain::SetEntryList(TEntryList *elist, Option_t *opt)
02215 {
02216    //Set the input entry list (processing the entries of the chain will then be
02217    //limited to the entries in the list)
02218    //This function finds correspondance between the sub-lists of the TEntryList
02219    //and the trees of the TChain
02220    //By default (opt=""), both the file names of the chain elements and
02221    //the file names of the TEntryList sublists are expanded to full path name.
02222    //If opt = "ne", the file names are taken as they are and not expanded
02223 
02224    if (fEntryList){
02225       //check, if the chain is the owner of the previous entry list
02226       //(it happens, if the previous entry list was created from a user-defined
02227       //TEventList in SetEventList() function)
02228       if (fEntryList->TestBit(kCanDelete)) {
02229          TEntryList *tmp = fEntryList;
02230          fEntryList = 0; // Avoid problem with RecursiveRemove.
02231          delete tmp;
02232       } else {
02233          fEntryList = 0;
02234       }
02235    }
02236    if (!elist){
02237       fEntryList = 0;
02238       fEventList = 0;
02239       return;
02240    }
02241    if (!elist->TestBit(kCanDelete)){
02242       //this is a direct call to SetEntryList, not via SetEventList
02243       fEventList = 0;
02244    }
02245    if (elist->GetN() == 0){
02246       fEntryList = elist;
02247       return;
02248    }
02249    if (fProofChain){
02250       //for processing on proof, event list and entry list can't be
02251       //set at the same time.
02252       fEventList = 0;
02253       fEntryList = elist;
02254       return;
02255    }
02256 
02257    Int_t ne = fFiles->GetEntries();
02258    Int_t listfound=0;
02259    TString treename, filename;
02260 
02261    TEntryList *templist = 0;
02262    for (Int_t ie = 0; ie<ne; ie++){
02263       treename = gSystem->BaseName( ((TChainElement*)fFiles->UncheckedAt(ie))->GetName() );
02264       filename = ((TChainElement*)fFiles->UncheckedAt(ie))->GetTitle();
02265       templist = elist->GetEntryList(treename.Data(), filename.Data(), opt);
02266       if (templist) {
02267          listfound++;
02268          templist->SetTreeNumber(ie);
02269       }
02270    }
02271 
02272    if (listfound == 0){
02273       Error("SetEntryList", "No list found for the trees in this chain");
02274       fEntryList = 0;
02275       return;
02276    }
02277    fEntryList = elist;
02278    TList *elists = elist->GetLists();
02279    Bool_t shift = kFALSE;
02280    TIter next(elists);
02281 
02282    //check, if there are sub-lists in the entry list, that don't
02283    //correspond to any trees in the chain
02284    while((templist = (TEntryList*)next())){
02285       if (templist->GetTreeNumber() < 0){
02286          shift = kTRUE;
02287          break;
02288       }
02289    }
02290    fEntryList->SetShift(shift);
02291 
02292 }
02293 
02294 //_______________________________________________________________________
02295 void TChain::SetEntryListFile(const char *filename, Option_t * /*opt*/)
02296 {
02297 // Set the input entry list (processing the entries of the chain will then be
02298 // limited to the entries in the list). This function creates a special kind
02299 // of entry list (TEntryListFromFile object) that loads lists, corresponding
02300 // to the chain elements, one by one, so that only one list is in memory at a time.
02301 //
02302 // If there is an error opening one of the files, this file is skipped and the
02303 // next file is loaded
02304 //
02305 // File naming convention:
02306 // - by default, filename_elist.root is used, where filename is the
02307 //   name of the chain element
02308 // - xxx$xxx.root - $ sign is replaced by the name of the chain element
02309 // If the list name is not specified (by passing filename_elist.root/listname to
02310 // the TChain::SetEntryList() function, the first object of class TEntryList
02311 // in the file is taken.
02312 //
02313 // It is assumed, that there are as many list files, as there are elements in
02314 // the chain and they are in the same order
02315 
02316 
02317    if (fEntryList){
02318       //check, if the chain is the owner of the previous entry list
02319       //(it happens, if the previous entry list was created from a user-defined
02320       //TEventList in SetEventList() function)
02321       if (fEntryList->TestBit(kCanDelete)) {
02322          TEntryList *tmp = fEntryList;
02323          fEntryList = 0; // Avoid problem with RecursiveRemove.
02324          delete tmp;
02325       } else {
02326          fEntryList = 0;
02327       }
02328    }
02329 
02330    fEventList = 0;
02331 
02332    TString basename(filename);
02333 
02334    Int_t dotslashpos = basename.Index(".root/");
02335    TString behind_dot_root = "";
02336    if (dotslashpos>=0) {
02337       // Copy the list name specification
02338       behind_dot_root = basename(dotslashpos+6,basename.Length()-dotslashpos+6);
02339       // and remove it from basename
02340       basename.Remove(dotslashpos+5);
02341    }
02342    fEntryList = new TEntryListFromFile(basename.Data(), behind_dot_root.Data(), fNtrees);
02343    fEntryList->SetBit(kCanDelete, kTRUE);
02344    fEntryList->SetDirectory(0);
02345    ((TEntryListFromFile*)fEntryList)->SetFileNames(fFiles);
02346 }
02347 
02348 
02349 //_______________________________________________________________________
02350 void TChain::SetEventList(TEventList *evlist)
02351 {
02352 //This function transfroms the given TEventList into a TEntryList
02353 //
02354 //NOTE, that this function loads all tree headers, because the entry numbers
02355 //in the TEventList are global and have to be recomputed, taking into account
02356 //the number of entries in each tree.
02357 //
02358 //The new TEntryList is owned by the TChain and gets deleted when the chain
02359 //is deleted. This TEntryList is returned by GetEntryList() function, and after
02360 //GetEntryList() function is called, the TEntryList is not owned by the chain
02361 //any more and will not be deleted with it.
02362 
02363    fEventList = evlist;
02364    if (fEntryList) {
02365       if (fEntryList->TestBit(kCanDelete)) {
02366          TEntryList *tmp = fEntryList;
02367          fEntryList = 0; // Avoid problem with RecursiveRemove.
02368          delete tmp;
02369       } else {
02370          fEntryList = 0;
02371       }
02372    }
02373 
02374    if (!evlist) {
02375       fEntryList = 0;
02376       fEventList = 0;
02377       return;
02378    }
02379 
02380    if(fProofChain) {
02381       //on proof, fEventList and fEntryList shouldn't be set at the same time
02382       if (fEntryList){
02383          //check, if the chain is the owner of the previous entry list
02384          //(it happens, if the previous entry list was created from a user-defined
02385          //TEventList in SetEventList() function)
02386          if (fEntryList->TestBit(kCanDelete)){
02387             TEntryList *tmp = fEntryList;
02388             fEntryList = 0; // Avoid problem with RecursiveRemove.
02389             delete tmp;
02390          } else {
02391             fEntryList = 0;
02392          }
02393       }
02394       return;
02395    }
02396 
02397    char enlistname[100];
02398    snprintf(enlistname,100, "%s_%s", evlist->GetName(), "entrylist");
02399    TEntryList *enlist = new TEntryList(enlistname, evlist->GetTitle());
02400    enlist->SetDirectory(0);
02401 
02402    Int_t nsel = evlist->GetN();
02403    Long64_t globalentry, localentry;
02404    const char *treename;
02405    const char *filename;
02406    if (fTreeOffset[fNtrees-1]==theBigNumber){
02407       //Load all the tree headers if the tree offsets are not known
02408       //It is assumed here, that loading the last tree will load all
02409       //previous ones
02410       printf("loading trees\n");
02411       (const_cast<TChain*>(this))->LoadTree(evlist->GetEntry(evlist->GetN()-1));
02412    }
02413    for (Int_t i=0; i<nsel; i++){
02414       globalentry = evlist->GetEntry(i);
02415       //add some protection from globalentry<0 here
02416       Int_t treenum = 0;
02417       while (globalentry>=fTreeOffset[treenum])
02418          treenum++;
02419       treenum--;
02420       localentry = globalentry - fTreeOffset[treenum];
02421       // printf("globalentry=%lld, treeoffset=%lld, localentry=%lld\n", globalentry, fTreeOffset[treenum], localentry);
02422       treename = ((TNamed*)fFiles->At(treenum))->GetName();
02423       filename = ((TNamed*)fFiles->At(treenum))->GetTitle();
02424       //printf("entering for tree %s %s\n", treename, filename);
02425       enlist->SetTree(treename, filename);
02426       enlist->Enter(localentry);
02427    }
02428    enlist->SetBit(kCanDelete, kTRUE);
02429    enlist->SetReapplyCut(evlist->GetReapplyCut());
02430    SetEntryList(enlist);
02431 }
02432 
02433 //_______________________________________________________________________
02434 void TChain::SetPacketSize(Int_t size)
02435 {
02436    // -- Set number of entries per packet for parallel root.
02437 
02438    fPacketSize = size;
02439    TIter next(fFiles);
02440    TChainElement *element;
02441    while ((element = (TChainElement*)next())) {
02442       element->SetPacketSize(size);
02443    }
02444 }
02445 
02446 //______________________________________________________________________________
02447 void TChain::SetProof(Bool_t on, Bool_t refresh, Bool_t gettreeheader)
02448 {
02449    // Enable/Disable PROOF processing on the current default Proof (gProof).
02450    //
02451    // "Draw" and "Processed" commands will be handled by PROOF.
02452    // The refresh and gettreeheader are meaningfull only if on == kTRUE.
02453    // If refresh is kTRUE the underlying fProofChain (chain proxy) is always
02454    // rebuilt (even if already existing).
02455    // If gettreeheader is kTRUE the header of the tree will be read from the
02456    // PROOF cluster: this is only needed for browsing and should be used with
02457    // care because it may take a long time to execute.
02458 
02459    if (!on) {
02460       // Disable
02461       SafeDelete(fProofChain);
02462       // Reset related bit
02463       ResetBit(kProofUptodate);
02464    } else {
02465       if (fProofChain && !refresh &&
02466          (!gettreeheader || (gettreeheader && fProofChain->GetTree()))) {
02467          return;
02468       }
02469       SafeDelete(fProofChain);
02470       ResetBit(kProofUptodate);
02471 
02472       // Make instance of TChainProof via the plugin manager
02473       TPluginHandler *h;
02474       if ((h = gROOT->GetPluginManager()->FindHandler("TChain", "proof"))) {
02475          if (h->LoadPlugin() == -1)
02476             return;
02477          if (!(fProofChain = reinterpret_cast<TChain *>(h->ExecPlugin(2, this, gettreeheader))))
02478             Error("SetProof", "creation of TProofChain failed");
02479          // Set related bits
02480          SetBit(kProofUptodate);
02481       }
02482    }
02483 }
02484 
02485 //______________________________________________________________________________
02486 void TChain::SetWeight(Double_t w, Option_t* option)
02487 {
02488    // -- Set chain weight.
02489    //
02490    // The weight is used by TTree::Draw to automatically weight each
02491    // selected entry in the resulting histogram.
02492    //  For example the equivalent of
02493    //     chain.Draw("x","w")
02494    //  is
02495    //     chain.SetWeight(w,"global");
02496    //     chain.Draw("x");
02497    //
02498    //  By default the weight used will be the weight
02499    //  of each Tree in the TChain. However, one can force the individual
02500    //  weights to be ignored by specifying the option "global".
02501    //  In this case, the TChain global weight will be used for all Trees.
02502 
02503    fWeight = w;
02504    TString opt = option;
02505    opt.ToLower();
02506    ResetBit(kGlobalWeight);
02507    if (opt.Contains("global")) {
02508       SetBit(kGlobalWeight);
02509    }
02510 }
02511 
02512 //______________________________________________________________________________
02513 void TChain::Streamer(TBuffer& b)
02514 {
02515    // -- Stream a class object.
02516 
02517    if (b.IsReading()) {
02518       UInt_t R__s, R__c;
02519       Version_t R__v = b.ReadVersion(&R__s, &R__c);
02520       if (R__v > 2) {
02521          b.ReadClassBuffer(TChain::Class(), this, R__v, R__s, R__c);
02522          return;
02523       }
02524       //====process old versions before automatic schema evolution
02525       TTree::Streamer(b);
02526       b >> fTreeOffsetLen;
02527       b >> fNtrees;
02528       fFiles->Streamer(b);
02529       if (R__v > 1) {
02530          fStatus->Streamer(b);
02531          fTreeOffset = new Long64_t[fTreeOffsetLen];
02532          b.ReadFastArray(fTreeOffset,fTreeOffsetLen);
02533       }
02534       b.CheckByteCount(R__s, R__c, TChain::IsA());
02535       //====end of old versions
02536 
02537    } else {
02538       b.WriteClassBuffer(TChain::Class(),this);
02539    }
02540 }
02541 
02542 //______________________________________________________________________________
02543 void TChain::UseCache(Int_t /* maxCacheSize */, Int_t /* pageSize */)
02544 {
02545    // -- Dummy function kept for back compatibility.
02546    // The cache is now activated automatically when processing TTrees/TChain.
02547 }

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