//*CMZ :  2.23/08 01/11/99  21.38.07  by  Fons Rademakers
//*CMZ :  2.23/07 27/10/99  10.47.27  by  Fons Rademakers
//*CMZ :  2.22/07 05/07/99  18.54.53  by  Fons Rademakers
//*CMZ :  2.22/01 18/05/99  12.44.57  by  Fons Rademakers
//*CMZ :  2.21/01 12/01/99  11.47.03  by  Fons Rademakers
//*CMZ :  2.20/00 12/11/98  13.57.41  by  Fons Rademakers
//*-- Author :    Fons Rademakers   10/08/95

//*KEEP,CopyRight,T=C.
/*************************************************************************
 * Copyright(c) 1995-1999, The ROOT System, All rights reserved.         *
 * Authors: Rene Brun and Fons Rademakers.                               *
 *                                                                       *
 * For the licensing terms see $ROOTSYS/AA_LICENSE.                      *
 * For the list of contributors see $ROOTSYS/AA_CREDITS.                 *
 *************************************************************************/
//*KEND.

//////////////////////////////////////////////////////////////////////////
//                                                                      //
// TList                                                                //
//                                                                      //
// A doubly linked list. All classes inheriting from TObject can be     //
// inserted in a TList. Before being inserted into the list the object  //
// pointer is wrapped in a TObjLink object which contains, besides      //
// the object pointer also a previous and next pointer.                 //
//                                                                      //
// There are basically four ways to iterate over a TList (in order      //
// of preference, if not forced by other constraints):                  //
//    1) Using the ForEach macro:                                       //
//         GetListOfPrimitives()->ForEach(TObject,Paint)(option);       //
//                                                                      //
//    2) Using the TList iterator TListIter (via the wrapper class      //
//       TIter):                                                        //
//         TIter next(GetListOfPrimitives());                           //
//         while (TObject *obj = next())                                //
//            obj->Draw(next.GetOption());                              //
//                                                                      //
//    3) Using the TObjLink list entries (that wrap the TObject*):      //
//         TObjLink *lnk = GetListOfPrimitives()->FirstLink();          //
//         while (lnk) {                                                //
//            lnk->GetObject()->Draw(lnk->GetOption());                 //
//            lnk = lnk->Next();                                        //
//         }                                                            //
//                                                                      //
//    4) Using the TList's After() and Before() member functions:       //
//         TFree *idcur = this;                                         //
//         while (idcur) {                                              //
//            ...                                                       //
//            ...                                                       //
//            idcur = (TFree*)GetListOfFree()->After(idcur);            //
//         }                                                            //
//                                                                      //
//   Methods 2, 3 and 4 can also easily iterate backwards using either  //
//   a backward TIter (using argument kIterBackward) or by using        //
//   LastLink() and lnk->Prev() or by using the Before() member.        //
//                             //
//                                                                      //
//////////////////////////////////////////////////////////////////////////

//*KEEP,TList.
#include "TList.h"
//*KEND.


ClassImp(TList)

//______________________________________________________________________________
 TList::~TList()
{
   // Delete the list. Objects are not deleted.

   Clear();
}

//______________________________________________________________________________
 void TList::AddFirst(TObject *obj)
{
   // Add object at the beginning of the list.

   if (IsArgNull("AddFirst", obj)) return;

   if (!fFirst) {
      fFirst = NewLink(obj);
      fLast = fFirst;
   } else {
      TObjLink *t = NewLink(obj);
      t->fNext = fFirst;
      fFirst->fPrev = t;
      fFirst = t;
   }
   fSize++;
}

//______________________________________________________________________________
 void TList::AddFirst(TObject *obj, Option_t *opt)
{
   // Add object at the beginning of the list and also store option.
   // Storing an option is useful when one wants to change the behaviour
   // of an object a little without having to create a complete new
   // copy of the object. This feature is used, for example, by the Draw()
   // method. It allows the same object to be drawn in different ways.

   if (IsArgNull("AddFirst", obj)) return;

   if (!fFirst) {
      fFirst = NewOptLink(obj, opt);
      fLast = fFirst;
   } else {
      TObjLink *t = NewOptLink(obj, opt);
      t->fNext = fFirst;
      fFirst->fPrev = t;
      fFirst = t;
   }
   fSize++;
}

