Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TEventList.cxx
Go to the documentation of this file.
1// @(#)root/tree:$Id$
2// Author: Rene Brun 11/02/97
3
4/*************************************************************************
5 * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers. *
6 * All rights reserved. *
7 * *
8 * For the licensing terms see $ROOTSYS/LICENSE. *
9 * For the list of contributors see $ROOTSYS/README/CREDITS. *
10 *************************************************************************/
11
12/** \class TEventList
13\ingroup tree
14
15\legacy{TEventList}
16
17A TEventList object is a list of selected events (entries) in a TTree.
18
19A TEventList is automatically generated by TTree::Draw: example
20~~~ {.cpp}
21 tree->Draw(">>elist1","x<0 && y> 0");
22~~~
23In this example, a TEventList object named "elist1" will be
24generated. (Previous contents are overwritten).
25~~~ {.cpp}
26 tree->Draw(">>+elist1","x<0 && y> 0");
27~~~
28In this example, selected entries are added to the list.
29
30The TEventList object is added to the list of objects in the current
31directory.
32
33Use TTree:SetEventList(TEventList *list) to inform TTree that you
34want to use the list as input. The following code gets a pointer to
35the TEventList object created in the above commands:
36~~~ {.cpp}
37 TEventList *list = (TEventList*)gDirectory->Get("elist1");
38~~~
39- Use function Enter to enter an element in the list
40- The function Add may be used to merge two lists.
41- The function Subtract may be used to subtract two lists.
42- The function Reset may be used to reset a list.
43- Use TEventList::Print(option) to print the contents.
44 (option "all" prints all the list entries).
45- Operators + and - correspond to functions Add and Subtract.
46- A TEventList object can be saved on a file via the Write function.
47*/
48
49#include "TEventList.h"
50#include "TBuffer.h"
51#include "TCut.h"
52#include "TDirectory.h"
53#include "TList.h"
54#include "TMath.h"
55#include "strlcpy.h"
56#include "snprintf.h"
57
58
59////////////////////////////////////////////////////////////////////////////////
60/// Default constructor for a EventList.
61
63{
64 fN = 0;
65 fSize = 100;
66 fDelta = 100;
67 fList = nullptr;
68 fDirectory = nullptr;
69 fReapply = false;
70}
71
72////////////////////////////////////////////////////////////////////////////////
73/// Create a EventList.
74///
75/// This Eventlist is added to the list of objects in current directory.
76
77TEventList::TEventList(const char *name, const char *title, Int_t initsize, Int_t delta)
78 :TNamed(name,title), fReapply(false)
79{
80 fN = 0;
81 if (initsize > 100) fSize = initsize;
82 else fSize = 100;
83 if (delta > 100) fDelta = delta;
84 else fDelta = 100;
85 fList = nullptr;
87 if (fDirectory) fDirectory->Append(this);
88}
89
90////////////////////////////////////////////////////////////////////////////////
91/// Copy constructor.
92
94{
95 fN = list.fN;
96 fSize = list.fSize;
97 fDelta = list.fDelta;
98 fList = new Long64_t[fSize];
99 for (Int_t i=0; i<fN; i++)
100 fList[i] = list.fList[i];
101 fReapply = list.fReapply;
102 fDirectory = nullptr;
103}
104
105////////////////////////////////////////////////////////////////////////////////
106/// Default destructor for a EventList.
107
109{
110 delete [] fList; fList = nullptr;
111 if (fDirectory) fDirectory->Remove(this);
112 fDirectory = nullptr;
113}
114
115////////////////////////////////////////////////////////////////////////////////
116/// Merge contents of alist with this list.
117///
118/// Both alist and this list are assumed to be sorted prior to this call
119
120void TEventList::Add(const TEventList *alist)
121{
122 Int_t i;
123 Int_t an = alist->GetN();
124 if (!an) return;
125 Long64_t *alst = alist->GetList();
126 if (!fList) {
127 fList = new Long64_t[an];
128 for (i=0;i<an;i++) fList[i] = alst[i];
129 fN = an;
130 fSize = an;
131 return;
132 }
133 Int_t newsize = fN + an;
136 newpos = alpos = 0;
137 for (i=0;i<fN;i++) {
138 while (alpos < an && fList[i] > alst[alpos]) {
140 newpos++;
141 alpos++;
142 }
143 if (alpos < an && fList[i] == alst[alpos]) alpos++;
144 newlist[newpos] = fList[i];
145 newpos++;
146 }
147 while (alpos < an) {
149 newpos++;
150 alpos++;
151 }
152 delete [] fList;
153 fN = newpos;
154 fSize = newsize;
155 fList = newlist;
156
157 TCut orig = GetTitle();
158 TCut added = alist->GetTitle();
159 TCut updated = orig || added;
160 SetTitle(updated.GetTitle());
161}
162
163////////////////////////////////////////////////////////////////////////////////
164/// Return TRUE if list contains entry.
165
167{
168 if (GetIndex(entry) < 0) return false;
169 return true;
170}
171
172////////////////////////////////////////////////////////////////////////////////
173/// Return TRUE if list contains entries from entrymin to entrymax included.
174
176{
178 //printf("ContainsRange: entrymin=%lld, entrymax=%lld,imax=%lld, fList[imax]=%lld\n",entrymin,entrymax,imax,fList[imax]);
179
180 if (fList[imax] < entrymin) return false;
181 return true;
182}
183
184////////////////////////////////////////////////////////////////////////////////
185/// Called by TKey and others to automatically add us to a directory when we are read from a file.
186
191
192////////////////////////////////////////////////////////////////////////////////
193/// Enter element entry into the list.
194
196{
197 if (!fList) {
198 fList = new Long64_t[fSize];
199 fList[0] = entry;
200 fN = 1;
201 return;
202 }
203 if (fN>0 && entry==fList[fN-1]) return;
204 if (fN >= fSize) {
207 }
208 if(fN==0 || entry>fList[fN-1]) {
209 fList[fN] = entry;
210 ++fN;
211 } else {
213 if(pos>=0 && entry==fList[pos])
214 return;
215 ++pos;
216 memmove( &(fList[pos+1]), &(fList[pos]), 8*(fN-pos));
217 fList[pos] = entry;
218 ++fN;
219 }
220}
221
222////////////////////////////////////////////////////////////////////////////////
223/// Return value of entry at index in the list.
224/// Return -1 if index is not in the list range.
225
227{
228 if (!fList) return -1;
229 if (index < 0 || index >= fN) return -1;
230 return fList[index];
231}
232
233////////////////////////////////////////////////////////////////////////////////
234/// Return index in the list of element with value entry
235/// array is supposed to be sorted prior to this call.
236/// If match is found, function returns position of element.
237/// If no match found, function returns -1.
238
240{
242 nabove = fN+1;
243 nbelow = 0;
244 while(nabove-nbelow > 1) {
245 middle = (nabove+nbelow)/2;
246 if (entry == fList[middle-1]) return middle-1;
247 if (entry < fList[middle-1]) nabove = middle;
248 else nbelow = middle;
249 }
250 return -1;
251}
252
253////////////////////////////////////////////////////////////////////////////////
254/// Remove elements from this list that are NOT present in alist.
255
257{
258 if (!alist) return;
259 if (!fList) return;
260
261 Long64_t *newlist = new Long64_t[fN];
262 Int_t newpos = 0;
263 Int_t i;
264 for (i=0;i<fN;i++) {
265 if (alist->GetIndex(fList[i]) >= 0) {
266 newlist[newpos] = fList[i];
267 newpos++;
268 }
269 }
270 delete [] fList;
271 fN = newpos;
272 fList = newlist;
273
274 TCut orig = GetTitle();
275 TCut removed = alist->GetTitle();
277 SetTitle(updated.GetTitle());
278}
279
280////////////////////////////////////////////////////////////////////////////////
281/// Merge entries in all the TEventList in the collection in this event list.
282
284{
285 if (!list) return -1;
286 TIter next(list);
287
288 //first loop to count the number of entries
289 TEventList *el;
290 Int_t nevents = 0;
291 while ((el = (TEventList*)next())) {
292 if (!el->InheritsFrom(TEventList::Class())) {
293 Error("Add","Attempt to add object of class: %s to a %s",el->ClassName(),this->ClassName());
294 return -1;
295 }
296 Add(el);
297 nevents += el->GetN();
298 }
299
300 return nevents;
301}
302
303////////////////////////////////////////////////////////////////////////////////
304/// Print contents of this list.
305
307{
308 printf("EventList:%s/%s, number of entries =%d, size=%d\n",GetName(),GetTitle(),fN,fSize);
309 if (!strstr(option,"all")) return;
310 Int_t i;
311 Int_t nbuf = 0;
312 char element[10];
313 char *line = new char[100];
314 snprintf(line,100,"%5d : ",0);
315 for (i=0;i<fN;i++) {
316 nbuf++;
317 if (nbuf > 10) {
318 printf("%s\n",line);
319 snprintf(line,100,"%5d : ",i);
320 nbuf = 1;
321 }
322 snprintf(element,10,"%7lld ",fList[i]);
323 strlcat(line,element,100);
324 }
325 if (nbuf) printf("%s\n",line);
326 delete [] line;
327}
328
329////////////////////////////////////////////////////////////////////////////////
330/// Reset number of entries in event list.
331
333{
334 fN = 0;
335}
336
337////////////////////////////////////////////////////////////////////////////////
338/// Resize list by delta entries.
339
341{
342 if (!delta) delta = fDelta;
343 fSize += delta;
345 for (Int_t i=0;i<fN;i++) newlist[i] = fList[i];
346 delete [] fList;
347 fList = newlist;
348}
349
350////////////////////////////////////////////////////////////////////////////////
351/// Remove reference to this EventList from current directory and add
352/// reference to new directory dir. dir can be 0 in which case the list
353/// does not belong to any directory.
354
356{
357 if (fDirectory == dir) return;
358 if (fDirectory) fDirectory->Remove(this);
359 fDirectory = dir;
360 if (fDirectory) fDirectory->Append(this);
361}
362
363////////////////////////////////////////////////////////////////////////////////
364/// Change the name of this TEventList.
365
366void TEventList::SetName(const char *name)
367{
368 // TEventLists are named objects in a THashList.
369 // We must update the hashlist if we change the name
370 if (fDirectory) fDirectory->Remove(this);
371 fName = name;
372 if (fDirectory) fDirectory->Append(this);
373}
374
375////////////////////////////////////////////////////////////////////////////////
376/// Sort list entries in increasing order
377
379{
380 Int_t *index = new Int_t[fN];
382 Int_t i,ind;
383 TMath::Sort(fN,fList,index); //sort in decreasing order
384 for (i=0;i<fN;i++) {
385 ind = index[fN-i-1];
386 newlist[i] = fList[ind];
387 }
388 for (i=fN;i<fSize;i++) {
389 newlist[i] = 0;
390 }
391 delete [] index;
392 delete [] fList;
393 fList = newlist;
394}
395
396////////////////////////////////////////////////////////////////////////////////
397/// Stream an object of class TEventList.
398
400{
401 if (b.IsReading()) {
403 Version_t R__v = b.ReadVersion(&R__s, &R__c);
404 fDirectory = nullptr;
405 if (R__v > 1) {
406 b.ReadClassBuffer(TEventList::Class(), this, R__v, R__s, R__c);
408 return;
409 }
410 //====process old versions before automatic schema evolution
412 b >> fN;
413 b >> fSize;
414 b >> fDelta;
415 if (fN) {
416 Int_t *tlist = new Int_t[fSize];
417 b.ReadFastArray(tlist,fN);
418 fList = new Long64_t[fSize];
419 for (Int_t i=0;i<fN;i++) fList[i] = tlist[i];
420 delete [] tlist;
421 }
423 b.CheckByteCount(R__s, R__c, TEventList::IsA());
424 //====end of old versions
425
426 } else {
427 b.WriteClassBuffer(TEventList::Class(), this);
428 }
429}
430
431////////////////////////////////////////////////////////////////////////////////
432/// Remove elements from this list that are present in alist.
433
435{
436 if (!alist) return;
437 if (!fList) return;
438
439 Long64_t *newlist = new Long64_t[fN];
440 Int_t newpos = 0;
441 Int_t i;
442 for (i=0;i<fN;i++) {
443 if (alist->GetIndex(fList[i]) < 0) {
444 newlist[newpos] = fList[i];
445 newpos++;
446 }
447 }
448 delete [] fList;
449 fN = newpos;
450 fList = newlist;
451
452 TCut orig = GetTitle();
453 TCut removed = alist->GetTitle();
454 TCut updated = orig && !removed;
455 SetTitle(updated.GetTitle());
456}
457
458////////////////////////////////////////////////////////////////////////////////
459/// Assingment.
460
462{
463 if (this != &list) {
464 TNamed::operator=(list);
465 if (fSize < list.fSize) {
466 delete [] fList;
467 fList = new Long64_t[list.fSize];
468 }
469 fN = list.fN;
470 fSize = list.fSize;
471 fDelta = list.fDelta;
472 for (Int_t i=0; i<fN; i++)
473 fList[i] = list.fList[i];
474 }
475 return *this;
476}
477
478////////////////////////////////////////////////////////////////////////////////
479/// Addition.
480
482{
483 TEventList newlist = list1;
484 newlist.Add(&list2);
485 return newlist;
486}
487
488////////////////////////////////////////////////////////////////////////////////
489/// Subtraction
490
492{
493 TEventList newlist = list1;
494 newlist.Subtract(&list2);
495 return newlist;
496}
497
498////////////////////////////////////////////////////////////////////////////////
499/// Intersection.
500
502{
503 TEventList newlist = list1;
504 newlist.Intersect(&list2);
505 return newlist;
506}
507
#define b(i)
Definition RSha256.hxx:100
short Version_t
Class version identifier (short)
Definition RtypesCore.h:79
long long Long64_t
Portable signed long integer 8 bytes.
Definition RtypesCore.h:83
const char Option_t
Option string (const char)
Definition RtypesCore.h:80
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
#define gDirectory
Definition TDirectory.h:385
TEventList operator-(const TEventList &list1, const TEventList &list2)
Subtraction.
TEventList operator+(const TEventList &list1, const TEventList &list2)
Addition.
TEventList operator*(const TEventList &list1, const TEventList &list2)
Intersection.
Option_t Option_t option
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t index
char name[80]
Definition TGX11.cxx:110
#define snprintf
Definition civetweb.c:1579
Buffer base class used for serializing objects.
Definition TBuffer.h:43
Collection abstract base class.
Definition TCollection.h:65
A specialized string object used for TTree selections.
Definition TCut.h:25
Describe directory structure in memory.
Definition TDirectory.h:45
virtual void Append(TObject *obj, Bool_t replace=kFALSE)
Append object to this directory.
virtual TObject * Remove(TObject *)
Remove an object from the in-memory list.
<div class="legacybox"><h2>Legacy Code</h2> TEventList is a legacy interface: there will be no bug fi...
Definition TEventList.h:31
Long64_t * fList
[fN]Array of elements
Definition TEventList.h:38
static TClass * Class()
~TEventList() override
Default destructor for a EventList.
void Streamer(TBuffer &) override
Stream an object of class TEventList.
bool fReapply
If true, TTree::Draw will 'reapply' the original cut.
Definition TEventList.h:37
TEventList()
Default constructor for a EventList.
virtual void Reset(Option_t *option="")
Reset number of entries in event list.
virtual bool ContainsRange(Long64_t entrymin, Long64_t entrymax)
Return TRUE if list contains entries from entrymin to entrymax included.
TDirectory * fDirectory
! Pointer to directory holding this tree
Definition TEventList.h:39
virtual Long64_t GetEntry(Int_t index) const
Return value of entry at index in the list.
virtual Int_t GetIndex(Long64_t entry) const
Return index in the list of element with value entry array is supposed to be sorted prior to this cal...
virtual Int_t GetN() const
Definition TEventList.h:56
virtual Long64_t * GetList() const
Definition TEventList.h:55
TClass * IsA() const override
Definition TEventList.h:77
virtual Int_t Merge(TCollection *list)
Merge entries in all the TEventList in the collection in this event list.
virtual void DirectoryAutoAdd(TDirectory *)
Called by TKey and others to automatically add us to a directory when we are read from a file.
virtual void Subtract(const TEventList *list)
Remove elements from this list that are present in alist.
void Print(Option_t *option="") const override
Print contents of this list.
virtual void SetDirectory(TDirectory *dir)
Remove reference to this EventList from current directory and add reference to new directory dir.
virtual void Add(const TEventList *list)
Merge contents of alist with this list.
Int_t fN
Number of elements in the list.
Definition TEventList.h:34
virtual void Resize(Int_t delta=0)
Resize list by delta entries.
virtual void Enter(Long64_t entry)
Enter element entry into the list.
Int_t fSize
Size of array.
Definition TEventList.h:35
virtual void Sort()
Sort list entries in increasing order.
virtual bool Contains(Long64_t entry)
Return TRUE if list contains entry.
TEventList & operator=(const TEventList &list)
Assingment.
virtual void Intersect(const TEventList *list)
Remove elements from this list that are NOT present in alist.
Int_t fDelta
Increment size.
Definition TEventList.h:36
void SetName(const char *name) override
Change the name of this TEventList.
The TNamed class is the base class for all named ROOT classes.
Definition TNamed.h:29
virtual void SetTitle(const char *title="")
Set the title of the TNamed.
Definition TNamed.cxx:173
const char * GetName() const override
Returns name of object.
Definition TNamed.h:49
void Streamer(TBuffer &) override
Stream an object of class TObject.
const char * GetTitle() const override
Returns title of object.
Definition TNamed.h:50
TString fName
Definition TNamed.h:32
TNamed & operator=(const TNamed &rhs)
TNamed assignment operator.
Definition TNamed.cxx:50
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:1071
void ResetBit(UInt_t f)
Definition TObject.h:201
@ kMustCleanup
if object destructor must call RecursiveRemove()
Definition TObject.h:70
TLine * line
Short_t Max(Short_t a, Short_t b)
Returns the largest of a and b.
Definition TMathBase.h:251
void Sort(Index n, const Element *a, Index *index, Bool_t down=kTRUE)
Sort the n elements of the array a of generic templated type Element.
Definition TMathBase.h:432
Long64_t BinarySearch(Long64_t n, const T *array, T value)
Binary search in an array of n values to locate value.
Definition TMathBase.h:348