Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TRefTable.cxx
Go to the documentation of this file.
1// @(#)root/cont:$Id$
2// Author: Rene Brun 28/09/2001
3
4/*************************************************************************
5 * Copyright (C) 1995-2004, Rene Brun and Fons Rademakers. *
6 * All rights reserved. *
7 * *
8 * For the licensing terms see $ROOTSYS/LICENSE. *
9 * For the list of contributors see $ROOTSYS/README/CREDITS. *
10 *************************************************************************/
11
12/** \class TRefTable
13\ingroup Containers
14A TRefTable maintains the association between a referenced object
15and the parent object supporting this referenced object.
16
17The parent object is typically a branch of a TTree. For each object
18referenced in a TTree entry, the corresponding entry in the TTree's
19TBranchRef::fRefTable contains the index of the branch that
20needs to be loaded to bring the object into memory.
21
22Persistency of a TRefTable is split into two parts:
23 - entry specific information is stored (read) by FillBuffer
24 (ReadBuffer). For each referenced object the object's fUniqueID
25 and the referencing TRef::fPID is stored (to allow the TRefTable
26 to autoload references created by different processes).
27 - non-entry specific, i.e. global information is stored (read) by
28 the Streamer function. This comprises all members marked as
29 persistent.
30
31As TObject::fUniqueID is only unique for a given TProcessID, a table
32of unique IDs is kept for each used TProcessID. There is no natural
33order of TProcessIDs, so TRefTable stores a vector of the TGUID of
34all known TProcessIDs in fProcessGUIDs; the index of a TProcessID in
35this vector defines the index of the auto-loading info in fParentIDs
36for that TProcessID. The mapping of TProcessID* to index is cached
37for quick non-persistent lookup.
38*/
39
40#include "TRefTable.h"
41#include "TBuffer.h"
42#include "TObjArray.h"
43#include "TProcessID.h"
44#include <algorithm>
45
47
49////////////////////////////////////////////////////////////////////////////////
50/// Default constructor for I/O.
51
52TRefTable::TRefTable() : fNumPIDs(0), fAllocSize(0), fN(0), fParentIDs(0), fParentID(-1),
53 fDefaultSize(10), fUID(0), fUIDContext(0), fSize(0), fParents(0), fOwner(0)
54{
55 fgRefTable = this;
56}
57
58////////////////////////////////////////////////////////////////////////////////
59/// Create a TRefTable with initial size.
60
62 fNumPIDs(0), fAllocSize(0), fN(0), fParentIDs(0), fParentID(-1),
63 fDefaultSize(size<10 ? 10 : size), fUID(0), fUIDContext(0), fSize(0), fParents(new TObjArray(1)), fOwner(owner)
64{
65 fgRefTable = this;
66}
67
68////////////////////////////////////////////////////////////////////////////////
69/// Destructor.
70
72{
73 delete [] fAllocSize;
74 delete [] fN;
75 for (Int_t pid = 0; pid < fNumPIDs; ++pid) {
76 delete [] fParentIDs[pid];
77 }
78 delete [] fParentIDs;
79 delete fParents;
80 if (fgRefTable == this) fgRefTable = 0;
81}
82
83////////////////////////////////////////////////////////////////////////////////
84/// Add a new uid to the table.
85/// we add a new pair (uid,fparent) to the map
86/// This function is called by TObject::Streamer or TStreamerInfo::WriteBuffer
87
89{
90 if (!context)
92 Int_t iid = GetInternalIdxForPID(context);
93
94 Int_t newsize = 0;
95 uid = uid & 0xffffff;
96 if (uid >= fAllocSize[iid]) {
97 newsize = uid + uid / 2;
98 if (newsize < fDefaultSize)
99 newsize = fDefaultSize;
100 newsize = ExpandForIID(iid, newsize);
101 }
102 if (newsize < 0) {
103 Error("Add", "Cannot allocate space to store uid=%d", uid);
104 return -1;
105 }
106 if (fParentID < 0) {
107 Error("Add", "SetParent must be called before adding uid=%d", uid);
108 return -1;
109 }
110 fParentIDs[iid][uid] = fParentID + 1;
111 if (uid >= fN[iid]) fN[iid] = uid + 1;
112 return uid;
113}
114
115
116////////////////////////////////////////////////////////////////////////////////
117/// Add the internal index for fProcessIDs, fAllocSize, etc given a PID.
118
120{
121 if (!procid)
123 Int_t pid = procid->GetUniqueID();
124 if (fMapPIDtoInternal.size() <= (size_t) pid)
126
127 Int_t iid = fMapPIDtoInternal[pid];
128 if (iid == -1) {
129 // need to update
130 iid = FindPIDGUID(procid->GetTitle());
131 if (iid == -1) {
132 fProcessGUIDs.push_back(procid->GetTitle());
133 iid = fProcessGUIDs.size() - 1;
134 }
135 fMapPIDtoInternal[pid] = iid;
136 }
137
138 ExpandPIDs(iid + 1);
139 return iid;
140}
141
142////////////////////////////////////////////////////////////////////////////////
143/// Clear all entries in the table.
144
145void TRefTable::Clear(Option_t * /*option*/ )
146{
147 for (Int_t iid = 0; iid < fNumPIDs; ++iid) {
148 memset(fParentIDs[iid], 0, sizeof(Int_t) * fN[iid]);
149 }
150 memset(fN, 0, sizeof(Int_t) * fNumPIDs);
151 fParentID = -1;
152}
153
154////////////////////////////////////////////////////////////////////////////////
155/// Expand fParentIDs to newsize for ProcessID pid.
156
158{
159 Int_t iid = GetInternalIdxForPID(pid);
160 if (iid < 0) return -1;
161 return ExpandForIID(iid, newsize);
162}
163
164////////////////////////////////////////////////////////////////////////////////
165/// Expand fParentIDs to newsize for internel ProcessID index iid.
166
168{
169 if (newsize < 0) return newsize;
170 if (newsize != fAllocSize[iid]) {
171 Int_t *temp = fParentIDs[iid];
172 if (newsize != 0) {
173 fParentIDs[iid] = new Int_t[newsize];
174 if (newsize < fAllocSize[iid])
175 memcpy(fParentIDs[iid], temp, newsize * sizeof(Int_t));
176 else {
177 memcpy(fParentIDs[iid], temp, fAllocSize[iid] * sizeof(Int_t));
178 memset(&fParentIDs[iid][fAllocSize[iid]], 0,
179 (newsize - fAllocSize[iid]) * sizeof(Int_t));
180 }
181 } else {
182 fParentIDs[iid] = 0;
183 }
184 if (fAllocSize[iid]) delete [] temp;
185 fAllocSize[iid] = newsize;
186 }
187 return newsize;
188}
189
190////////////////////////////////////////////////////////////////////////////////
191/// Expand the arrays of managed PIDs
192
194{
195 if (numpids <= fNumPIDs) return;
196
197 // else add to internal tables
198 Int_t oldNumPIDs = fNumPIDs;
199 fNumPIDs = numpids;
200
201 Int_t *temp = fAllocSize;
203 if (temp) memcpy(fAllocSize, temp, oldNumPIDs * sizeof(Int_t));
204 memset(&fAllocSize[oldNumPIDs], 0,
205 (fNumPIDs - oldNumPIDs) * sizeof(Int_t));
206 delete [] temp;
207
208 temp = fN;
209 fN = new Int_t[fNumPIDs];
210 if (temp) memcpy(fN, temp, oldNumPIDs * sizeof(Int_t));
211 memset(&fN[oldNumPIDs], 0, (fNumPIDs - oldNumPIDs) * sizeof(Int_t));
212 delete [] temp;
213
214 Int_t **temp2 = fParentIDs;
215 fParentIDs = new Int_t *[fNumPIDs];
216 if (temp2) memcpy(fParentIDs, temp2, oldNumPIDs * sizeof(Int_t *));
217 memset(&fParentIDs[oldNumPIDs], 0,
218 (fNumPIDs - oldNumPIDs) * sizeof(Int_t*));
219}
220
221////////////////////////////////////////////////////////////////////////////////
222/// Fill buffer b with the fN elements in fParentdIDs.
223/// This function is called by TBranchRef::FillLeaves.
224
226{
227 b << -fNumPIDs; // write out "-" to signal new TRefTable buffer format using PID table
228 for (Int_t iid = 0; iid < fNumPIDs; ++iid) {
229 b << fN[iid];
230 b.WriteFastArray(fParentIDs[iid], fN[iid]);
231 }
232}
233
234
235////////////////////////////////////////////////////////////////////////////////
236/// Get fProcessGUIDs' index of the TProcessID with GUID guid
237
238Int_t TRefTable::FindPIDGUID(const char *guid) const
239{
240 std::vector<std::string>::const_iterator posPID
241 = std::find(fProcessGUIDs.begin(), fProcessGUIDs.end(), guid);
242 if (posPID == fProcessGUIDs.end()) return -1;
243 return posPID - fProcessGUIDs.begin();
244}
245
246////////////////////////////////////////////////////////////////////////////////
247/// Return object corresponding to uid.
248
249TObject *TRefTable::GetParent(Int_t uid, TProcessID *context /* =0 */ ) const
250{
251 if (!fParents) return 0;
252
253 Int_t iid = -1;
254 if (!context) context = TProcessID::GetSessionProcessID();
255 iid = GetInternalIdxForPID(context);
256
257 uid = uid & 0xFFFFFF;
258 if (uid < 0 || uid >= fN[iid]) return 0;
259 Int_t pnumber = fParentIDs[iid][uid] - 1;
260 Int_t nparents = fParents->GetEntriesFast();
261 if (pnumber < 0 || pnumber >= nparents) return 0;
262 return fParents->UncheckedAt(pnumber);
263}
264
265////////////////////////////////////////////////////////////////////////////////
266/// Get the index for fProcessIDs, fAllocSize, etc given a PID.
267/// Uses fMapPIDtoInternal and the pid's GUID / fProcessGUID
268
270{
271 return const_cast <TRefTable*>(this)->AddInternalIdxForPID(procid);
272}
273
274////////////////////////////////////////////////////////////////////////////////
275/// Get the index for fProcessIDs, fAllocSize, etc given a PID.
276/// Uses fMapPIDtoInternal and the pid's GUID / fProcessGUID
277
279{
281}
282
283
284////////////////////////////////////////////////////////////////////////////////
285/// Static function returning the current TRefTable.
286
288{
289 return fgRefTable;
290}
291
292////////////////////////////////////////////////////////////////////////////////
293/// This function is called by TRef::Streamer or TStreamerInfo::ReadBuffer
294/// when reading a reference.
295/// This function, in turns, notifies the TRefTable owner for action.
296/// eg, when the owner is a TBranchRef, TBranchRef::Notify is called
297/// to read the branch containing the referenced object.
298
300{
301 return fOwner->Notify();
302}
303
304////////////////////////////////////////////////////////////////////////////////
305/// Fill buffer b with the fN elements in fParentdIDs.
306/// This function is called by TBranchRef::ReadLeaves
307
309{
310 Int_t firstInt = 0; // we don't know yet what it means
311 b >> firstInt;
312
313 Int_t numIids = -1;
314 Int_t startIid = 0;
315 if (firstInt < 0) numIids = -firstInt; // new format
316 else {
317 // old format, only one PID
318 numIids = 1;
319
320 TProcessID *fileProcessID = b.GetLastProcessID(this);
321
322 startIid = GetInternalIdxForPID(fileProcessID);
323 if (startIid == -1) {
324 fProcessGUIDs.push_back(fileProcessID->GetTitle());
325 startIid = fProcessGUIDs.size() - 1;
326 }
327 numIids += startIid;
328 }
329
330 ExpandPIDs(numIids);
331 for (Int_t iid = startIid; iid < numIids; ++iid) {
332 Int_t newN = 0;
333 if (firstInt < 0) b >> newN;
334 else newN = firstInt;
335 if (newN > fAllocSize[iid])
336 ExpandForIID(iid, newN + newN / 2);
337 fN[iid] = newN;
338 b.ReadFastArray(fParentIDs[iid], fN[iid]);
339 }
340}
341
342////////////////////////////////////////////////////////////////////////////////
343/// Clear all entries in the table.
344
345void TRefTable::Reset(Option_t * /*option*/ )
346{
347 Clear();
348 if (fParents) fParents->Clear();
349}
350
351////////////////////////////////////////////////////////////////////////////////
352/// -- Set current parent object, typically a branch of a tree.
353///
354/// This function is called by TBranchElement::Fill() and by
355/// TBranchElement::GetEntry().
356
357Int_t TRefTable::SetParent(const TObject* parent, Int_t branchID)
358{
359 if (!fParents) {
360 return -1;
361 }
362 Int_t nparents = fParents->GetEntriesFast();
363 if (branchID != -1) {
364 // -- The branch already has an id cached, just use it.
365 fParentID = branchID;
366 }
367 else {
368 // -- The branch does *not* have an id cached, find it or generate one.
369 // Lookup the branch.
370 fParentID = fParents->IndexOf(parent);
371 if (fParentID < 0) {
372 // -- The branch is not known, generate an id number.
373 fParents->AddAtAndExpand(const_cast<TObject*>(parent), nparents);
374 fParentID = nparents;
375 }
376 }
377 return fParentID;
378}
379
380////////////////////////////////////////////////////////////////////////////////
381/// Static function setting the current TRefTable.
382
384{
385 fgRefTable = table;
386}
387
388////////////////////////////////////////////////////////////////////////////////
389/// Stream an object of class TRefTable.
390
391void TRefTable::Streamer(TBuffer &R__b)
392{
393 if (R__b.IsReading()) {
394 R__b.ReadClassBuffer(TRefTable::Class(),this);
395 } else {
396 R__b.WriteClassBuffer(TRefTable::Class(),this);
397 //make sure that all TProcessIDs referenced in the Tree are put to the buffer
398 //this is important in case the buffer is a TMessage to be sent through a TSocket
399#if 0
401 Int_t npids = pids->GetEntries();
402 Int_t npid2 = fProcessGUIDs.size();
403 for (Int_t i = 0; i < npid2; i++) {
404 TProcessID *pid;
405 for (Int_t ipid = 0;ipid<npids;ipid++) {
406 pid = (TProcessID*)pids->At(ipid);
407 if (!pid) continue;
408 if (!strcmp(pid->GetTitle(),fProcessGUIDs[i].c_str()))
409 R__b.WriteProcessID(pid);
410 }
411 }
412#endif
413 }
414}
size_t fSize
#define b(i)
Definition RSha256.hxx:100
size_t size(const MatrixT &matrix)
retrieve the size of a square matrix
int Int_t
Definition RtypesCore.h:45
const char Option_t
Definition RtypesCore.h:66
#define ClassImp(name)
Definition Rtypes.h:364
Buffer base class used for serializing objects.
Definition TBuffer.h:43
virtual UShort_t WriteProcessID(TProcessID *pid)=0
Always return 0 (current processID).
Definition TBuffer.cxx:353
virtual Int_t ReadClassBuffer(const TClass *cl, void *pointer, const TClass *onfile_class=0)=0
Bool_t IsReading() const
Definition TBuffer.h:86
virtual Int_t WriteClassBuffer(const TClass *cl, void *pointer)=0
virtual const char * GetTitle() const
Returns title of object.
Definition TNamed.h:48
An array of TObjects.
Definition TObjArray.h:31
Int_t IndexOf(const TObject *obj) const
Int_t GetEntriesFast() const
Definition TObjArray.h:58
virtual void AddAtAndExpand(TObject *obj, Int_t idx)
Add object at position idx.
Int_t GetEntries() const
Return the number of objects in array (i.e.
virtual void Clear(Option_t *option="")
Remove all objects from the array.
TObject * UncheckedAt(Int_t i) const
Definition TObjArray.h:84
TObject * At(Int_t idx) const
Definition TObjArray.h:164
Mother of all ROOT objects.
Definition TObject.h:41
virtual Bool_t Notify()
This method must be overridden to handle object notification.
Definition TObject.cxx:578
virtual UInt_t GetUniqueID() const
Return the unique object id.
Definition TObject.cxx:447
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:963
A TProcessID identifies a ROOT job in a unique way in time and space.
Definition TProcessID.h:74
static UInt_t GetNProcessIDs()
Return the (static) number of process IDs.
static TObjArray * GetPIDs()
static: returns array of TProcessIDs
static TProcessID * GetSessionProcessID()
static function returning the pointer to the session TProcessID
static TProcessID * GetProcessID(UShort_t pid)
static function returning a pointer to TProcessID number pid in fgPIDs
A TRefTable maintains the association between a referenced object and the parent object supporting th...
Definition TRefTable.h:35
Int_t * fN
[fNumPIDs] allocated size of array fParentIDs for each ProcessID
Definition TRefTable.h:40
virtual Int_t SetParent(const TObject *parent, Int_t branchID)
– Set current parent object, typically a branch of a tree.
std::vector< Int_t > fMapPIDtoInternal
Definition TRefTable.h:50
TRefTable()
Default constructor for I/O.
Definition TRefTable.cxx:52
virtual Int_t ExpandForIID(Int_t iid, Int_t newsize)
Expand fParentIDs to newsize for internel ProcessID index iid.
Int_t GetInternalIdxForPID(TProcessID *procid) const
Get the index for fProcessIDs, fAllocSize, etc given a PID.
Int_t fNumPIDs
Definition TRefTable.h:38
virtual Bool_t Notify()
This function is called by TRef::Streamer or TStreamerInfo::ReadBuffer when reading a reference.
TObjArray * fParents
Definition TRefTable.h:47
virtual void FillBuffer(TBuffer &b)
Fill buffer b with the fN elements in fParentdIDs.
Int_t FindPIDGUID(const char *guid) const
Get fProcessGUIDs' index of the TProcessID with GUID guid.
Int_t * fAllocSize
number of known ProcessIDs
Definition TRefTable.h:39
Int_t AddInternalIdxForPID(TProcessID *procid)
Add the internal index for fProcessIDs, fAllocSize, etc given a PID.
static TRefTable * fgRefTable
cache of pid to index in fProcessGUIDs
Definition TRefTable.h:51
static void SetRefTable(TRefTable *table)
Static function setting the current TRefTable.
Int_t fDefaultSize
current parent ID in fParents (latest call to SetParent)
Definition TRefTable.h:43
static TRefTable * GetRefTable()
Static function returning the current TRefTable.
virtual void Reset(Option_t *="")
Clear all entries in the table.
void ExpandPIDs(Int_t numpids)
Expand the arrays of managed PIDs.
TObject * GetParent(Int_t uid, TProcessID *context=0) const
Return object corresponding to uid.
std::vector< std::string > fProcessGUIDs
Definition TRefTable.h:49
virtual Int_t Expand(Int_t pid, Int_t newsize)
Expand fParentIDs to newsize for ProcessID pid.
virtual ~TRefTable()
Destructor.
Definition TRefTable.cxx:71
Int_t ** fParentIDs
[fNumPIDs] current maximum number of IDs in array fParentIDs for each ProcessID
Definition TRefTable.h:41
virtual void Clear(Option_t *="")
Clear all entries in the table.
virtual Int_t Add(Int_t uid, TProcessID *context=0)
Add a new uid to the table.
Definition TRefTable.cxx:88
Int_t fParentID
[fNumPIDs][fAllocSize] array of Parent IDs
Definition TRefTable.h:42
TObject * fOwner
Definition TRefTable.h:48
virtual void ReadBuffer(TBuffer &b)
Fill buffer b with the fN elements in fParentdIDs.