//______________________________________________________________________________
 void TList::AddLast(TObject *obj)
{
   // Add object at the end of the list.

   if (IsArgNull("AddLast", obj)) return;

   if (!fFirst) {
      fFirst = NewLink(obj);
      fLast  = fFirst;
   } else
      fLast = NewLink(obj, fLast);
   fSize++;
}

//______________________________________________________________________________
 void TList::AddLast(TObject *obj, Option_t *opt)
{
   // Add object at the end of the list and also store option.
   // Storing an option is useful when one wants to change the behaviour
   // of an object a little without having to create a complete new
   // copy of the object. This feature is used, for example, by the Draw()
   // method. It allows the same object to be drawn in different ways.

   if (IsArgNull("AddLast", obj)) return;

   if (!fFirst) {
      fFirst = NewOptLink(obj, opt);
      fLast  = fFirst;
   } else
      fLast = NewOptLink(obj, opt, fLast);
   fSize++;
}

//______________________________________________________________________________
 void TList::AddBefore(TObject *before, TObject *obj)
{
   // Insert object before object before in the list.

   if (IsArgNull("AddBefore", obj)) return;

   if (!before)
      TList::AddFirst(obj);
   else {
      Int_t    idx;
      TObjLink *t = FindLink(before, idx);
      if (!t) {
         Error("AddBefore", "before not found, object not added");
         return;
      }
      if (t == fFirst)
         TList::AddFirst(obj);
      else {
         NewLink(obj, t->Prev());
         fSize++;
      }
   }
}

//______________________________________________________________________________
 void TList::AddBefore(TObjLink *before, TObject *obj)
{
   // Insert object before the specified ObjLink object. If before = 0 then add
   // to the head of the list. An ObjLink can be obtained by looping over a list
   // using the above describe iterator method 3.

   if (IsArgNull("AddBefore", obj)) return;

   if (!before)
      TList::AddFirst(obj);
   else {
      if (before == fFirst)
         TList::AddFirst(obj);
      else {
         NewLink(obj, before->Prev());
         fSize++;
      }
   }
}

//______________________________________________________________________________
 void TList::AddAfter(TObject *after, TObject *obj)
{
   // Insert object after object after in the list.

   if (IsArgNull("AddAfter", obj)) return;

   if (!after)
      TList::AddLast(obj);
   else {
      Int_t    idx;
      TObjLink *t = FindLink(after, idx);
      if (!t) {
         Error("AddAfter", "after not found, object not added");
         return;
      }
      if (t == fLast)
         TList::AddLast(obj);
      else {
         NewLink(obj, t);
         fSize++;
      }
   }
}

//______________________________________________________________________________
 void TList::AddAfter(TObjLink *after, TObject *obj)
{
   // Insert object after the specified ObjLink object. If after = 0 then add
   // to the tail of the list. An ObjLink can be obtained by looping over a list
   // using the above describe iterator method 3.

   if (IsArgNull("AddAfter", obj)) return;

   if (!after)
      TList::AddLast(obj);
   else {
      if (after == fLast)
         TList::AddLast(obj);
      else {
         NewLink(obj, after);
         fSize++;
      }
   }
}

//______________________________________________________________________________
 void TList::AddAt(TObject *obj, Int_t idx)
{
   // Insert object at position idx in the list.

   if (IsArgNull("AddAt", obj)) return;

   TObjLink *lnk = LinkAt(idx);
   if (!lnk)
      TList::AddLast(obj);
   else if (lnk == fFirst)
      TList::AddFirst(obj);
   else {
      NewLink(obj, lnk->Prev());
      fSize++;
   }
}

//______________________________________________________________________________
 TObject *TList::After(TObject *obj) const
{
   // Returns the object after object obj. Returns 0 if obj is last in list.

   TObjLink *t;

   if (fCache && fCache->GetObject() == obj) {
      t = fCache;
      ((TList*)this)->fCache = fCache->Next();  //cast const away, fCache should be mutable
   } else {
      Int_t idx;
      t = FindLink(obj, idx);
      if (t) ((TList*)this)->fCache = t->Next();
   }

   if (t && t->Next())
      return t->Next()->GetObject();
   else
      return 0;
}

//______________________________________________________________________________
 TObject *TList::At(Int_t idx) const
{
   // Returns the object at position idx. Returns 0 if idx is out of range.

   TObjLink *lnk = LinkAt(idx);
   if (lnk) return lnk->GetObject();
   return 0;
}

