TRefTable.cxx

Go to the documentation of this file.
00001 // @(#)root/cont:$Id: TRefTable.cxx 37274 2010-12-04 21:31:29Z pcanal $
00002 // Author: Rene Brun   28/09/2001
00003 
00004 /*************************************************************************
00005  * Copyright (C) 1995-2004, 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 // A TRefTable maintains the association between a referenced object    //
00015 // and the parent object supporting this referenced object.             //
00016 //                                                                      //
00017 // The parent object is typically a branch of a TTree. For each object  //
00018 // referenced in a TTree entry, the corresponding entry in the TTree's  //
00019 // TBranchRef::fRefTable contains the index of the branch that          //
00020 // needs to be loaded to bring the object into memory.                  //
00021 //                                                                      //
00022 // Persistency of a TRefTable is split into two parts:                  //
00023 // * entry specific information is stored (read) by FillBuffer          //
00024 //   (ReadBuffer). For each referenced object the object's fUniqueID    //
00025 //   and the referencing TRef::fPID is stored (to allow the TRefTable   //
00026 //   to autoload references created by different processes).            //
00027 // * non-entry specific, i.e. global information is stored (read) by    //
00028 //   the Streamer function. This comprises all members marked as        //
00029 //   persistent.                                                        //
00030 //                                                                      //
00031 // As TObject::fUniqueID is only unique for a given TProcessID, a table //
00032 // of unique IDs is kept for each used TProcessID. There is no natural  //
00033 // order of TProcessIDs, so TRefTable stores a vector of the TGUID of   //
00034 // all known TProcessIDs in fProcessGUIDs; the index of a TProcessID in //
00035 // this vector defines the index of the auto-loading info in fParentIDs //
00036 // for that TProcessID. The mapping of TProcessID* to index is cached   //
00037 // for quick non-persistent lookup.                                     //
00038 //                                                                      //
00039 //////////////////////////////////////////////////////////////////////////
00040 
00041 #include "TRefTable.h"
00042 #include "TObjArray.h"
00043 #include "TProcessID.h"
00044 #include <algorithm>
00045 
00046 TRefTable *TRefTable::fgRefTable = 0;
00047 
00048 ClassImp(TRefTable)
00049 //______________________________________________________________________________
00050 TRefTable::TRefTable() : fNumPIDs(0), fAllocSize(0), fN(0), fParentIDs(0), fParentID(-1),
00051                          fDefaultSize(10), fUID(0), fUIDContext(0), fSize(0), fParents(0), fOwner(0)
00052 {
00053    // Default constructor for I/O.
00054 
00055    fgRefTable   = this;
00056 }
00057 
00058 //______________________________________________________________________________
00059 TRefTable::TRefTable(TObject *owner, Int_t size) :
00060      fNumPIDs(0), fAllocSize(0), fN(0), fParentIDs(0), fParentID(-1),
00061      fDefaultSize(size<10 ? 10 : size), fUID(0), fUIDContext(0), fSize(0), fParents(new TObjArray(1)), fOwner(owner)
00062 {
00063    // Create a TRefTable with initial size.
00064 
00065    fgRefTable   = this;
00066 }
00067 
00068 //______________________________________________________________________________
00069 TRefTable::~TRefTable()
00070 {
00071    // Destructor.
00072 
00073    delete [] fAllocSize;
00074    delete [] fN;
00075    for (Int_t pid = 0; pid < fNumPIDs; ++pid) {
00076       delete [] fParentIDs[pid];
00077    }
00078    delete [] fParentIDs;
00079    delete fParents;
00080    if (fgRefTable == this) fgRefTable = 0;
00081 }
00082 
00083 //______________________________________________________________________________
00084 Int_t TRefTable::Add(Int_t uid, TProcessID *context)
00085 {
00086    // Add a new uid to the table.
00087    // we add a new pair (uid,fparent) to the map
00088    // This function is called by TObject::Streamer or TStreamerInfo::WriteBuffer
00089 
00090    if (!context)
00091       context = TProcessID::GetSessionProcessID();
00092    Int_t iid = GetInternalIdxForPID(context);
00093 
00094    Int_t newsize = 0;
00095    uid = uid & 0xffffff;
00096    if (uid >= fAllocSize[iid]) {
00097       newsize = uid + uid / 2;
00098       if (newsize < fDefaultSize)
00099          newsize = fDefaultSize;
00100       newsize = ExpandForIID(iid, newsize);
00101    }
00102    if (newsize < 0) {
00103       Error("Add", "Cannot allocate space to store uid=%d", uid);
00104       return -1;
00105    }
00106    if (fParentID < 0) {
00107       Error("Add", "SetParent must be called before adding uid=%d", uid);
00108       return -1;
00109    }
00110    fParentIDs[iid][uid] = fParentID + 1;
00111    if (uid >= fN[iid]) fN[iid] = uid + 1;
00112    return uid;
00113 }
00114 
00115 
00116 //______________________________________________________________________________
00117 Int_t TRefTable::AddInternalIdxForPID(TProcessID *procid)
00118 {
00119    // Add the internal index for fProcessIDs, fAllocSize, etc given a PID.
00120 
00121    if (!procid)
00122       procid = TProcessID::GetSessionProcessID();
00123    Int_t pid = procid->GetUniqueID();
00124    if (fMapPIDtoInternal.size() <= (size_t) pid)
00125       fMapPIDtoInternal.resize(TProcessID::GetNProcessIDs(), -1);
00126 
00127    Int_t iid = fMapPIDtoInternal[pid];
00128    if (iid == -1) {
00129       // need to update
00130       iid = FindPIDGUID(procid->GetTitle());
00131       if (iid == -1) {
00132          fProcessGUIDs.push_back(procid->GetTitle());
00133          iid = fProcessGUIDs.size() - 1;
00134       }
00135       fMapPIDtoInternal[pid] = iid;
00136    }
00137 
00138    ExpandPIDs(iid + 1);
00139    return iid;
00140 }
00141 
00142 //______________________________________________________________________________
00143 void TRefTable::Clear(Option_t * /*option*/ )
00144 {
00145    // Clear all entries in the table.
00146 
00147    for (Int_t iid = 0; iid < fNumPIDs; ++iid) {
00148       memset(fParentIDs[iid], 0, sizeof(Int_t) * fN[iid]);
00149    }
00150    memset(fN, 0, sizeof(Int_t) * fNumPIDs);
00151    fParentID = -1;
00152 }
00153 
00154 //______________________________________________________________________________
00155 Int_t TRefTable::Expand(Int_t pid, Int_t newsize)
00156 {
00157    // Expand fParentIDs to newsize for ProcessID pid.
00158 
00159    Int_t iid = GetInternalIdxForPID(pid);
00160    if (iid < 0) return -1;
00161    return ExpandForIID(iid, newsize);
00162 }
00163 
00164 //______________________________________________________________________________
00165 Int_t TRefTable::ExpandForIID(Int_t iid, Int_t newsize)
00166 {
00167    // Expand fParentIDs to newsize for internel ProcessID index iid.
00168 
00169    if (newsize < 0)  return newsize;
00170    if (newsize != fAllocSize[iid]) {
00171       Int_t *temp = fParentIDs[iid];
00172       if (newsize != 0) {
00173          fParentIDs[iid] = new Int_t[newsize];
00174          if (newsize < fAllocSize[iid])
00175             memcpy(fParentIDs[iid], temp, newsize * sizeof(Int_t));
00176          else {
00177             memcpy(fParentIDs[iid], temp, fAllocSize[iid] * sizeof(Int_t));
00178             memset(&fParentIDs[iid][fAllocSize[iid]], 0,
00179                    (newsize - fAllocSize[iid]) * sizeof(Int_t));
00180          }
00181       } else {
00182          fParentIDs[iid] = 0;
00183       }
00184       if (fAllocSize[iid]) delete [] temp;
00185       fAllocSize[iid] = newsize;
00186    }
00187    return newsize;
00188 }
00189 
00190 //______________________________________________________________________________
00191 void TRefTable::ExpandPIDs(Int_t numpids)
00192 {
00193    // Expand the arrays of managed PIDs
00194 
00195    if (numpids <= fNumPIDs) return;
00196 
00197    // else add to internal tables
00198    Int_t oldNumPIDs = fNumPIDs;
00199    fNumPIDs  = numpids;
00200 
00201    Int_t *temp = fAllocSize;
00202    fAllocSize = new Int_t[fNumPIDs];
00203    if (temp) memcpy(fAllocSize, temp, oldNumPIDs * sizeof(Int_t));
00204    memset(&fAllocSize[oldNumPIDs], 0,
00205           (fNumPIDs - oldNumPIDs) * sizeof(Int_t));
00206    delete [] temp;
00207 
00208    temp = fN;
00209    fN = new Int_t[fNumPIDs];
00210    if (temp) memcpy(fN, temp, oldNumPIDs * sizeof(Int_t));
00211    memset(&fN[oldNumPIDs], 0, (fNumPIDs - oldNumPIDs) * sizeof(Int_t));
00212    delete [] temp;
00213 
00214    Int_t **temp2 = fParentIDs;
00215    fParentIDs = new Int_t *[fNumPIDs];
00216    if (temp2) memcpy(fParentIDs, temp2, oldNumPIDs * sizeof(Int_t *));
00217    memset(&fParentIDs[oldNumPIDs], 0,
00218           (fNumPIDs - oldNumPIDs) * sizeof(Int_t*));
00219 }
00220 
00221 //______________________________________________________________________________
00222 void TRefTable::FillBuffer(TBuffer & b)
00223 {
00224    // Fill buffer b with the fN elements in fParentdIDs.
00225    // This function is called by TBranchRef::FillLeaves.
00226 
00227    b << -fNumPIDs; // write out "-" to signal new TRefTable buffer format using PID table
00228    for (Int_t iid = 0; iid < fNumPIDs; ++iid) {
00229       b << fN[iid];
00230       b.WriteFastArray(fParentIDs[iid], fN[iid]);
00231    }
00232 }
00233 
00234 
00235 //______________________________________________________________________________
00236 Int_t TRefTable::FindPIDGUID(const char *guid) const
00237 {
00238    // Get fProcessGUIDs' index of the TProcessID with GUID guid
00239    std::vector<std::string>::const_iterator posPID
00240       = std::find(fProcessGUIDs.begin(), fProcessGUIDs.end(), guid);
00241    if (posPID == fProcessGUIDs.end()) return -1;
00242    return posPID - fProcessGUIDs.begin();
00243 }
00244 
00245 //______________________________________________________________________________
00246 TObject *TRefTable::GetParent(Int_t uid, TProcessID *context /* =0 */ ) const 
00247 {
00248    // Return object corresponding to uid.
00249    if (!fParents) return 0;
00250 
00251    Int_t iid = -1;
00252    if (!context) context = TProcessID::GetSessionProcessID();
00253    iid = GetInternalIdxForPID(context);
00254 
00255    uid = uid & 0xFFFFFF;
00256    if (uid < 0 || uid >= fN[iid]) return 0;
00257    Int_t pnumber = fParentIDs[iid][uid] - 1;
00258    Int_t nparents = fParents->GetEntriesFast();
00259    if (pnumber < 0 || pnumber >= nparents) return 0;
00260    return fParents->UncheckedAt(pnumber);
00261 }
00262 
00263 //______________________________________________________________________________
00264 Int_t TRefTable::GetInternalIdxForPID(TProcessID *procid) const 
00265 {
00266    // Get the index for fProcessIDs, fAllocSize, etc given a PID.
00267    // Uses fMapPIDtoInternal and the pid's GUID / fProcessGUID 
00268 
00269    return const_cast <TRefTable*>(this)->AddInternalIdxForPID(procid);
00270 }
00271 
00272 //______________________________________________________________________________
00273 Int_t TRefTable::GetInternalIdxForPID(Int_t pid) const 
00274 {
00275    // Get the index for fProcessIDs, fAllocSize, etc given a PID.
00276    // Uses fMapPIDtoInternal and the pid's GUID / fProcessGUID
00277 
00278    return GetInternalIdxForPID(TProcessID::GetProcessID(pid));
00279 }
00280 
00281 
00282 //______________________________________________________________________________
00283 TRefTable *TRefTable::GetRefTable()
00284 {
00285    // Static function returning the current TRefTable.
00286 
00287    return fgRefTable;
00288 }
00289 
00290 //______________________________________________________________________________
00291 Bool_t TRefTable::Notify()
00292 {
00293    // This function is called by TRef::Streamer or TStreamerInfo::ReadBuffer
00294    // when reading a reference.
00295    // This function, in turns, notifies the TRefTable owner for action.
00296    // eg, when the owner is a TBranchRef, TBranchRef::Notify is called
00297    // to read the branch containing the referenced object.
00298 
00299    return fOwner->Notify();
00300 }
00301 
00302 //______________________________________________________________________________
00303 void TRefTable::ReadBuffer(TBuffer &b)
00304 {
00305    // Fill buffer b with the fN elements in fParentdIDs.
00306    // This function is called by TBranchRef::ReadLeaves
00307 
00308    Int_t firstInt = 0;          // we don't know yet what it means
00309    b >> firstInt;
00310 
00311    Int_t numIids = -1;
00312    Int_t startIid = 0;
00313    if (firstInt < 0) numIids = -firstInt; // new format
00314    else {
00315       // old format, only one PID
00316       numIids = 1;
00317 
00318       TProcessID *fileProcessID = b.GetLastProcessID(this);
00319 
00320       startIid = GetInternalIdxForPID(fileProcessID);
00321       if (startIid == -1) {
00322          fProcessGUIDs.push_back(fileProcessID->GetTitle());
00323          startIid = fProcessGUIDs.size() - 1;
00324       }
00325       numIids += startIid;
00326    }
00327 
00328    ExpandPIDs(numIids);
00329    for (Int_t iid = startIid; iid < numIids; ++iid) {
00330       Int_t newN = 0;
00331       if (firstInt < 0) b >> newN;
00332       else newN = firstInt;
00333       if (newN > fAllocSize[iid])
00334          ExpandForIID(iid, newN + newN / 2);
00335       fN[iid] = newN;
00336       b.ReadFastArray(fParentIDs[iid], fN[iid]);
00337    }
00338 }
00339 
00340 //______________________________________________________________________________
00341 void TRefTable::Reset(Option_t * /*option*/ )
00342 {
00343    // Clear all entries in the table.
00344    Clear();
00345    if (fParents) fParents->Clear();
00346 }
00347 
00348 //______________________________________________________________________________
00349 Int_t TRefTable::SetParent(const TObject* parent, Int_t branchID)
00350 {
00351    // -- Set current parent object, typically a branch of a tree.
00352    //
00353    // This function is called by TBranchElement::Fill() and by
00354    // TBranchElement::GetEntry().
00355    //
00356    if (!fParents) {
00357       return -1;
00358    }
00359    Int_t nparents = fParents->GetEntriesFast();
00360    if (branchID != -1) {
00361       // -- The branch already has an id cached, just use it.
00362       fParentID = branchID;
00363    }
00364    else {
00365       // -- The branch does *not* have an id cached, find it or generate one.
00366       // Lookup the branch.
00367       fParentID = fParents->IndexOf(parent);
00368       if (fParentID < 0) {
00369          // -- The branch is not known, generate an id number.
00370          fParents->AddAtAndExpand(const_cast<TObject*>(parent), nparents);
00371          fParentID = nparents;
00372       }
00373    }
00374    return fParentID;
00375 }
00376 
00377 //______________________________________________________________________________
00378 void TRefTable::SetRefTable(TRefTable *table)
00379 {
00380    // Static function setting the current TRefTable.
00381 
00382    fgRefTable = table;
00383 }
00384 
00385 //______________________________________________________________________________
00386 void TRefTable::Streamer(TBuffer &R__b)
00387 {
00388    // Stream an object of class TRefTable.
00389 
00390    if (R__b.IsReading()) {
00391       R__b.ReadClassBuffer(TRefTable::Class(),this);
00392    } else {
00393       R__b.WriteClassBuffer(TRefTable::Class(),this);
00394       //make sure that all TProcessIDs referenced in the Tree are put to the buffer
00395       //this is important in case the buffer is a TMessage to be sent through a TSocket
00396 #if 0
00397       TObjArray *pids = TProcessID::GetPIDs();
00398       Int_t npids = pids->GetEntries();
00399       Int_t npid2 = fProcessGUIDs.size();
00400       for (Int_t i = 0; i < npid2; i++) {
00401          TProcessID *pid;
00402          for (Int_t ipid = 0;ipid<npids;ipid++) {
00403             pid = (TProcessID*)pids->At(ipid);
00404             if (!pid) continue;
00405             if (!strcmp(pid->GetTitle(),fProcessGUIDs[i].c_str()))
00406                R__b.WriteProcessID(pid);
00407          }
00408       }
00409 #endif
00410    }
00411 }

Generated on Tue Jul 5 14:11:37 2011 for ROOT_528-00b_version by  doxygen 1.5.1