//*CMZ :  2.22/01 26/04/99  11.51.22  by  Rene Brun
//*CMZ :  2.21/08 17/03/99  18.09.58  by  Rene Brun
//*CMZ :  2.21/01 03/03/99  14.18.31  by  Rene Brun
//*CMZ :  2.20/04 09/12/98  19.25.22  by  Rene Brun
//*CMZ :  2.20/00 04/11/98  12.54.44  by  Rene Brun
//*CMZ :  2.00/11 23/08/98  18.00.10  by  Rene Brun
//*CMZ :  2.00/10 30/06/98  08.02.36  by  Rene Brun
//*CMZ :  2.00/08 02/06/98  13.36.31  by  Rene Brun
//*CMZ :  2.00/00 11/02/98  12.19.22  by  Rene Brun
//*CMZ :  1.03/09 08/12/97  13.41.52  by  Fons Rademakers
//*-- Author :    Rene Brun   11/02/97

//*KEEP,CopyRight,T=C.
/*************************************************************************
 * Copyright(c) 1995-1999, The ROOT System, All rights reserved.         *
 * Authors: Rene Brun, Fons Rademakers.                                  *
 * For list of contributors see $ROOTSYS/AA_CREDITS.                     *
 *                                                                       *
 * Permission to use, copy, modify and distribute this software and its  *
 * documentation for non-commercial purposes is hereby granted without   *
 * fee, provided that the above copyright notice appears in all copies   *
 * and that both the copyright notice and this permission notice appear  *
 * in the supporting documentation. The authors make no claims about the *
 * suitability of this software for any purpose. It is provided "as is"  *
 * without express or implied warranty.                                  *
 *************************************************************************/
//*KEND.

//////////////////////////////////////////////////////////////////////////
//                                                                      //
// A TEventList object is a list of selected events (entries) in a TTree//
//                                                                      //
// A TEventList is automatically generated by TTree::Draw: example      //
//      tree->Draw(">>elist1","x<0 && y> 0");                           //
//         In this example, a TEventList object named "elist1" will be  //
//         generated. (Previous contents are overwritten).              //
//      tree->Draw(">>+elist1","x<0 && y> 0");                          //
//         In this example, selected entries are added to the list.     //
//                                                                      //
// The TEventList object is added to the list of objects in the current //
// directory.                                                           //
// Use TTree:SetEventList(TEventList *list) to inform TTree that you    //
// want to use the list as input. The following code gets a pointer to  //
// the TEventList object created in the above commands:                 //
//     TEventList *list = (TEventList*)gDirectory->Get("elist1");       //
//                                                                      //
// Use function Enter to enter an element in the list                   //
// The function Add may be used to merge two lists.                     //
// The function Subtract may be used to subtract two lists.             //
// The function Reset may be used to reset a list.                      //
// Use TEventList::Print(option) to print the contents.                 //
//      (option "all" prints all the list entries).                     //
//   Operators + and - correspond to functions Add and Subtract.        //
// A TEventList object can be saved on a file via the Write function.   //
//                                                                      //
//////////////////////////////////////////////////////////////////////////

//*KEEP,TEventList,T=C++.
#include "TEventList.h"
//*KEEP,TTree.
#include "TTree.h"
//*KEEP,TFile.
#include "TFile.h"
//*KEEP,TMath.
#include "TMath.h"
//*KEND.

ClassImp(TEventList)

//______________________________________________________________________________
 TEventList::TEventList(): TNamed()
{
//*-*-*-*-*-*Default constructor for a EventList*-*-*-*-*-*-*-*-*-*-*-*-*
//*-*        ==================================

   fN          = 0;
   fSize       = 100;
   fDelta      = 100;
   fList       = 0;
   fDirectory  = 0;
}

//______________________________________________________________________________
 TEventList::TEventList(const Text_t *name, const Text_t *title, Int_t initsize, Int_t delta)
    :TNamed(name,title)
{
//*-*-*-*-*-*-*-*-*-*-*-*-*Create a EventList*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
//*-*                      =================
//
//  This Eventlist is added to the list of objects in current directory

   fN = 0;
   if (initsize > 100) fSize  = initsize;
   else                fSize  = 100;
   if (delta > 100)    fDelta = delta;
   else                fDelta = 100;
   fList       = 0;
   fDirectory  = gDirectory;
   gDirectory->Append(this);
}