//______________________________________________________________________________
 TObject *TList::Before(TObject *obj) const
{
   // Returns the object before object obj. Returns 0 if obj is first in list.

   TObjLink *t;

   if (fCache && fCache->GetObject() == obj) {
      t = fCache;
      ((TList*)this)->fCache = fCache->Prev();  //cast const away, fCache should be mutable
   } else {
      Int_t idx;
      t = FindLink(obj, idx);
      if (t) ((TList*)this)->fCache = t->Prev();
   }

   if (t && t->Prev())
      return t->Prev()->GetObject();
   else
      return 0;
}

//______________________________________________________________________________
 void TList::Clear(Option_t *option)
{
   // Remove all objects from the list. Does not delete the objects.

   Bool_t nodel = !strcmp(option, "nodelete") ? kTRUE : kFALSE;

   while (fFirst) {
      TObjLink *tlk = fFirst;
      fFirst = fFirst->Next();
      fSize--;
      // delete only heap objects marked OK to clear
      if (!nodel && tlk->GetObject() && tlk->GetObject()->IsOnHeap()) {
         if (tlk->GetObject()->TestBit(kCanDelete))
            TCollection::GarbageCollect(tlk->GetObject());
      }
      delete tlk;
   }
   fFirst = fLast = fCache = 0;
   fSize  = 0;
}

//______________________________________________________________________________
 void TList::Delete(Option_t *)
{
   // Remove all objects from the list AND delete all heap based objects.

   while (fFirst) {
      TObjLink *tlk = fFirst;
      fFirst = fFirst->Next();
      fSize--;
      // delete only heap objects
      if (tlk->GetObject() && tlk->GetObject()->IsOnHeap()) {
         TCollection::GarbageCollect(tlk->GetObject());
      }
      delete tlk;
   }
   fFirst = fLast = fCache = 0;
   fSize  = 0;
}

//______________________________________________________________________________
 void TList::DeleteLink(TObjLink *lnk)
{
   // Delete a TObjLink object.

   lnk->fNext = lnk->fPrev = 0;
   lnk->fObject = 0;
   delete lnk;
}

//______________________________________________________________________________
 TObject *TList::FindObject(const Text_t *name) const
{
   // Find an object in this list using its name. Requires a sequential
   // scan till the object has been found. Returns 0 if object with specified
   // name is not found. This method overrides the generic FindObject()
   // of TCollection for efficiency reasons.

   if (!name) return 0;
   TObjLink *lnk = FirstLink();
   while (lnk) {
      TObject *obj = lnk->GetObject();
      if (obj->GetName() && !strcmp(name, obj->GetName())) return obj;
      lnk = lnk->Next();
   }
   return 0;
}

//______________________________________________________________________________
 TObject *TList::FindObject(TObject *obj) const
{
   // Find an object in this list using the object's IsEqual()
   // member function. Requires a sequential scan till the object has
   // been found. Returns 0 if object is not found.
   // This method overrides the generic FindObject() of TCollection for
   // efficiency reasons.

   TObjLink *lnk = FirstLink();

   while (lnk) {
      TObject *ob = lnk->GetObject();
      if (ob->IsEqual(obj)) return ob;
      lnk = lnk->Next();
   }
   return 0;
}

//______________________________________________________________________________
 TObjLink *TList::FindLink(TObject *obj, Int_t &idx) const
{
   // Returns the TObjLink object that contains object obj. In idx it returns
   // the position of the object in the list.

   if (!fFirst) return 0;

   TObject *object;
   TObjLink *lnk = fFirst;
   idx = 0;

   while (lnk) {
      object = lnk->GetObject();
      if (object) {
         if (object->TestBit(kNotDeleted)) {
            if (object->IsEqual(obj)) return lnk;
         }
      }
      lnk = lnk->Next();
      idx++;
   }
   return 0;
}

//______________________________________________________________________________
 TObject *TList::First() const
{
   // Return the first object in the list. Returns 0 when list is empty.

   if (fFirst) return fFirst->GetObject();
   return 0;
}

//______________________________________________________________________________
 TObject *TList::Last() const
{
   // Return the last object in the list. Returns 0 when list is empty.

   if (fLast) return fLast->GetObject();
   return 0;
}

