TList.cxx

Go to the documentation of this file.
00001 // @(#)root/cont:$Id: TList.cxx 38184 2011-02-21 14:53:55Z rdm $
00002 // Author: Fons Rademakers   10/08/95
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 // TList                                                                //
00015 //                                                                      //
00016 // A doubly linked list. All classes inheriting from TObject can be     //
00017 // inserted in a TList. Before being inserted into the list the object  //
00018 // pointer is wrapped in a TObjLink object which contains, besides      //
00019 // the object pointer also a previous and next pointer.                 //
00020 //                                                                      //
00021 // There are basically four ways to iterate over a TList (in order      //
00022 // of preference, if not forced by other constraints):                  //
00023 //    1) Using the R__FOR_EACH macro:                                   //
00024 //         GetListOfPrimitives()->R__FOR_EACH(TObject,Paint)(option);   //
00025 //                                                                      //
00026 //    2) Using the TList iterator TListIter (via the wrapper class      //
00027 //       TIter):                                                        //
00028 //         TIter next(GetListOfPrimitives());                           //
00029 //         while ((TObject *obj = next()))                              //
00030 //            obj->Draw(next.GetOption());                              //
00031 //                                                                      //
00032 //    3) Using the TList iterator TListIter and std::for_each           //
00033 //       algorithm:                                                     //
00034 //         // A function object, which will be applied to each element  //
00035 //         // of the given range.                                       //
00036 //         struct STestFunctor {                                        //
00037 //            bool operator()(TObject *aObj) {                          //
00038 //               ...                                                    //
00039 //               return true;                                           //
00040 //            }                                                         //
00041 //        }                                                             //
00042 //        ...                                                           //
00043 //        ...                                                           //
00044 //        TIter iter(mylist);                                           //
00045 //        for_each( iter.Begin(), TIter::End(), STestFunctor() );       //
00046 //                                                                      //
00047 //    4) Using the TObjLink list entries (that wrap the TObject*):      //
00048 //         TObjLink *lnk = GetListOfPrimitives()->FirstLink();          //
00049 //         while (lnk) {                                                //
00050 //            lnk->GetObject()->Draw(lnk->GetOption());                 //
00051 //            lnk = lnk->Next();                                        //
00052 //         }                                                            //
00053 //                                                                      //
00054 //    5) Using the TList's After() and Before() member functions:       //
00055 //         TFree *idcur = this;                                         //
00056 //         while (idcur) {                                              //
00057 //            ...                                                       //
00058 //            ...                                                       //
00059 //            idcur = (TFree*)GetListOfFree()->After(idcur);            //
00060 //         }                                                            //
00061 //                                                                      //
00062 //   Methods 2, 3 and 4 can also easily iterate backwards using either  //
00063 //   a backward TIter (using argument kIterBackward) or by using        //
00064 //   LastLink() and lnk->Prev() or by using the Before() member.        //
00065 //Begin_Html <img src=gif/tlist.gif> End_Html                           //
00066 //                                                                      //
00067 //////////////////////////////////////////////////////////////////////////
00068 
00069 #include "TList.h"
00070 #include "TClass.h"
00071 
00072 #include <string>
00073 namespace std {} using namespace std;
00074 
00075 ClassImp(TList)
00076 
00077 //______________________________________________________________________________
00078 TList::~TList()
00079 {
00080    // Delete the list. Objects are not deleted unless the TList is the
00081    // owner (set via SetOwner()).
00082 
00083    Clear();
00084 }
00085 
00086 //______________________________________________________________________________
00087 void TList::AddFirst(TObject *obj)
00088 {
00089    // Add object at the beginning of the list.
00090 
00091    if (IsArgNull("AddFirst", obj)) return;
00092 
00093    if (!fFirst) {
00094       fFirst = NewLink(obj);
00095       fLast = fFirst;
00096    } else {
00097       TObjLink *t = NewLink(obj);
00098       t->fNext = fFirst;
00099       fFirst->fPrev = t;
00100       fFirst = t;
00101    }
00102    fSize++;
00103    Changed();
00104 }
00105 
00106 //______________________________________________________________________________
00107 void TList::AddFirst(TObject *obj, Option_t *opt)
00108 {
00109    // Add object at the beginning of the list and also store option.
00110    // Storing an option is useful when one wants to change the behaviour
00111    // of an object a little without having to create a complete new
00112    // copy of the object. This feature is used, for example, by the Draw()
00113    // method. It allows the same object to be drawn in different ways.
00114 
00115    if (IsArgNull("AddFirst", obj)) return;
00116 
00117    if (!fFirst) {
00118       fFirst = NewOptLink(obj, opt);
00119       fLast = fFirst;
00120    } else {
00121       TObjLink *t = NewOptLink(obj, opt);
00122       t->fNext = fFirst;
00123       fFirst->fPrev = t;
00124       fFirst = t;
00125    }
00126    fSize++;
00127    Changed();
00128 }
00129 
00130 //______________________________________________________________________________
00131 void TList::AddLast(TObject *obj)
00132 {
00133    // Add object at the end of the list.
00134 
00135    if (IsArgNull("AddLast", obj)) return;
00136 
00137    if (!fFirst) {
00138       fFirst = NewLink(obj);
00139       fLast  = fFirst;
00140    } else
00141       fLast = NewLink(obj, fLast);
00142    fSize++;
00143    Changed();
00144 }
00145 
00146 //______________________________________________________________________________
00147 void TList::AddLast(TObject *obj, Option_t *opt)
00148 {
00149    // Add object at the end of the list and also store option.
00150    // Storing an option is useful when one wants to change the behaviour
00151    // of an object a little without having to create a complete new
00152    // copy of the object. This feature is used, for example, by the Draw()
00153    // method. It allows the same object to be drawn in different ways.
00154 
00155    if (IsArgNull("AddLast", obj)) return;
00156 
00157    if (!fFirst) {
00158       fFirst = NewOptLink(obj, opt);
00159       fLast  = fFirst;
00160    } else
00161       fLast = NewOptLink(obj, opt, fLast);
00162    fSize++;
00163    Changed();
00164 }
00165 
00166 //______________________________________________________________________________
00167 void TList::AddBefore(const TObject *before, TObject *obj)
00168 {
00169    // Insert object before object before in the list.
00170 
00171    if (IsArgNull("AddBefore", obj)) return;
00172 
00173    if (!before)
00174       TList::AddFirst(obj);
00175    else {
00176       Int_t    idx;
00177       TObjLink *t = FindLink(before, idx);
00178       if (!t) {
00179          Error("AddBefore", "before not found, object not added");
00180          return;
00181       }
00182       if (t == fFirst)
00183          TList::AddFirst(obj);
00184       else {
00185          NewLink(obj, t->Prev());
00186          fSize++;
00187          Changed();
00188       }
00189    }
00190 }
00191 
00192 //______________________________________________________________________________
00193 void TList::AddBefore(TObjLink *before, TObject *obj)
00194 {
00195    // Insert object before the specified ObjLink object. If before = 0 then add
00196    // to the head of the list. An ObjLink can be obtained by looping over a list
00197    // using the above describe iterator method 3.
00198 
00199    if (IsArgNull("AddBefore", obj)) return;
00200 
00201    if (!before)
00202       TList::AddFirst(obj);
00203    else {
00204       if (before == fFirst)
00205          TList::AddFirst(obj);
00206       else {
00207          NewLink(obj, before->Prev());
00208          fSize++;
00209          Changed();
00210       }
00211    }
00212 }
00213 
00214 //______________________________________________________________________________
00215 void TList::AddAfter(const TObject *after, TObject *obj)
00216 {
00217    // Insert object after object after in the list.
00218 
00219    if (IsArgNull("AddAfter", obj)) return;
00220 
00221    if (!after)
00222       TList::AddLast(obj);
00223    else {
00224       Int_t    idx;
00225       TObjLink *t = FindLink(after, idx);
00226       if (!t) {
00227          Error("AddAfter", "after not found, object not added");
00228          return;
00229       }
00230       if (t == fLast)
00231          TList::AddLast(obj);
00232       else {
00233          NewLink(obj, t);
00234          fSize++;
00235          Changed();
00236       }
00237    }
00238 }
00239 
00240 //______________________________________________________________________________
00241 void TList::AddAfter(TObjLink *after, TObject *obj)
00242 {
00243    // Insert object after the specified ObjLink object. If after = 0 then add
00244    // to the tail of the list. An ObjLink can be obtained by looping over a list
00245    // using the above describe iterator method 3.
00246 
00247    if (IsArgNull("AddAfter", obj)) return;
00248 
00249    if (!after)
00250       TList::AddLast(obj);
00251    else {
00252       if (after == fLast)
00253          TList::AddLast(obj);
00254       else {
00255          NewLink(obj, after);
00256          fSize++;
00257          Changed();
00258       }
00259    }
00260 }
00261 
00262 //______________________________________________________________________________
00263 void TList::AddAt(TObject *obj, Int_t idx)
00264 {
00265    // Insert object at position idx in the list.
00266 
00267    if (IsArgNull("AddAt", obj)) return;
00268 
00269    TObjLink *lnk = LinkAt(idx);
00270    if (!lnk)
00271       TList::AddLast(obj);
00272    else if (lnk == fFirst)
00273       TList::AddFirst(obj);
00274    else {
00275       NewLink(obj, lnk->Prev());
00276       fSize++;
00277       Changed();
00278    }
00279 }
00280 
00281 //______________________________________________________________________________
00282 TObject *TList::After(const TObject *obj) const
00283 {
00284    // Returns the object after object obj. Obj is found using the
00285    // object's IsEqual() method.  Returns 0 if obj is last in list.
00286 
00287    TObjLink *t;
00288 
00289    if (fCache && fCache->GetObject() && fCache->GetObject()->IsEqual(obj)) {
00290       t = fCache;
00291       ((TList*)this)->fCache = fCache->Next();  //cast const away, fCache should be mutable
00292    } else {
00293       Int_t idx;
00294       t = FindLink(obj, idx);
00295       if (t) ((TList*)this)->fCache = t->Next();
00296    }
00297 
00298    if (t && t->Next())
00299       return t->Next()->GetObject();
00300    else
00301       return 0;
00302 }
00303 
00304 //______________________________________________________________________________
00305 TObject *TList::At(Int_t idx) const
00306 {
00307    // Returns the object at position idx. Returns 0 if idx is out of range.
00308 
00309    TObjLink *lnk = LinkAt(idx);
00310    if (lnk) return lnk->GetObject();
00311    return 0;
00312 }
00313 
00314 //______________________________________________________________________________
00315 TObject *TList::Before(const TObject *obj) const
00316 {
00317    // Returns the object before object obj. Obj is found using the
00318    // object's IsEqual() method.  Returns 0 if obj is first in list.
00319 
00320    TObjLink *t;
00321 
00322    if (fCache && fCache->GetObject() && fCache->GetObject()->IsEqual(obj)) {
00323       t = fCache;
00324       ((TList*)this)->fCache = fCache->Prev();  //cast const away, fCache should be mutable
00325    } else {
00326       Int_t idx;
00327       t = FindLink(obj, idx);
00328       if (t) ((TList*)this)->fCache = t->Prev();
00329    }
00330 
00331    if (t && t->Prev())
00332       return t->Prev()->GetObject();
00333    else
00334       return 0;
00335 }
00336 
00337 //______________________________________________________________________________
00338 void TList::Clear(Option_t *option)
00339 {
00340    // Remove all objects from the list. Does not delete the objects
00341    // unless the TList is the owner (set via SetOwner()) and option
00342    // "nodelete" is not set.
00343    // If option="nodelete" then don't delete any heap objects that were
00344    // marked with the kCanDelete bit, otherwise these objects will be
00345    // deleted (this option is used by THashTable::Clear()).
00346 
00347    Bool_t nodel = option ? (!strcmp(option, "nodelete") ? kTRUE : kFALSE) : kFALSE;
00348 
00349    if (!nodel && IsOwner()) {
00350       Delete(option);
00351       return;
00352    }
00353 
00354    while (fFirst) {
00355       TObjLink *tlk = fFirst;
00356       fFirst = fFirst->Next();
00357       fSize--;
00358       // delete only heap objects marked OK to clear
00359       if (!nodel && tlk->GetObject() && tlk->GetObject()->IsOnHeap()) {
00360          if (tlk->GetObject()->TestBit(kCanDelete)) {
00361             if(tlk->GetObject()->TestBit(kNotDeleted)) {
00362                TCollection::GarbageCollect(tlk->GetObject());
00363             }
00364          }
00365       }
00366       delete tlk;
00367    }
00368    fFirst = fLast = fCache = 0;
00369    fSize  = 0;
00370    Changed();
00371 }
00372 
00373 //______________________________________________________________________________
00374 void TList::Delete(Option_t *option)
00375 {
00376    // Remove all objects from the list AND delete all heap based objects.
00377    // If option="slow" then keep list consistent during delete. This allows
00378    // recursive list operations during the delete (e.g. during the dtor
00379    // of an object in this list one can still access the list to search for
00380    // other not yet deleted objects).
00381 
00382    Bool_t slow = option ? (!strcmp(option, "slow") ? kTRUE : kFALSE) : kFALSE;
00383 
00384    TList removeDirectory; // need to deregistere these from their directory
00385 
00386    if (slow) {
00387 
00388       while (fFirst) {
00389          TObjLink *tlk = fFirst;
00390          fFirst = fFirst->Next();
00391          fSize--;
00392          // delete only heap objects
00393          if (tlk->GetObject() && tlk->GetObject()->IsOnHeap())
00394             TCollection::GarbageCollect(tlk->GetObject());
00395          else if (tlk->GetObject() && tlk->GetObject()->IsA()->GetDirectoryAutoAdd())
00396             removeDirectory.Add(tlk->GetObject());
00397 
00398          delete tlk;
00399       }
00400       fFirst = fLast = fCache = 0;
00401       fSize  = 0;
00402 
00403    } else {
00404 
00405       TObjLink *first = fFirst;    //pointer to first entry in linked list
00406       fFirst = fLast = fCache = 0;
00407       fSize  = 0;
00408       while (first) {
00409          TObjLink *tlk = first;
00410          first = first->Next();
00411          // delete only heap objects
00412          if (tlk->GetObject() && tlk->GetObject()->IsOnHeap())
00413             TCollection::GarbageCollect(tlk->GetObject());
00414          else if (tlk->GetObject() && tlk->GetObject()->IsA()->GetDirectoryAutoAdd())
00415             removeDirectory.Add(tlk->GetObject());
00416 
00417          delete tlk;
00418       }
00419    }
00420 
00421    // These objects cannot expect to have a valid TDirectory anymore;
00422    // e.g. because *this is the TDirectory's list of objects. Even if
00423    // not, they are supposed to be deleted, so we can as well unregister
00424    // them from their directory, even if they are stack-based:
00425    TIter iRemDir(&removeDirectory);
00426    TObject* dirRem = 0;
00427    while ((dirRem = iRemDir())) {
00428       (*dirRem->IsA()->GetDirectoryAutoAdd())(dirRem, 0);
00429    }
00430    Changed();
00431 }
00432 
00433 //______________________________________________________________________________
00434 void TList::DeleteLink(TObjLink *lnk)
00435 {
00436    // Delete a TObjLink object.
00437 
00438    lnk->fNext = lnk->fPrev = 0;
00439    lnk->fObject = 0;
00440    delete lnk;
00441 }
00442 
00443 //______________________________________________________________________________
00444 TObject *TList::FindObject(const char *name) const
00445 {
00446    // Find an object in this list using its name. Requires a sequential
00447    // scan till the object has been found. Returns 0 if object with specified
00448    // name is not found. This method overrides the generic FindObject()
00449    // of TCollection for efficiency reasons.
00450 
00451    if (!name) return 0;
00452    TObjLink *lnk = FirstLink();
00453    while (lnk) {
00454       TObject *obj = lnk->GetObject();
00455       const char *objname = obj->GetName();
00456       if (objname && !strcmp(name, objname)) return obj;
00457       lnk = lnk->Next();
00458    }
00459    return 0;
00460 }
00461 
00462 //______________________________________________________________________________
00463 TObject *TList::FindObject(const TObject *obj) const
00464 {
00465    // Find an object in this list using the object's IsEqual()
00466    // member function. Requires a sequential scan till the object has
00467    // been found. Returns 0 if object is not found.
00468    // This method overrides the generic FindObject() of TCollection for
00469    // efficiency reasons.
00470 
00471    TObjLink *lnk = FirstLink();
00472 
00473    while (lnk) {
00474       TObject *ob = lnk->GetObject();
00475       if (ob->IsEqual(obj)) return ob;
00476       lnk = lnk->Next();
00477    }
00478    return 0;
00479 }
00480 
00481 //______________________________________________________________________________
00482 TObjLink *TList::FindLink(const TObject *obj, Int_t &idx) const
00483 {
00484    // Returns the TObjLink object that contains object obj. In idx it returns
00485    // the position of the object in the list.
00486 
00487    if (!fFirst) return 0;
00488 
00489    TObject *object;
00490    TObjLink *lnk = fFirst;
00491    idx = 0;
00492 
00493    while (lnk) {
00494       object = lnk->GetObject();
00495       if (object) {
00496          if (object->TestBit(kNotDeleted)) {
00497             if (object->IsEqual(obj)) return lnk;
00498          }
00499       }
00500       lnk = lnk->Next();
00501       idx++;
00502    }
00503    return 0;
00504 }
00505 
00506 //______________________________________________________________________________
00507 TObject *TList::First() const
00508 {
00509    // Return the first object in the list. Returns 0 when list is empty.
00510 
00511    if (fFirst) return fFirst->GetObject();
00512    return 0;
00513 }
00514 
00515 //______________________________________________________________________________
00516 TObject **TList::GetObjectRef(const TObject *obj) const
00517 {
00518    // Return address of pointer to obj
00519 
00520    TObjLink *lnk = FirstLink();
00521 
00522    while (lnk) {
00523       TObject *ob = lnk->GetObject();
00524       if (ob->IsEqual(obj)) return lnk->GetObjectRef();
00525       lnk = lnk->Next();
00526    }
00527    return 0;
00528 }
00529 
00530 //______________________________________________________________________________
00531 TObject *TList::Last() const
00532 {
00533    // Return the last object in the list. Returns 0 when list is empty.
00534 
00535    if (fLast) return fLast->GetObject();
00536    return 0;
00537 }
00538 
00539 //______________________________________________________________________________
00540 TObjLink *TList::LinkAt(Int_t idx) const
00541 {
00542    // Return the TObjLink object at index idx.
00543 
00544    Int_t    i = 0;
00545    TObjLink *lnk = fFirst;
00546    while (i < idx && lnk) {
00547       i++;
00548       lnk = lnk->Next();
00549    }
00550    return lnk;
00551 }
00552 
00553 //______________________________________________________________________________
00554 TIterator *TList::MakeIterator(Bool_t dir) const
00555 {
00556    // Return a list iterator.
00557 
00558    return new TListIter(this, dir);
00559 }
00560 
00561 //______________________________________________________________________________
00562 TObjLink *TList::NewLink(TObject *obj, TObjLink *prev)
00563 {
00564    // Return a new TObjLink.
00565 
00566    if (prev)
00567       return new TObjLink(obj, prev);
00568    else
00569       return new TObjLink(obj);
00570 }
00571 
00572 //______________________________________________________________________________
00573 TObjLink *TList::NewOptLink(TObject *obj, Option_t *opt, TObjLink *prev)
00574 {
00575    // Return a new TObjOptLink (a TObjLink that also stores the option).
00576 
00577    if (prev)
00578       return new TObjOptLink(obj, prev, opt);
00579    else
00580       return new TObjOptLink(obj, opt);
00581 }
00582 
00583 //______________________________________________________________________________
00584 void TList::RecursiveRemove(TObject *obj)
00585 {
00586    // Remove object from this collection and recursively remove the object
00587    // from all other objects (and collections).
00588 
00589    if (!obj) return;
00590 
00591    TObjLink *lnk  = fFirst;
00592    TObjLink *next = 0;
00593    while (lnk) {
00594       next = lnk->Next();
00595       TObject *ob = lnk->GetObject();
00596       if (ob->TestBit(kNotDeleted)) {
00597          if (ob->IsEqual(obj)) {
00598             if (lnk == fFirst) {
00599                fFirst = next;
00600                if (lnk == fLast)
00601                   fLast = fFirst;
00602                else
00603                   fFirst->fPrev = 0;
00604                DeleteLink(lnk);
00605             } else if (lnk == fLast) {
00606                fLast = lnk->Prev();
00607                fLast->fNext = 0;
00608                DeleteLink(lnk);
00609             } else {
00610                lnk->Prev()->fNext = next;
00611                lnk->Next()->fPrev = lnk->Prev();
00612                DeleteLink(lnk);
00613             }
00614             fSize--;
00615             fCache = 0;
00616             Changed();
00617          } else
00618             ob->RecursiveRemove(obj);
00619       }
00620       lnk = next;
00621    }
00622 }
00623 
00624 //______________________________________________________________________________
00625 TObject *TList::Remove(TObject *obj)
00626 {
00627    // Remove object from the list.
00628 
00629    if (!obj) return 0;
00630 
00631    Int_t    idx;
00632    TObjLink *lnk = FindLink(obj, idx);
00633 
00634    if (!lnk) return 0;
00635 
00636    // return object found, which may be (pointer wise) different than the
00637    // input object (depending on what IsEqual() is doing)
00638    TObject *ob = lnk->GetObject();
00639 
00640    if (lnk == fFirst) {
00641       fFirst = lnk->Next();
00642       if (lnk == fLast)
00643          fLast = fFirst;
00644       else
00645          fFirst->fPrev = 0;
00646       DeleteLink(lnk);
00647    } else if (lnk == fLast) {
00648       fLast = lnk->Prev();
00649       fLast->fNext = 0;
00650       DeleteLink(lnk);
00651    } else {
00652       lnk->Prev()->fNext = lnk->Next();
00653       lnk->Next()->fPrev = lnk->Prev();
00654       DeleteLink(lnk);
00655    }
00656    fSize--;
00657    fCache = 0;
00658    Changed();
00659 
00660    return ob;
00661 }
00662 
00663 //______________________________________________________________________________
00664 TObject *TList::Remove(TObjLink *lnk)
00665 {
00666    // Remove object link (and therefore the object it contains)
00667    // from the list.
00668 
00669    if (!lnk) return 0;
00670 
00671    TObject *obj = lnk->GetObject();
00672 
00673    if (lnk == fFirst) {
00674       fFirst = lnk->Next();
00675       if (lnk == fLast)
00676          fLast = fFirst;
00677       else
00678          fFirst->fPrev = 0;
00679       DeleteLink(lnk);
00680    } else if (lnk == fLast) {
00681       fLast = lnk->Prev();
00682       fLast->fNext = 0;
00683       DeleteLink(lnk);
00684    } else {
00685       lnk->Prev()->fNext = lnk->Next();
00686       lnk->Next()->fPrev = lnk->Prev();
00687       DeleteLink(lnk);
00688    }
00689    fSize--;
00690    fCache = 0;
00691    Changed();
00692 
00693    return obj;
00694 }
00695 
00696 //______________________________________________________________________________
00697 void TList::RemoveLast()
00698 {
00699    // Remove the last object of the list.
00700 
00701    TObjLink *lnk = fLast;
00702    if (!lnk) return;
00703    
00704    if (lnk == fFirst) {
00705       fFirst = 0;
00706       fLast = 0;
00707    } else {
00708       fLast = lnk->Prev();
00709       fLast->fNext = 0;
00710    }
00711    DeleteLink(lnk);
00712    
00713    fSize--;
00714    fCache = 0;
00715    Changed();   
00716 }
00717    
00718 //______________________________________________________________________________
00719 void TList::Sort(Bool_t order)
00720 {
00721    // Sort linked list. Real sorting is done in private function DoSort().
00722    // The list can only be sorted when is contains objects of a sortable
00723    // class.
00724 
00725    if (!fFirst) return;
00726 
00727    fAscending = order;
00728 
00729    if (!fFirst->GetObject()->IsSortable()) {
00730       Error("Sort", "objects in list are not sortable");
00731       return;
00732    }
00733 
00734    DoSort(&fFirst, fSize);
00735 
00736    // correct back links
00737    TObjLink *ol, *lnk = fFirst;
00738 
00739    if (lnk) lnk->fPrev = 0;
00740    while ((ol = lnk)) {
00741       lnk = lnk->fNext;
00742       if (lnk)
00743          lnk->fPrev = ol;
00744       else
00745          fLast = ol;
00746    }
00747    fSorted = kTRUE;
00748 }
00749 
00750 //______________________________________________________________________________
00751 Bool_t TList::LnkCompare(TObjLink *l1, TObjLink *l2)
00752 {
00753    // Compares the objects stored in the TObjLink objects.
00754    // Depending on the flag IsAscending() the function returns
00755    // true if the object in l1 <= l2 (ascending) or l2 <= l1 (descending).
00756 
00757    Int_t cmp = l1->GetObject()->Compare(l2->GetObject());
00758 
00759    if ((IsAscending() && cmp <=0) || (!IsAscending() && cmp > 0))
00760       return kTRUE;
00761    return kFALSE;
00762 }
00763 
00764 //______________________________________________________________________________
00765 TObjLink **TList::DoSort(TObjLink **head, Int_t n)
00766 {
00767    // Sort linked list.
00768 
00769    TObjLink *p1, *p2, **h2, **t2;
00770 
00771    switch (n) {
00772       case 0:
00773          return head;
00774 
00775       case 1:
00776          return &((*head)->fNext);
00777 
00778       case 2:
00779          p2 = (p1 = *head)->fNext;
00780          if (LnkCompare(p1, p2)) return &(p2->fNext);
00781          p1->fNext = (*head = p2)->fNext;
00782          return &((p2->fNext = p1)->fNext);
00783    }
00784 
00785    int m;
00786    n -= m = n / 2;
00787 
00788    t2 = DoSort(h2 = DoSort(head, n), m);
00789 
00790    if (LnkCompare((p1 = *head), (p2 = *h2))) {
00791       do {
00792          if (!--n) return *h2 = p2, t2;
00793       } while (LnkCompare((p1 = *(head = &(p1->fNext))), p2));
00794    }
00795 
00796    while (1) {
00797       *head = p2;
00798       do {
00799          if (!--m) return *h2 = *t2, *t2 = p1, h2;
00800       } while (!LnkCompare(p1, (p2 = *(head = &(p2->fNext)))));
00801       *head = p1;
00802       do {
00803          if (!--n) return *h2 = p2, t2;
00804       } while (LnkCompare((p1 = *(head = &(p1->fNext))), p2));
00805    }
00806 }
00807 
00808 //////////////////////////////////////////////////////////////////////////
00809 //                                                                      //
00810 // TObjLink                                                             //
00811 //                                                                      //
00812 // Wrapper around a TObject so it can be stored in a TList.             //
00813 //                                                                      //
00814 //////////////////////////////////////////////////////////////////////////
00815 
00816 //______________________________________________________________________________
00817 TObjLink::TObjLink(TObject *obj, TObjLink *prev)
00818           : fNext(prev->fNext), fPrev(prev), fObject(obj)
00819 {
00820    // Create a new TObjLink.
00821 
00822    fPrev->fNext = this;
00823    if (fNext) fNext->fPrev = this;
00824 }
00825 
00826 //////////////////////////////////////////////////////////////////////////
00827 //                                                                      //
00828 // TListIter                                                            //
00829 //                                                                      //
00830 // Iterator of linked list.                                             //
00831 //                                                                      //
00832 //////////////////////////////////////////////////////////////////////////
00833 
00834 ClassImp(TListIter)
00835 
00836 //______________________________________________________________________________
00837 TListIter::TListIter(const TList *l, Bool_t dir)
00838         : fList(l), fCurCursor(0), fCursor(0), fDirection(dir), fStarted(kFALSE)
00839 {
00840    // Create a new list iterator. By default the iteration direction
00841    // is kIterForward. To go backward use kIterBackward.
00842 }
00843 
00844 //______________________________________________________________________________
00845 TListIter::TListIter(const TListIter &iter) : TIterator(iter)
00846 {
00847    // Copy ctor.
00848 
00849    fList      = iter.fList;
00850    fCurCursor = iter.fCurCursor;
00851    fCursor    = iter.fCursor;
00852    fDirection = iter.fDirection;
00853    fStarted   = iter.fStarted;
00854 }
00855 
00856 //______________________________________________________________________________
00857 TIterator &TListIter::operator=(const TIterator &rhs)
00858 {
00859    // Overridden assignment operator.
00860 
00861    if (this != &rhs && rhs.IsA() == TListIter::Class()) {
00862       const TListIter &rhs1 = (const TListIter &)rhs;
00863       fList      = rhs1.fList;
00864       fCurCursor = rhs1.fCurCursor;
00865       fCursor    = rhs1.fCursor;
00866       fDirection = rhs1.fDirection;
00867       fStarted   = rhs1.fStarted;
00868    }
00869    return *this;
00870 }
00871 
00872 //______________________________________________________________________________
00873 TListIter &TListIter::operator=(const TListIter &rhs)
00874 {
00875    // Overloaded assignment operator.
00876 
00877    if (this != &rhs) {
00878       fList      = rhs.fList;
00879       fCurCursor = rhs.fCurCursor;
00880       fCursor    = rhs.fCursor;
00881       fDirection = rhs.fDirection;
00882       fStarted   = rhs.fStarted;
00883    }
00884    return *this;
00885 }
00886 
00887 //______________________________________________________________________________
00888 TObject *TListIter::Next()
00889 {
00890    // Return next object in the list. Returns 0 when no more objects in list.
00891 
00892    if (!fList) return 0;
00893 
00894    if (fDirection == kIterForward) {
00895       if (!fStarted) {
00896          fCursor = fList->fFirst;
00897          fStarted = kTRUE;
00898       }
00899       fCurCursor = fCursor;
00900       if (fCursor) fCursor = fCursor->Next();
00901    } else {
00902       if (!fStarted) {
00903          fCursor = fList->fLast;
00904          fStarted = kTRUE;
00905       }
00906       fCurCursor = fCursor;
00907       if (fCursor) fCursor = fCursor->Prev();
00908    }
00909 
00910    if (fCurCursor) return fCurCursor->GetObject();
00911    return 0;
00912 }
00913 
00914 //______________________________________________________________________________
00915 Option_t *TListIter::GetOption() const
00916 {
00917    // Returns the object option stored in the list.
00918 
00919    if (fCurCursor) return fCurCursor->GetOption();
00920    return "";
00921 }
00922 
00923 //______________________________________________________________________________
00924 void TListIter::SetOption(Option_t *option)
00925 {
00926    // Sets the object option stored in the list.
00927 
00928    if (fCurCursor) fCurCursor->SetOption(option);
00929 }
00930 
00931 //______________________________________________________________________________
00932 void TListIter::Reset()
00933 {
00934    // Reset list iterator.
00935 
00936    fStarted = kFALSE;
00937 }
00938 
00939 //______________________________________________________________________________
00940 bool TListIter::operator!=(const TIterator &aIter) const
00941 {
00942    // This operator compares two TIterator objects.
00943 
00944    if (nullptr == (&aIter))
00945       return fCurCursor;
00946 
00947    if ((aIter.IsA() == TListIter::Class())) {
00948       const TListIter &iter(dynamic_cast<const TListIter &>(aIter));
00949       return (fCurCursor != iter.fCurCursor);
00950    }
00951    return false; // for base class we don't implement a comparison
00952 }
00953 
00954 //______________________________________________________________________________
00955 bool TListIter::operator!=(const TListIter &aIter) const
00956 {
00957    // This operator compares two TListIter objects.
00958 
00959    if (nullptr == (&aIter))
00960       return fCurCursor;
00961 
00962    return (fCurCursor != aIter.fCurCursor);
00963 }
00964 
00965 //_______________________________________________________________________
00966 void TList::Streamer(TBuffer &b)
00967 {
00968    // Stream all objects in the collection to or from the I/O buffer.
00969 
00970    Int_t nobjects;
00971    UChar_t nch;
00972    Int_t nbig;
00973    TObject *obj;
00974    UInt_t R__s, R__c;
00975 
00976    if (b.IsReading()) {
00977       Version_t v = b.ReadVersion(&R__s, &R__c);
00978       if (v > 3) {
00979          TObject::Streamer(b);
00980          fName.Streamer(b);
00981          b >> nobjects;
00982          string readOption;
00983          for (Int_t i = 0; i < nobjects; i++) {
00984             b >> obj;
00985             b >> nch;
00986             if (v > 4 && nch == 255)  {
00987                b >> nbig;
00988             } else {
00989                nbig = nch;
00990             }
00991             readOption.resize(nbig,'\0');
00992             b.ReadFastArray((char*) readOption.data(),nbig);
00993             if (obj) { // obj can be null if the class had a custom streamer and we do not have the shared library nor a streamerInfo.
00994                if (nch) {
00995                   Add(obj,readOption.c_str());
00996                } else {
00997                   Add(obj);
00998                }
00999             }
01000          }
01001          b.CheckByteCount(R__s, R__c,TList::IsA());
01002          return;
01003       }
01004 
01005       //  process old versions when TList::Streamer was in TCollection::Streamer
01006       if (v > 2)
01007          TObject::Streamer(b);
01008       if (v > 1)
01009          fName.Streamer(b);
01010       b >> nobjects;
01011       for (Int_t i = 0; i < nobjects; i++) {
01012          b >> obj;
01013          Add(obj);
01014       }
01015       b.CheckByteCount(R__s, R__c,TList::IsA());
01016 
01017    } else {
01018       R__c = b.WriteVersion(TList::IsA(), kTRUE);
01019       TObject::Streamer(b);
01020       fName.Streamer(b);
01021       nobjects = GetSize();
01022       b << nobjects;
01023 
01024       TObjLink *lnk = fFirst;
01025       while (lnk) {
01026          obj = lnk->GetObject();
01027          b << obj;
01028 
01029          nbig = strlen(lnk->GetAddOption());
01030          if (nbig > 254) {
01031             nch = 255;
01032             b << nch;
01033             b << nbig;
01034          } else {
01035             nch = UChar_t(nbig);
01036             b << nch;
01037          }
01038          b.WriteFastArray(lnk->GetAddOption(),nbig);
01039 
01040          lnk = lnk->Next();
01041       }
01042       b.SetByteCount(R__c, kTRUE);
01043    }
01044 }

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