//______________________________________________________________________________
 TEventList::TEventList(const TEventList &list)
{
//*-*-*-*-*-*-*-*-*-*-*-*-*Copy constructor*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
//*-*                      ================

   fN     = list.fN;
   fSize  = list.fSize;
   fDelta = list.fDelta;
   fList  = new Int_t[fSize];
   for (Int_t i=0; i<fN; i++)
      fList[i] = list.fList[i];
}

//______________________________________________________________________________
 TEventList::~TEventList()
{
//*-*-*-*-*-*Default destructor for a EventList*-*-*-*-*-*-*-*-*-*-*-*
//*-*        =================================

   delete [] fList;
   if (fDirectory) fDirectory->GetList()->Remove(this);
   fDirectory  = 0;
}

//______________________________________________________________________________
 void TEventList::Add(const TEventList *alist)
{
//          Merge contents of alist with this list.
//   Both alist and this list are assumed to be sorted prior to this call

   Int_t i;
   Int_t an = alist->GetN();
   if (!an) return;
   Int_t *alst = alist->GetList();
   if (!fList) {
      fList = new Int_t[an];
      for (i=0;i<an;i++) fList[i] = alst[i];
      fN = an;
      fSize = an;
      return;
   }
   Int_t newsize = fN + an;
   Int_t *newlist = new Int_t[newsize];
   Int_t newpos, alpos;
   newpos = alpos = 0;
   for (i=0;i<fN;i++) {
      while (alpos < an && fList[i] > alst[alpos]) {
         newlist[newpos] = alst[alpos];
         newpos++;
         alpos++;
      }
      if (fList[i] == alst[alpos]) alpos++;
      newlist[newpos] = fList[i];
      newpos++;
   }
   while (alpos < an) {
      newlist[newpos] = alst[alpos];
      newpos++;
      alpos++;
   }
   delete [] fList;
   fN    = newpos;
   fSize = newsize;
   fList = newlist;
}

//______________________________________________________________________________
 Bool_t TEventList::Contains(Int_t entry)
{
//          Return TRUE if list contains entry.

   if (GetIndex(entry) < 0) return kFALSE;
   return kTRUE;
}

//______________________________________________________________________________
 void TEventList::Enter(Int_t entry)
{
//          Enter element entry into the list

   if (!fList) {
      fList = new Int_t[fSize];
      fList[0] = entry;
      fN = 1;
      return;
   }
   if (fN >= fSize) Resize();
   fList[fN] = entry;
   fN++;
}

//______________________________________________________________________________
 Int_t TEventList::GetEntry(Int_t index) const
{
//       Return value of entry at index in the list.
//       Return -1 if index is not in the list range

   if (!fList)   return -1;
   if (index < 0 || index >= fN)   return -1;
   return fList[index];
}

//______________________________________________________________________________
 Int_t TEventList::GetIndex(Int_t entry) const
{
   // Return index in the list of element with value entry
   // array is supposed  to be sorted prior to this call.
   // If match is found, function returns position of element.
   // If no match found, function returns -1.

   Int_t nabove, nbelow, middle;
   nabove = fN+1;
   nbelow = 0;
   while(nabove-nbelow > 1) {
      middle = (nabove+nbelow)/2;
      if (entry == fList[middle-1]) return middle-1;
      if (entry  < fList[middle-1]) nabove = middle;
      else                          nbelow = middle;
   }
   return -1;
}