//______________________________________________________________________________
 TObjLink *TList::LinkAt(Int_t idx) const
{
   // Return the TObjLink object at index idx.

   Int_t    i = 0;
   TObjLink *lnk = fFirst;
   while (i < idx && lnk) {
      i++;
      lnk = lnk->Next();
   }
   return lnk;
}

//______________________________________________________________________________
 TIterator *TList::MakeIterator(Bool_t dir) const
{
   // Return a list iterator.

   return new TListIter(this, dir);
}

//______________________________________________________________________________
 TObjLink *TList::NewLink(TObject *obj, TObjLink *prev)
{
   // Return a new TObjLink.

   if (prev)
      return new TObjLink(obj, prev);
   else
      return new TObjLink(obj);
}

//______________________________________________________________________________
 TObjLink *TList::NewOptLink(TObject *obj, Option_t *opt, TObjLink *prev)
{
   // Return a new TObjOptLink (a TObjLink that also stores the option).

   if (prev)
      return new TObjOptLink(obj, prev, opt);
   else
      return new TObjOptLink(obj, opt);
}

//______________________________________________________________________________
 TObject *TList::Remove(TObject *obj)
{
   // Remove object from the list.

   if (!obj) return 0;

   Int_t    idx;
   TObjLink *lnk = FindLink(obj, idx);

   if (!lnk) return 0;

   // return object found, which may be (pointer wise) different than the
   // input object (depending on what IsEqual() is doing)
   TObject *ob = lnk->GetObject();

   if (lnk == fFirst) {
      fFirst = lnk->Next();
      if (lnk == fLast) fLast = fFirst;
      else              fFirst->fPrev = 0;
      DeleteLink(lnk);
   } else if (lnk == fLast) {
      fLast = lnk->Prev();
      fLast->fNext = 0;
      DeleteLink(lnk);
   } else {
      lnk->Prev()->fNext = lnk->Next();
      lnk->Next()->fPrev = lnk->Prev();
      DeleteLink(lnk);
   }
   fSize--;
   fCache = 0;

   return ob;
}

//______________________________________________________________________________
 TObject *TList::Remove(TObjLink *lnk)
{
   // Remove object link (and therefore the object it contains)
   // from the list.

   if (!lnk) return 0;

   TObject *obj = lnk->GetObject();

   if (lnk == fFirst) {
      fFirst = lnk->Next();
      if (lnk == fLast) fLast = fFirst;
      else              fFirst->fPrev = 0;
      DeleteLink(lnk);
   } else if (lnk == fLast) {
      fLast = lnk->Prev();
      fLast->fNext = 0;
      DeleteLink(lnk);
   } else {
      lnk->Prev()->fNext = lnk->Next();
      lnk->Next()->fPrev = lnk->Prev();
      DeleteLink(lnk);
   }
   fSize--;
   fCache = 0;

   return obj;
}


//______________________________________________________________________________
 void TList::Sort(Bool_t order)
{
   // Sort linked list. Real sorting is done in private function DoSort().
   // The list can only be sorted when is contains objects of a sortable
   // class.

   if (!fFirst) return;

   fAscending = order;

   if (!fFirst->GetObject()->IsSortable()) {
      Error("Sort", "objects in list are not sortable");
      return;
   }

   DoSort(&fFirst, fSize);

   // correct back links
   TObjLink *ol, *lnk = fFirst;

   if (lnk) lnk->fPrev = 0;
   while ((ol = lnk)) {
      lnk = lnk->fNext;
      if (lnk)
         lnk->fPrev = ol;
      else
         fLast = ol;
   }
}

//______________________________________________________________________________
 Bool_t TList::LnkCompare(TObjLink *l1, TObjLink *l2)
{
   // Compares the objects stored in the TObjLink objects.
   // Depending on the flag IsAscending() the function returns
   // true if the object in l1 <= l2 (ascending) or l2 <= l1 (descending).

   Int_t cmp = l1->GetObject()->Compare(l2->GetObject());

   if ((IsAscending() && cmp <=0) || (!IsAscending() && cmp > 0))
      return kTRUE;
   return kFALSE;
}

