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
59
60////////////////////////////////////////////////////////////////////////////////
61/// Default constructor for a EventList.
62
64{
65 fN = 0;
66 fSize = 100;
67 fDelta = 100;
68 fList = nullptr;
69 fDirectory = nullptr;
70 fReapply = false;
71}
72
73////////////////////////////////////////////////////////////////////////////////
74/// Create a EventList.
75///
76/// This Eventlist is added to the list of objects in current directory.
77
78TEventList::TEventList(const char *name, const char *title, Int_t initsize, Int_t delta)
79 :TNamed(name,title), fReapply(false)
80{
81 fN = 0;
82 if (initsize > 100) fSize = initsize;
83 else fSize = 100;
84 if (delta > 100) fDelta = delta;
85 else fDelta = 100;
86 fList = nullptr;
88 if (fDirectory) fDirectory->Append(this);
89}
90
91////////////////////////////////////////////////////////////////////////////////
92/// Copy constructor.
93
95{
96 fN = list.fN;
97 fSize = list.fSize;
98 fDelta = list.fDelta;
99 fList = new Long64_t[fSize];
100 for (Int_t i=0; i<fN; i++)
101 fList[i] = list.fList[i];
102 fReapply = list.fReapply;
103 fDirectory = nullptr;
104}
105
106////////////////////////////////////////////////////////////////////////////////
107/// Default destructor for a EventList.
108
110{
111 delete [] fList; fList = nullptr;
112 if (fDirectory) fDirectory->Remove(this);
113 fDirectory = nullptr;
114}
115
116////////////////////////////////////////////////////////////////////////////////
117/// Merge contents of alist with this list.
118///
119/// Both alist and this list are assumed to be sorted prior to this call
120
121void TEventList::Add(const TEventList *alist)
122{
123 Int_t i;
124 Int_t an = alist->GetN();
125 if (!an) return;
126 Long64_t *alst = alist->GetList();
127 if (!fList) {
128 fList = new Long64_t[an];
129 for (i=0;i<an;i++) fList[i] = alst[i];
130 fN = an;
131 fSize = an;
132 return;
133 }
134 Int_t newsize = fN + an;
135 Long64_t *newlist = new Long64_t[newsize];
136 Int_t newpos, alpos;
137 newpos = alpos = 0;
138 for (i=0;i<fN;i++) {
139 while (alpos < an && fList[i] > alst[alpos]) {
140 newlist[newpos] = alst[alpos];
141 newpos++;
142 alpos++;
143 }
144 if (alpos < an && fList[i] == alst[alpos]) alpos++;
145 newlist[newpos] = fList[i];
146 newpos++;
147 }
148 while (alpos < an) {
149 newlist[newpos] = alst[alpos];
150 newpos++;
151 alpos++;
152 }
153 delete [] fList;
154 fN = newpos;
155 fSize = newsize;
156 fList = newlist;
157
158 TCut orig = GetTitle();
159 TCut added = alist->GetTitle();
160 TCut updated = orig || added;
161 SetTitle(updated.GetTitle());
162}
163
164////////////////////////////////////////////////////////////////////////////////
165/// Return TRUE if list contains entry.
166
168{
169 if (GetIndex(entry) < 0) return false;
170 return true;
171}
172
173////////////////////////////////////////////////////////////////////////////////
174/// Return TRUE if list contains entries from entrymin to entrymax included.
175
177{
178 Long64_t imax = TMath::BinarySearch(fN,fList,entrymax);
179 //printf("ContainsRange: entrymin=%lld, entrymax=%lld,imax=%lld, fList[imax]=%lld\n",entrymin,entrymax,imax,fList[imax]);
180
181 if (fList[imax] < entrymin) return false;
182 return true;
183}
184
185////////////////////////////////////////////////////////////////////////////////
186/// Called by TKey and others to automatically add us to a directory when we are read from a file.
187
189{
190 SetDirectory(dir);
191}
192
193////////////////////////////////////////////////////////////////////////////////
194/// Enter element entry into the list.
195
197{
198 if (!fList) {
199 fList = new Long64_t[fSize];
200 fList[0] = entry;
201 fN = 1;
202 return;
203 }
204 if (fN>0 && entry==fList[fN-1]) return;
205 if (fN >= fSize) {
206 Int_t newsize = TMath::Max(2*fSize,fN+fDelta);
207 Resize(newsize-fSize);
208 }
209 if(fN==0 || entry>fList[fN-1]) {
210 fList[fN] = entry;
211 ++fN;
212 } else {
213 Int_t pos = TMath::BinarySearch(fN, fList, entry);
214 if(pos>=0 && entry==fList[pos])
215 return;
216 ++pos;
217 memmove( &(fList[pos+1]), &(fList[pos]), 8*(fN-pos));
218 fList[pos] = entry;
219 ++fN;
220 }
221}
222
223////////////////////////////////////////////////////////////////////////////////
224/// Return value of entry at index in the list.
225/// Return -1 if index is not in the list range.
226
228{
229 if (!fList) return -1;
230 if (index < 0 || index >= fN) return -1;
231 return fList[index];
232}
233
234////////////////////////////////////////////////////////////////////////////////
235/// Return index in the list of element with value entry
236/// array is supposed to be sorted prior to this call.
237/// If match is found, function returns position of element.
238/// If no match found, function returns -1.
239
241{
242 Long64_t nabove, nbelow, middle;
243 nabove = fN+1;
244 nbelow = 0;
245 while(nabove-nbelow > 1) {
246 middle = (nabove+nbelow)/2;
247 if (entry == fList[middle-1]) return middle-1;
248 if (entry < fList[middle-1]) nabove = middle;
249 else nbelow = middle;
250 }
251 return -1;
252}
253
254////////////////////////////////////////////////////////////////////////////////
255/// Remove elements from this list that are NOT present in alist.
256
258{
259 if (!alist) return;
260 if (!fList) return;
261
262 Long64_t *newlist = new Long64_t[fN];
263 Int_t newpos = 0;
264 Int_t i;
265 for (i=0;i<fN;i++) {
266 if (alist->GetIndex(fList[i]) >= 0) {
267 newlist[newpos] = fList[i];
268 newpos++;
269 }
270 }
271 delete [] fList;
272 fN = newpos;
273 fList = newlist;
274
275 TCut orig = GetTitle();
276 TCut removed = alist->GetTitle();
277 TCut updated = orig && removed;
278 SetTitle(updated.GetTitle());
279}
280
281////////////////////////////////////////////////////////////////////////////////
282/// Merge entries in all the TEventList in the collection in this event list.
283
285{
286 if (!list) return -1;
287 TIter next(list);
288
289 //first loop to count the number of entries
290 TEventList *el;
291 Int_t nevents = 0;
292 while ((el = (TEventList*)next())) {
293 if (!el->InheritsFrom(TEventList::Class())) {
294 Error("Add","Attempt to add object of class: %s to a %s",el->ClassName(),this->ClassName());
295 return -1;
296 }
297 Add(el);
298 nevents += el->GetN();
299 }
300
301 return nevents;
302}
303
304////////////////////////////////////////////////////////////////////////////////
305/// Print contents of this list.
306
308{
309 printf("EventList:%s/%s, number of entries =%d, size=%d\n",GetName(),GetTitle(),fN,fSize);
310 if (!strstr(option,"all")) return;
311 Int_t i;
312 Int_t nbuf = 0;
313 char element[10];
314 char *line = new char[100];
315 snprintf(line,100,"%5d : ",0);
316 for (i=0;i<fN;i++) {
317 nbuf++;
318 if (nbuf > 10) {
319 printf("%s\n",line);
320 snprintf(line,100,"%5d : ",i);
321 nbuf = 1;
322 }
323 snprintf(element,10,"%7lld ",fList[i]);
324 strlcat(line,element,100);
325 }
326 if (nbuf) printf("%s\n",line);
327 delete [] line;
328}
329
330////////////////////////////////////////////////////////////////////////////////
331/// Reset number of entries in event list.
332
334{
335 fN = 0;
336}
337
338////////////////////////////////////////////////////////////////////////////////
339/// Resize list by delta entries.
340
342{
343 if (!delta) delta = fDelta;
344 fSize += delta;
345 Long64_t *newlist = new Long64_t[fSize];
346 for (Int_t i=0;i<fN;i++) newlist[i] = fList[i];
347 delete [] fList;
348 fList = newlist;
349}
350
351////////////////////////////////////////////////////////////////////////////////
352/// Remove reference to this EventList from current directory and add
353/// reference to new directory dir. dir can be 0 in which case the list
354/// does not belong to any directory.
355
357{
358 if (fDirectory == dir) return;
359 if (fDirectory) fDirectory->Remove(this);
360 fDirectory = dir;
361 if (fDirectory) fDirectory->Append(this);
362}
363
364////////////////////////////////////////////////////////////////////////////////
365/// Change the name of this TEventList.
366
367void TEventList::SetName(const char *name)
368{
369 // TEventLists are named objects in a THashList.
370 // We must update the hashlist if we change the name
371 if (fDirectory) fDirectory->Remove(this);
372 fName = name;
373 if (fDirectory) fDirectory->Append(this);
374}
375
376////////////////////////////////////////////////////////////////////////////////
377/// Sort list entries in increasing order
378
380{
381 Int_t *index = new Int_t[fN];
382 Long64_t *newlist = new Long64_t[fSize];
383 Int_t i,ind;
384 TMath::Sort(fN,fList,index); //sort in decreasing order
385 for (i=0;i<fN;i++) {
386 ind = index[fN-i-1];
387 newlist[i] = fList[ind];
388 }
389 for (i=fN;i<fSize;i++) {
390 newlist[i] = 0;
391 }
392 delete [] index;
393 delete [] fList;
394 fList = newlist;
395}
396
397////////////////////////////////////////////////////////////////////////////////
398/// Stream an object of class TEventList.
399
401{
402 if (b.IsReading()) {
403 UInt_t R__s, R__c;
404 Version_t R__v = b.ReadVersion(&R__s, &R__c);
405 fDirectory = nullptr;
406 if (R__v > 1) {
407 b.ReadClassBuffer(TEventList::Class(), this, R__v, R__s, R__c);
409 return;
410 }
411 //====process old versions before automatic schema evolution
413 b >> fN;
414 b >> fSize;
415 b >> fDelta;
416 if (fN) {
417 Int_t *tlist = new Int_t[fSize];
418 b.ReadFastArray(tlist,fN);
419 fList = new Long64_t[fSize];
420 for (Int_t i=0;i<fN;i++) fList[i] = tlist[i];
421 delete [] tlist;
422 }
424 b.CheckByteCount(R__s, R__c, TEventList::IsA());
425 //====end of old versions
426
427 } else {
428 b.WriteClassBuffer(TEventList::Class(), this);
429 }
430}
431
432////////////////////////////////////////////////////////////////////////////////
433/// Remove elements from this list that are present in alist.
434
436{
437 if (!alist) return;
438 if (!fList) return;
439
440 Long64_t *newlist = new Long64_t[fN];
441 Int_t newpos = 0;
442 Int_t i;
443 for (i=0;i<fN;i++) {
444 if (alist->GetIndex(fList[i]) < 0) {
445 newlist[newpos] = fList[i];
446 newpos++;
447 }
448 }
449 delete [] fList;
450 fN = newpos;
451 fList = newlist;
452
453 TCut orig = GetTitle();
454 TCut removed = alist->GetTitle();
455 TCut updated = orig && !removed;
456 SetTitle(updated.GetTitle());
457}
458
459////////////////////////////////////////////////////////////////////////////////
460/// Assingment.
461
463{
464 if (this != &list) {
465 TNamed::operator=(list);
466 if (fSize < list.fSize) {
467 delete [] fList;
468 fList = new Long64_t[list.fSize];
469 }
470 fN = list.fN;
471 fSize = list.fSize;
472 fDelta = list.fDelta;
473 for (Int_t i=0; i<fN; i++)
474 fList[i] = list.fList[i];
475 }
476 return *this;
477}
478
479////////////////////////////////////////////////////////////////////////////////
480/// Addition.
481
482TEventList operator+(const TEventList &list1, const TEventList &list2)
483{
484 TEventList newlist = list1;
485 newlist.Add(&list2);
486 return newlist;
487}
488
489////////////////////////////////////////////////////////////////////////////////
490/// Subtraction
491
492TEventList operator-(const TEventList &list1, const TEventList &list2)
493{
494 TEventList newlist = list1;
495 newlist.Subtract(&list2);
496 return newlist;
497}
498
499////////////////////////////////////////////////////////////////////////////////
500/// Intersection.
501
502TEventList operator*(const TEventList &list1, const TEventList &list2)
503{
504 TEventList newlist = list1;
505 newlist.Intersect(&list2);
506 return newlist;
507}
508
#define b(i)
Definition RSha256.hxx:100
short Version_t
Definition RtypesCore.h:65
long long Long64_t
Definition RtypesCore.h:69
const char Option_t
Definition RtypesCore.h:66
#define ClassImp(name)
Definition Rtypes.h:382
#define gDirectory
Definition TDirectory.h:384
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:1540
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:164
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
void Streamer(TBuffer &) override
Stream an object of class TObject.
const char * GetTitle() const override
Returns title of object.
Definition TNamed.h:48
TString fName
Definition TNamed.h:32
TNamed & operator=(const TNamed &rhs)
TNamed assignment operator.
Definition TNamed.cxx:51
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition TObject.cxx:213
virtual Bool_t InheritsFrom(const char *classname) const
Returns kTRUE if object inherits from class "classname".
Definition TObject.cxx:530
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:993
void ResetBit(UInt_t f)
Definition TObject.h:198
@ kMustCleanup
if object destructor must call RecursiveRemove()
Definition TObject.h:64
TLine * line
Short_t Max(Short_t a, Short_t b)
Returns the largest of a and b.
Definition TMathBase.h:250
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:431
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:347