//______________________________________________________________________________
 void TEventList::Print(Option_t *option)
{
//          Print contents of this list

   printf("EventList:%s/%s, number of entries =%d, size=%dn",GetName(),GetTitle(),fN,fSize);
   if (!strstr(option,"all")) return;
   Int_t i;
   Int_t nbuf = 0;
   char element[10];
   char *line = new char[100];
   sprintf(line,"%-5d ",0);
   for (i=0;i<fN;i++) {
      nbuf++;
      if (nbuf > 10) {
         printf("%sn",line);
         sprintf(line,"%-5d ",i);
         nbuf = 1;
      }
      sprintf(element,"%-7d ",fList[i]);
      strcat(line,element);
   }
   if (nbuf) printf("%sn",line);
   delete [] line;
}

//______________________________________________________________________________
 void TEventList::Reset(Option_t *)
{
//          Reset number of entries in event list

   fN = 0;
}

//______________________________________________________________________________
 void TEventList::Resize(Int_t delta)
{
//          Resize list by delta entries

   if (!delta) delta = fDelta;
   fSize += delta;
   Int_t *newlist = new Int_t[fSize];
   for (Int_t i=0;i<fN;i++) newlist[i] = fList[i];
   delete [] fList;
   fList = newlist;
}

//______________________________________________________________________________
 void TEventList::SetDirectory(TDirectory *dir)
{
   // Remove reference to this EventList from current directory and add
   // reference to new directory dir. dir can be 0 in which case the list
   // does not belong to any directory.

   if (fDirectory == dir) return;
   if (fDirectory) fDirectory->GetList()->Remove(this);
   fDirectory = dir;
   if (fDirectory) fDirectory->GetList()->Add(this);
}

//______________________________________________________________________________
 void TEventList::Sort()
{
//          Sort list entries in increasing order

   Int_t *index   = new Int_t[fN];
   Int_t *newlist = new Int_t[fSize];
   Int_t i,ind;
   TMath::Sort(fN,fList,index); //sort in decreasing order
   for (i=0;i<fN;i++) {
      ind = index[fN-i-1];
      newlist[i] = fList[ind];
   }
   for (i=fN;i<fSize;i++) {
      newlist[i] = 0;
   }
   delete [] index;
   delete [] fList;
   fList = newlist;
}

//______________________________________________________________________________
 void TEventList::Streamer(TBuffer &b)
{
   // Stream an object of class TEventList.

   if (b.IsReading()) {
      b.ReadVersion();  //Version_t v = b.ReadVersion();
      TNamed::Streamer(b);
      b >> fN;
      b >> fSize;
      b >> fDelta;
      if (fN) {
         fList = new Int_t[fSize];
         b.ReadFastArray(fList,fN);
      }
      fDirectory = gDirectory;
      gDirectory->Append(this);
   } else {
      b.WriteVersion(TEventList::IsA());
      TNamed::Streamer(b);
      b << fN;
      b << fSize;
      b << fDelta;
      if (fN) {
         b.WriteFastArray(fList,fN);
      }
   }
}

//______________________________________________________________________________
 void TEventList::Subtract(const TEventList *alist)
{
   // Remove elements from this list that are present in alist.

   if (!alist) return;
   if (!fList) return;

   Int_t *newlist = new Int_t[fN];
   Int_t newpos = 0;
   Int_t i;
   for (i=0;i<fN;i++) {
      if (alist->GetIndex(fList[i]) < 0) {
         newlist[newpos] = fList[i];
         newpos++;
      }
   }
   delete [] fList;
   fN    = newpos;
   fList = newlist;
}

//______________________________________________________________________________
TEventList& TEventList::operator=(const TEventList &list)
{
   if (this != &list) {
      TNamed::operator=(list);
      if (fSize < list.fSize) {
         delete [] fList;
         fList  = new Int_t[list.fSize];
      }
      fN     = list.fN;
      fSize  = list.fSize;
      fDelta = list.fDelta;
      for (Int_t i=0; i<fN; i++)
         fList[i] = list.fList[i];
   }
   return *this;
}

//______________________________________________________________________________
TEventList operator+(const TEventList &list1, const TEventList &list2)
{
   TEventList newlist = list1;
   newlist.Add(&list2);
   return newlist;
}

//______________________________________________________________________________
TEventList operator-(const TEventList &list1, const TEventList &list2)
{
   TEventList newlist = list1;
   newlist.Subtract(&list2);
   return newlist;
}



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.