//______________________________________________________________________________
 TObjLink **TList::DoSort(TObjLink **head, Int_t n)
{
   // Sort linked list.

   TObjLink *p1, *p2, **h2, **t2;

   switch (n) {
      case 0:
         return head;

      case 1:
         return &((*head)->fNext);

      case 2:
         p2 = (p1 = *head)->fNext;
         if (LnkCompare(p1, p2)) return &(p2->fNext);
         p1->fNext = (*head = p2)->fNext;
         return &((p2->fNext = p1)->fNext);
   }

   int m;
   n -= m = n / 2;

   t2 = DoSort(h2 = DoSort(head, n), m);

   if (LnkCompare((p1 = *head), (p2 = *h2))) {
      do {
         if (!--n) return *h2 = p2, t2;
      } while (LnkCompare((p1 = *(head = &(p1->fNext))), p2));
   }

   while (1) {
      *head = p2;
      do {
         if (!--m) return *h2 = *t2, *t2 = p1, h2;
      } while (!LnkCompare(p1, (p2 = *(head = &(p2->fNext)))));
      *head = p1;
      do {
         if (!--n) return *h2 = p2, t2;
      } while (LnkCompare((p1 = *(head = &(p1->fNext))), p2));
   }
}

//////////////////////////////////////////////////////////////////////////
//                                                                      //
// TObjLink                                                             //
//                                                                      //
// Wrapper around a TObject so it can be stored in a TList.             //
//                                                                      //
//////////////////////////////////////////////////////////////////////////

//______________________________________________________________________________
TObjLink::TObjLink(TObject *obj, TObjLink *prev)
          : fNext(prev->fNext), fPrev(prev), fObject(obj)
{
   // Create a new TObjLink.

   fPrev->fNext = this;
   if (fNext) fNext->fPrev = this;
}

//////////////////////////////////////////////////////////////////////////
//                                                                      //
// TListIter                                                            //
//                                                                      //
// Iterator of linked list.                                             //
//                                                                      //
//////////////////////////////////////////////////////////////////////////

ClassImp(TListIter)

//______________________________________________________________________________
TListIter::TListIter(const TList *l, Bool_t dir)
        : fList(l), fCurCursor(0), fCursor(0), fDirection(dir), fStarted(kFALSE)
{
   // Create a new list iterator. By default the iteration direction
   // is kIterForward. To go backward use kIterBackward.
}

//______________________________________________________________________________
TListIter::TListIter(const TListIter &iter)
{
   // Copy ctor.

   fList      = iter.fList;
   fCurCursor = iter.fCurCursor;
   fCursor    = iter.fCursor;
   fDirection = iter.fDirection;
   fStarted   = iter.fStarted;
}

//______________________________________________________________________________
TIterator &TListIter::operator=(const TIterator &rhs)
{
   // Overridden assignment operator.

   if (this != &rhs && rhs.IsA() == TListIter::Class()) {
      const TListIter &rhs1 = (const TListIter &)rhs;
      fList      = rhs1.fList;
      fCurCursor = rhs1.fCurCursor;
      fCursor    = rhs1.fCursor;
      fDirection = rhs1.fDirection;
      fStarted   = rhs1.fStarted;
   }
   return *this;
}

//______________________________________________________________________________
TListIter &TListIter::operator=(const TListIter &rhs)
{
   // Overloaded assignment operator.

   if (this != &rhs) {
      fList      = rhs.fList;
      fCurCursor = rhs.fCurCursor;
      fCursor    = rhs.fCursor;
      fDirection = rhs.fDirection;
      fStarted   = rhs.fStarted;
   }
   return *this;
}

//______________________________________________________________________________
TObject *TListIter::Next()
{
   // Return next object in the list. Returns 0 when no more objects in list.

   if (!fList) return 0;

   if (fDirection == kIterForward) {
      if (!fStarted) {
         fCursor = fList->fFirst;
         fStarted = kTRUE;
      }
      fCurCursor = fCursor;
      if (fCursor) fCursor = fCursor->Next();
   } else {
      if (!fStarted) {
         fCursor = fList->fLast;
         fStarted = kTRUE;
      }
      fCurCursor = fCursor;
      if (fCursor) fCursor = fCursor->Prev();
   }

   if (fCurCursor) return fCurCursor->GetObject();
   return 0;
}

//______________________________________________________________________________
Option_t *TListIter::GetOption() const
{
   // Returns the object option stored in the list.

   if (fCurCursor) return fCurCursor->GetOption();
   return "";
}

//______________________________________________________________________________
void TListIter::SetOption(Option_t *option)
{
   // Sets the object option stored in the list.

   if (fCurCursor) fCurCursor->SetOption(option);
}


ROOT page - Class index - 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.