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