library: libCore
#include "TClonesArray.h"

TClonesArray


class description - source file - inheritance tree (.pdf)

class TClonesArray : public TObjArray

Inheritance Chart:
TObject
<-
TCollection
<-
TSeqCollection
<-
TObjArray
<-
TClonesArray

    public:
TClonesArray() TClonesArray(const char* classname, Int_t size = 1000, Bool_t call_dtor = kFALSE) TClonesArray(const TClass* cl, Int_t size = 1000, Bool_t call_dtor = kFALSE) virtual ~TClonesArray() virtual void AddAfter(const TObject*, TObject*) virtual void AddAt(TObject*, Int_t) virtual void AddAtAndExpand(TObject*, Int_t) virtual Int_t AddAtFree(TObject*) virtual void AddBefore(const TObject*, TObject*) virtual void AddFirst(TObject*) virtual void AddLast(TObject*) TObject* AddrAt(Int_t idx) void BypassStreamer(Bool_t bypass = kTRUE) Bool_t CanBypassStreamer() const static TClass* Class() virtual void Clear(Option_t* option = "") virtual void Compress() virtual void Delete(Option_t* option = "") virtual void Expand(Int_t newSize) virtual void ExpandCreate(Int_t n) virtual void ExpandCreateFast(Int_t n) TClass* GetClass() const virtual TClass* IsA() const TObject* New(Int_t idx) virtual TObject*& operator[](Int_t idx) virtual TObject* operator[](Int_t idx) const virtual TObject* Remove(TObject* obj) virtual TObject* RemoveAt(Int_t idx) virtual void ShowMembers(TMemberInspector& insp, char* parent) virtual void Sort(Int_t upto = kMaxInt) virtual void Streamer(TBuffer& b) void StreamerNVirtual(TBuffer& b)

Data Members


    protected:
TClass* fClass !Pointer to the class TObjArray* fKeep !Saved copies of pointers to objects public:
static const enum TClonesArray:: kForgetBits static const enum TClonesArray:: kNoSplit static const enum TClonesArray:: kBypassStreamer

Class Description

                                                                      
 An array of clone (identical) objects. Memory for the objects        
 stored in the array is allocated only once in the lifetime of the    
 clones array. All objects must be of the same class. For the rest    
 this class has the same properties as TObjArray.                     
                                                                      
 To reduce the very large number of new and delete calls in large     
 loops like this (O(100000) x O(10000) times new/delete):             
                                                                      
   TObjArray a(10000);                                                
   while (TEvent *ev = (TEvent *)next()) {      // O(100000) events   
      for (int i = 0; i < ev->Ntracks; i++) {   // O(10000) tracks    
         a[i] = new TTrack(x,y,z,...);                                
         ...                                                          
         ...                                                          
      }                                                               
      ...                                                             
      a.Delete();                                                     
   }                                                                  
                                                                      
 One better uses a TClonesArray which reduces the number of           
 new/delete calls to only O(10000):                                   
                                                                      
   TCloneArray a("TTrack", 10000);                                    
   while (TEvent *ev = (TEvent *)next()) {      // O(100000) events   
      for (int i = 0; i < ev->Ntracks; i++) {   // O(10000) tracks    
         new(a[i]) TTrack(x,y,z,...);                                 
         ...                                                          
         ...                                                          
      }                                                               
      ...                                                             
      a.Delete();                                                     
   }                                                                  
                                                                      
 Note: the only supported way to add objects to a TClonesArray is     
 via the new with placement method. The diffrent Add() methods of     
 TObjArray and its base classes are not allowed.                      
                                                                      
 Considering that a new/delete costs about 70 mus on a 300 MHz HP,    
 O(10^9) new/deletes will save about 19 hours.                        
                                                                      


TClonesArray() : TObjArray()

TClonesArray(const char *classname, Int_t s, Bool_t) : TObjArray(s)
 Create an array of clone objects of classname. The class must inherit from
 TObject. If the class defines an own operator delete(), make sure that
 it looks like this:

    void MyClass::operator delete(void *vp)
    {
       if ((Long_t) vp != TObject::GetDtorOnly())
          ::operator delete(vp);       // delete space
       else
          TObject::SetDtorOnly(0);
    }

 The second argument s indicates an approximate number of objects
 that will be entered in the array. If more than s objects are entered,
 the array will be automatically expanded.

 The third argument is not used anymore and only there for backward
 compatibility reasons.

 In case you want to send a TClonesArray (or object containing a
 TClonesArray) via a TMessage over a TSocket don't forget to call
 BypassStreamer(kFALSE). See TClonesArray::BypassStreamer().

TClonesArray(const TClass *cl, Int_t s, Bool_t) : TObjArray(s)
 Create an array of clone objects of class cl. The class must inherit from
 TObject. If the class defines an own operator delete(), make sure that
 it looks like this:

    void MyClass::operator delete(void *vp)
    {
       if ((Long_t) vp != TObject::GetDtorOnly())
          ::operator delete(vp);       // delete space
       else
          TObject::SetDtorOnly(0);
    }

 The second argument s indicates an approximate number of objects
 that will be entered in the array. If more than s objects are entered,
 the array will be automatically expanded.

 The third argument is not used anymore and only there for backward
 compatibility reasons.

 In case you want to send a TClonesArray (or object containing a
 TClonesArray) via a TMessage over a TSocket don't forget to call
 BypassStreamer(kFALSE). See TClonesArray::BypassStreamer().

