#include "TEventList.h"
#include "TCut.h"
#include "TClass.h"
#include "TFile.h"
#include "TMath.h"
ClassImp(TEventList)
TEventList::TEventList(): TNamed()
{
fN = 0;
fSize = 100;
fDelta = 100;
fList = 0;
fDirectory = 0;
fReapply = kFALSE;
}
TEventList::TEventList(const char *name, const char *title, Int_t initsize, Int_t delta)
:TNamed(name,title), fReapply(kFALSE)
{
fN = 0;
if (initsize > 100) fSize = initsize;
else fSize = 100;
if (delta > 100) fDelta = delta;
else fDelta = 100;
fList = 0;
fDirectory = gDirectory;
if (fDirectory) fDirectory->Append(this);
}
TEventList::TEventList(const TEventList &list) : TNamed(list)
{
fN = list.fN;
fSize = list.fSize;
fDelta = list.fDelta;
fList = new Long64_t[fSize];
for (Int_t i=0; i<fN; i++)
fList[i] = list.fList[i];
fReapply = list.fReapply;
fDirectory = 0;
}
TEventList::~TEventList()
{
delete [] fList; fList = 0;
if (fDirectory) fDirectory->Remove(this);
fDirectory = 0;
}
void TEventList::Add(const TEventList *alist)
{
Int_t i;
Int_t an = alist->GetN();
if (!an) return;
Long64_t *alst = alist->GetList();
if (!fList) {
fList = new Long64_t[an];
for (i=0;i<an;i++) fList[i] = alst[i];
fN = an;
fSize = an;
return;
}
Int_t newsize = fN + an;
Long64_t *newlist = new Long64_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;
TCut orig = GetTitle();
TCut added = alist->GetTitle();
TCut updated = orig || added;
SetTitle(updated.GetTitle());
}
Bool_t TEventList::Contains(Long64_t entry)
{
if (GetIndex(entry) < 0) return kFALSE;
return kTRUE;
}
Bool_t TEventList::ContainsRange(Long64_t entrymin, Long64_t entrymax)
{
Long64_t imax = TMath::BinarySearch(fN,fList,entrymax);
if (fList[imax] < entrymin) return kFALSE;
return kTRUE;
}
void TEventList::DirectoryAutoAdd(TDirectory* dir)
{
SetDirectory(dir);
}
void TEventList::Enter(Long64_t entry)
{
if (!fList) {
fList = new Long64_t[fSize];
fList[0] = entry;
fN = 1;
return;
}
if (entry==fList[fN-1]) return;
if (fN >= fSize) {
Int_t newsize = TMath::Max(2*fSize,fN+fDelta);
Resize(newsize-fSize);
}
if(entry>fList[fN-1]) {
fList[fN] = entry;
++fN;
} else {
Int_t pos = TMath::BinarySearch(fN, fList, entry);
if(pos>=0 && entry==fList[pos])
return;
++pos;
memmove( &(fList[pos+1]), &(fList[pos]), 8*(fN-pos));
fList[pos] = entry;
++fN;
}
}
Long64_t TEventList::GetEntry(Int_t index) const
{
if (!fList) return -1;
if (index < 0 || index >= fN) return -1;
return fList[index];
}
Int_t TEventList::GetIndex(Long64_t entry) const
{
Long64_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::Intersect(const TEventList *alist)
{
if (!alist) return;
if (!fList) return;
Long64_t *newlist = new Long64_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;
TCut orig = GetTitle();
TCut removed = alist->GetTitle();
TCut updated = orig && removed;
SetTitle(updated.GetTitle());
}
Int_t TEventList::Merge(TCollection *list)
{
if (!list) return -1;
TIter next(list);
TEventList *el;
Int_t nevents = 0;
while ((el = (TEventList*)next())) {
if (!el->InheritsFrom(TEventList::Class())) {
Error("Add","Attempt to add object of class: %s to a %s",el->ClassName(),this->ClassName());
return -1;
}
Add(el);
nevents += el->GetN();
}
return nevents;
}
void TEventList::Print(Option_t *option) const
{
printf("EventList:%s/%s, number of entries =%d, size=%d\n",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("%s\n",line);
sprintf(line,"%5d : ",i);
nbuf = 1;
}
sprintf(element,"%7lld ",fList[i]);
strcat(line,element);
}
if (nbuf) printf("%s\n",line);
delete [] line;
}
void TEventList::Reset(Option_t *)
{
fN = 0;
}
void TEventList::Resize(Int_t delta)
{
if (!delta) delta = fDelta;
fSize += delta;
Long64_t *newlist = new Long64_t[fSize];
for (Int_t i=0;i<fN;i++) newlist[i] = fList[i];
delete [] fList;
fList = newlist;
}
void TEventList::SetDirectory(TDirectory *dir)
{
if (fDirectory == dir) return;
if (fDirectory) fDirectory->Remove(this);
fDirectory = dir;
if (fDirectory) fDirectory->Append(this);
}
void TEventList::SetName(const char *name)
{
if (fDirectory) fDirectory->Remove(this);
fName = name;
if (fDirectory) fDirectory->Append(this);
}
void TEventList::Sort()
{
Int_t *index = new Int_t[fN];
Long64_t *newlist = new Long64_t[fSize];
Int_t i,ind;
TMath::Sort(fN,fList,index);
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)
{
if (b.IsReading()) {
UInt_t R__s, R__c;
Version_t R__v = b.ReadVersion(&R__s, &R__c);
fDirectory = 0;
if (R__v > 1) {
b.ReadClassBuffer(TEventList::Class(), this, R__v, R__s, R__c);
ResetBit(kMustCleanup);
return;
}
TNamed::Streamer(b);
b >> fN;
b >> fSize;
b >> fDelta;
if (fN) {
Int_t *tlist = new Int_t[fSize];
b.ReadFastArray(tlist,fN);
fList = new Long64_t[fSize];
for (Int_t i=0;i<fN;i++) fList[i] = tlist[i];
delete [] tlist;
}
ResetBit(kMustCleanup);
b.CheckByteCount(R__s, R__c, TEventList::IsA());
} else {
b.WriteClassBuffer(TEventList::Class(), this);
}
}
void TEventList::Subtract(const TEventList *alist)
{
if (!alist) return;
if (!fList) return;
Long64_t *newlist = new Long64_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;
TCut orig = GetTitle();
TCut removed = alist->GetTitle();
TCut updated = orig && !removed;
SetTitle(updated.GetTitle());
}
TEventList& TEventList::operator=(const TEventList &list)
{
if (this != &list) {
TNamed::operator=(list);
if (fSize < list.fSize) {
delete [] fList;
fList = new Long64_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;
}
TEventList operator*(const TEventList &list1, const TEventList &list2)
{
TEventList newlist = list1;
newlist.Intersect(&list2);
return newlist;
}
Last change: Wed Jun 25 08:38:44 2008
Last generated: 2008-06-25 08:38
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.