~TClonesArray()
 Delete a clones array.

void BypassStreamer(Bool_t bypass)
 When the kBypassStreamer bit is set, the automatically
 generated Streamer can call directly TClass::WriteBuffer.
 Bypassing the Streamer improves the performance when writing/reading
 the objects in the TClonesArray. However there is a drawback:
 When a TClonesArray is written with split=0 bypassing the Streamer,
 the StreamerInfo of the class in the array being optimized,
 one cannot use later the TClonesArray with split>0. For example,
 there is a problem with the following scenario:
  1- A class Foo has a TClonesArray of Bar objects
  2- The Foo object is written with split=0 to Tree T1.
     In this case the StreamerInfo for the class Bar is created
     in optimized mode in such a way that data members of the same type
     are written as an array improving the I/O performance.
  3- In a new program, T1 is read and a new Tree T2 is created
     with the object Foo in split>1
  4- When the T2 branch is created, the StreamerInfo for the class Bar
     is created with no optimization (mandatory for the split mode).
     The optimized Bar StreamerInfo is going to be used to read
     the TClonesArray in T1. The result will be Bar objects with
     data member values not in the right sequence.
 The solution to this problem is to call BypassStreamer(kFALSE)
 for the TClonesArray. In this case, the normal Bar::Streamer function
 will be called. The Bar::Streamer function works OK independently
 if the Bar StreamerInfo had been generated in optimized mode or not.
 In case you want to send Foo via a TMessage over a TSocket you also
 need to disable the streamer bypass.

void Compress()
 Remove empty slots from array.

void Clear(Option_t *option)
 Clear the clones array. Only use this routine when your objects don't
 allocate memory since it will not call the object dtors.
 However, if the class in the TClonesArray implements the function
 Clear(Option_t *option) and if option = "C" the function Clear()
 is called for all objects in the array. In the function Clear(), one
 can delete objects or dynamic arrays allocated in the class.
 This procedure is much faster than calling TClonesArray::Delete().

void Delete(Option_t *)
 Clear the clones array. Use this routine when your objects allocate
 memory (e.g. objects inheriting from TNamed or containing TStrings
 allocate memory). If not you better use Clear() since if is faster.

void Expand(Int_t newSize)
 Expand or shrink the array to newSize elements.

void ExpandCreate(Int_t n)
 Expand or shrink the array to n elements and create the clone
 objects by calling their default ctor. If n is less than the current size
 the array is shrinked and the allocated space is freed.
 This routine is typically used to create a clonesarray into which
 one can directly copy object data without going via the
 "new (arr[i]) MyObj()" (i.e. the vtbl is already set correctly).

void ExpandCreateFast(Int_t n)
 Expand or shrink the array to n elements and create the clone
 objects by calling their default ctor. If n is less than the current size
 the array is shrinked and the allocated space is freed.
 This routine is typically used to create a clonesarray into which
 one can directly copy object data without going via the
 "new (arr[i]) MyObj()" (i.e. the vtbl is already set correctly).
 This is a simplified version of ExpandCreate used in the TTree mechanism.

TObject* RemoveAt(Int_t idx)
 Remove object at index idx.

TObject* Remove(TObject *obj)
 Remove object from array.

void Sort(Int_t upto)
 If objects in array are sortable (i.e. IsSortable() returns true
 for all objects) then sort array.

void Streamer(TBuffer &b)
 Write all objects in array to the I/O buffer. ATTENTION: empty slots
 are also stored (using one byte per slot). If you don't want this
 use a TOrdCollection or TList.

TObject* New(Int_t idx)
 Create an object of type fClass with the default ctor at the specified
 index. Returns 0 in case of error.



Inline Functions


            TClass* GetClass() const
               void AddFirst(TObject*)
               void AddLast(TObject*)
               void AddAt(TObject*, Int_t)
               void AddAtAndExpand(TObject*, Int_t)
              Int_t AddAtFree(TObject*)
               void AddAfter(const TObject*, TObject*)
               void AddBefore(const TObject*, TObject*)
             Bool_t CanBypassStreamer() const
           TObject* AddrAt(Int_t idx)
          TObject*& operator[](Int_t idx)
           TObject* operator[](Int_t idx) const
            TClass* Class()
            TClass* IsA() const
               void ShowMembers(TMemberInspector& insp, char* parent)
               void StreamerNVirtual(TBuffer& b)


Author: Rene Brun 11/02/96
Last update: root/cont:$Name: $:$Id: TClonesArray.cxx,v 1.46 2005/06/09 18:20:02 pcanal Exp $
Copyright (C) 1995-2000, Rene Brun and Fons Rademakers. *


ROOT page - Class index - Class Hierarchy - Top of the page

This page has been automatically generated. If you have any comments or suggestions about the page layout send a mail to ROOT support, or contact the developers with any questions or problems regarding ROOT.