ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/MitAna/TAM/src/TAMSelector.cxx
Revision: 1.21
Committed: Wed Mar 28 12:15:37 2012 UTC (13 years, 1 month ago) by paus
Content type: text/plain
Branch: MAIN
CVS Tags: Mit_025e
Changes since 1.20: +12 -2 lines
Log Message:
Enable skimming.

File Contents

# User Rev Content
1 loizides 1.1 //
2 paus 1.21 // $Id: TAMSelector.cxx,v 1.20 2011/03/21 15:58:37 paus Exp $
3 loizides 1.1 //
4    
5 loizides 1.17 #include "MitAna/TAM/interface/TAMSelector.h"
6 loizides 1.1
7    
8     #ifndef ROOT_RVersion
9     #include <RVersion.h>
10     #endif
11     #ifndef ROOT_TFile
12     #include "TFile.h"
13     #endif
14     #ifndef ROOT_TTree
15     #include "TTree.h"
16     #endif
17     #ifndef ROOT_TError
18     #include "TError.h"
19     #endif
20     #ifndef ROOT_TROOT
21     #include "TROOT.h"
22     #endif
23     #ifndef ROOT_TList
24     #include "TList.h"
25     #endif
26     #ifndef ROOT_TClass
27     #include "TClass.h"
28     #endif
29     #ifndef ROOT_TProcessID
30     #include "TProcessID.h"
31     #endif
32     #if ROOT_VERSION_CODE >= ROOT_VERSION(5,0,0) && \
33     ROOT_VERSION_CODE < ROOT_VERSION(5,11,2)
34     #ifndef ROOT_TUrl
35     #include "TUrl.h"
36     #endif
37     #endif
38 loizides 1.6 #ifndef ROOT_TRefTable
39     #include <TRefTable.h>
40     #endif
41     #ifndef ROOT_TBranchRef
42     #include <TBranchRef.h>
43     #endif
44     #ifndef ROOT_TFriendElement
45     #include "TFriendElement.h"
46     #endif
47 loizides 1.1 #ifndef TAM_TAModule
48 loizides 1.17 #include "MitAna/TAM/interface/TAModule.h"
49 loizides 1.1 #endif
50     #ifndef TAM_TAMOutput
51 loizides 1.17 #include "MitAna/TAM/interface/TAMOutput.h"
52 loizides 1.1 #endif
53     #ifndef TAM_TAMVirtualLoader
54 loizides 1.17 #include "MitAna/TAM/interface/TAMVirtualLoader.h"
55 loizides 1.1 #endif
56     #ifndef TAM_TAMVirtualBranchLoader
57 loizides 1.17 #include "MitAna/TAM/interface/TAMVirtualBranchLoader.h"
58 loizides 1.1 #endif
59     #ifndef TAM_TAMTreeLoader
60 loizides 1.17 #include "MitAna/TAM/interface/TAMTreeLoader.h"
61 loizides 1.1 #endif
62 bendavid 1.10 #include <map>
63 loizides 1.1
64     //////////////////////////////////////////////////////////////////////////
65     // //
66     // TAMSelector //
67     // //
68     // A selector class for modular processing of a tree (or chain). //
69     // Processing is done by TAModule(s). //
70     // //
71     // Each TAModule contains pointers to the objects that will be read //
72     // in from the tree. For each such pointer, it will call //
73     // ReqBranch(name, pointer) to request that the branch be available for //
74     // processing by the module. The module's pointer would then be //
75     // automatically set to point to the object for the current event by //
76     // the TAMSelector, whenever the module calls LoadBranch. //
77     // //
78     // TAMSelector provides (and hides) all interaction with the tree //
79     // for the TAModule(s). //
80     // //
81     // Author : Corey Reed 07/20/2004 //
82     // Maarten Ballintijn 12/06/2005 //
83     // Constantin Loizides 12/06/2005 //
84     // //
85     //////////////////////////////////////////////////////////////////////////
86    
87     ClassImp(TAMSelector)
88    
89    
90     #if ROOT_VERSION_CODE < ROOT_VERSION(5,11,3)
91     #define R__ASSERT(e) \
92     if (!(e)) Fatal("", kAssertMsg, _QUOTE_(e), __LINE__, __FILE__)
93     #endif
94    
95 loizides 1.6 //______________________________________________________________________________
96 loizides 1.14 TAMSelector::BranchProxy::BranchProxy(TAMSelector *sel, Bool_t e) :
97     fSel(sel),
98     fOrig(0),
99 loizides 1.15 fFake(new TRefTable(this, 10)),
100 loizides 1.14 fCurEntry(-1)
101     {
102     // Default constructor.
103    
104     memset(fBrRead, 0, 1024*sizeof(Bool_t));
105 loizides 1.15
106 loizides 1.14 if (e)
107     Enable();
108     }
109    
110    
111     //______________________________________________________________________________
112     TAMSelector::BranchProxy::~BranchProxy()
113     {
114     // Destructor.
115    
116     Disable();
117     delete fFake;
118     fFake = 0;
119     }
120    
121    
122     //______________________________________________________________________________
123     void TAMSelector::BranchProxy::Disable()
124     {
125     // Disable the auto branch loading via TAM.
126    
127     if (fOrig)
128     TRefTable::SetRefTable(fOrig);
129     }
130    
131    
132     //______________________________________________________________________________
133     void TAMSelector::BranchProxy::Enable()
134     {
135     // Enable the auto branch loading via TAM.
136    
137     TBranchRef *bref = fSel->GetTree()->GetTree()->GetBranchRef();
138     if (bref) {
139     fOrig = bref->GetRefTable();
140     TRefTable::SetRefTable(fFake);
141     }
142     }
143    
144    
145     //______________________________________________________________________________
146     TObject *TAMSelector::BranchProxy::GetObjectWithID(UInt_t uid, TProcessID *pid)
147     {
148     // This function should be called by a reference class which is designed to
149     // use TAM. Branch loading via TAM is handled through the BranchProxy class.
150     // Optimized loading is used if the object table cleaning has been
151     // enabled in TAM.
152 paus 1.21
153     // printf("enter TAMSelector::BranchProxy::GetObjectWithID: %d\n",uid);
154 loizides 1.14
155     if (fSel->GetObjTabClean()) {
156     // can trust that object table was cleaned
157     TObject *obj = pid->GetObjectWithID(uid);
158     if (obj)
159     return obj;
160     }
161    
162     TBranchRef *br = fSel->GetTree()->GetTree()->GetBranchRef();
163     if (!br)
164     return pid->GetObjectWithID(uid);
165 paus 1.21
166     //printf(" br table loaded %p\n",(void*)br);
167    
168 loizides 1.14 TRefTable *table = br->GetRefTable();
169     if (!table)
170     return pid->GetObjectWithID(uid);
171 paus 1.21
172     //printf(" ref table loaded %p\n",(void*)table);
173    
174 loizides 1.14 table->SetUID(uid,pid);
175    
176 paus 1.21
177     //printf("loading\n");
178 loizides 1.14 Load(uid,pid,br,table);
179 paus 1.21 //printf(" loaded %p\n",(void*)pid->GetObjectWithID(uid));
180 loizides 1.14
181     return pid->GetObjectWithID(uid);
182     }
183    
184    
185     //______________________________________________________________________________
186     Bool_t TAMSelector::BranchProxy::Load(UInt_t uid, TProcessID *pid,
187     TBranchRef *br, TRefTable *table)
188     {
189     // Find and load branch corresponding to given uid/pid.
190    
191     // read entry for this event
192     Long64_t readentry = br->GetReadEntry();
193     Long64_t tamentry = fSel->fCurEvt;
194     if (readentry!=tamentry) {
195     Fatal("Load",
196 bendavid 1.19 "Entries from BranchRef (%lld) differs from TAM current entry (%lld)",
197 loizides 1.14 readentry, tamentry);
198     }
199    
200 bendavid 1.18 Long64_t chainentry = readentry + fSel->GetTree()->GetChainOffset();
201    
202 loizides 1.14 // read branchref if needed
203 bendavid 1.18 if (fCurEntry != chainentry) {
204 loizides 1.14 Int_t bytes = br->GetEntry(readentry);
205     if (bytes<0) {
206 bendavid 1.19 Fatal("Load", "Could not get entry %lld from %s branch",
207 loizides 1.14 readentry, br->GetName());
208     }
209     }
210    
211     // get branch from RefTable
212     TBranch *branch = dynamic_cast<TBranch*>(table->GetParent(uid, pid));
213    
214     if (!branch) { //scan for possible friend Trees
215     TList *friends = fSel->GetTree()->GetListOfFriends();
216     if (friends) {
217    
218     // reset branch read flags if new entry to be read
219 bendavid 1.18 if (fCurEntry != chainentry)
220 loizides 1.14 memset(fBrRead, 0, 1024*sizeof(Bool_t));
221    
222     TObjLink *lnk = friends->FirstLink();
223     Int_t nfriend = 0;
224     while (lnk) {
225     TFriendElement *elem =
226     static_cast<TFriendElement*>(lnk->GetObject());
227     TBranchRef *bref = elem->GetTree()->GetBranchRef();
228     if (bref) {
229     if (!fBrRead[nfriend]) {
230     bref->GetEntry(readentry);
231     fBrRead[nfriend] = kTRUE;
232     }
233     TBranch *branch = dynamic_cast<TBranch*>
234     (bref->GetRefTable()->GetParent(uid, pid));
235     if (branch)
236     break; // found branch
237     }
238     lnk = lnk->Next();
239     ++nfriend;
240     R__ASSERT(nfriend<1024);
241     }
242     }
243     }
244    
245     // cache last branch read attempt
246 bendavid 1.18 fCurEntry = chainentry;
247 loizides 1.14
248     if (!branch) {
249     return kFALSE;
250     }
251    
252     // access the main branch
253     TBranch *readbranch = branch->GetMother();
254     if (!readbranch) {
255     return kFALSE;
256     }
257    
258     // check if TAMBranchInfo already exists
259     // and if not add it
260     const char *brname = readbranch->GetName();
261     TObject *brInfo = fSel->fBranchTable.FindObject(brname);
262     TAMBranchInfo *branchInfo;
263     if (brInfo==0) {
264     branchInfo = new TAMBranchInfo(brname);
265     fSel->fBranchTable.Add(branchInfo);
266     }
267     else
268     branchInfo = dynamic_cast<TAMBranchInfo*>(brInfo);
269    
270     // load the branch
271     fSel->LoadBranch(branchInfo);
272    
273     return kTRUE;
274     }
275    
276    
277     //______________________________________________________________________________
278     Bool_t TAMSelector::BranchProxy::Notify()
279     {
280     // Notification received from TRefTable::Notify. Originally this would
281     // have been address to the TBranchRef owner of the TRefTable, but
282     // we intercept this call.
283    
284 loizides 1.15 if (!fOrig)
285     return kFALSE;
286    
287 loizides 1.14 TBranchRef *br = dynamic_cast<TBranchRef*>(fOrig->GetOwner());
288     if (!br)
289     return kFALSE;
290    
291     // get the desired uid/pid pair
292     UInt_t uid = fFake->GetUID();
293     TProcessID *pid = fFake->GetUIDContext();
294     fOrig->SetUID(uid,pid);
295     if (0) // this is what ROOT would have done
296     return br->Notify();
297    
298     return Load(uid, pid, br, fOrig);
299     }
300    
301     //______________________________________________________________________________
302 loizides 1.1 TAMSelector::TAMSelector() :
303     fTree(0),
304 paus 1.20 fTreeCache(0),
305     fCacheSize(0),
306 loizides 1.5 fBranchTable(TCollection::kInitHashTableCapacity, 1),
307     fEventObjs(TCollection::kInitHashTableCapacity, 1),
308 loizides 1.1 fAModules(new TAModule("TAMTopModule",
309     "Top-most module containing all other modules.")),
310     fCurEvt(-1),
311     fOwnInput(0),
312     fAnalysisAborted(kFALSE),
313     fModAborted(kFALSE),
314     fEventAborted(kFALSE),
315     fActNotify(kFALSE),
316     fObjCounter(0),
317     fVerbosity(0),
318 loizides 1.14 fProxy(this),
319 loizides 1.6 fDoProxy(kFALSE),
320 loizides 1.14 fDoObjTabClean(kFALSE),
321     fObjTabClean(kFALSE),
322     fLoaders()
323 loizides 1.1 {
324     // Default constructor.
325    
326     fBranchTable.SetOwner(kTRUE);
327     fEventObjs.SetOwner(kTRUE);
328     gROOT->GetListOfBrowsables()->Add(fAModules,"Analysis Modules");
329     fLoaders.SetName("TAM_LOADERS");
330     }
331    
332    
333     //______________________________________________________________________________
334     TAMSelector::~TAMSelector()
335     {
336     // Destructor.
337    
338     gROOT->GetListOfBrowsables()->Remove(fAModules);
339     TList *submods = fAModules->GetListOfTasks();
340     if (submods!=0) submods->Clear("nodelete");
341     delete fAModules;
342     delete fOwnInput;
343     }
344    
345    
346     //______________________________________________________________________________
347     void TAMSelector::AbortAnalysis()
348     {
349     // Aborts the analysis by aborting all modules and never executing
350     // further module funtions.
351    
352     // no way known to propagate the stop of analysis in PROOF from
353     // one slave to the other slaves (and master).
354     Fatal("AbortAnalysis", "Analysis aborted!");
355    
356     fAModules->DeactivateAll();
357     fAnalysisAborted = kTRUE;
358     }
359    
360    
361     //______________________________________________________________________________
362     void TAMSelector::AbortEvent()
363     {
364     // Aborts the current event by setting all modules inactive.
365    
366     fAModules->DeactivateAll();
367     fEventAborted = kTRUE;
368     }
369    
370    
371     //______________________________________________________________________________
372     void TAMSelector::AbortModule(TAModule *mod)
373     {
374     // Abort the specified module by setting it (and its sub modules) inactive.
375    
376     mod->DeactivateAll();
377     fModAborted = kTRUE;
378     }
379    
380    
381     //______________________________________________________________________________
382     void TAMSelector::AddInput(TAModule *mod)
383     {
384     // Adds the module to the top-most module in the hierarchy.
385     // The "misnomer" on the function name is intentional:
386     // This allows users to call gProof->AddInput(myMod) when
387     // using proof and mySelector->AddInput(myMod) when not
388     // using proof.
389    
390     fAModules->Add(mod);
391     }
392    
393     //______________________________________________________________________________
394     Bool_t TAMSelector::AddObjThisEvt(TObject *obj)
395     {
396     // Add this object to the list of objects stored for this event.
397     // See further description below.
398    
399 paus 1.4 if (obj)
400 loizides 1.1 return AddObjThisEvt(obj,obj->GetName());
401     else {
402 paus 1.4 Error("AddObjThisEvt","Can not add null object to event.");
403 loizides 1.1 return kFALSE;
404     }
405     }
406    
407    
408     //______________________________________________________________________________
409     Bool_t TAMSelector::AddObjThisEvt(TObject *obj, const char *name)
410     {
411     // Add this object to the list of objects stored for this event.
412     // NOTE:
413     // - The object must have a unique name.
414     // - The object must be on the heap.
415     // - The object will be owned by the selector and deleted at the
416     // end of the processing of the current event.
417     //
418     // Returns true iff the object meets the requirements and is added to
419     // the list successfully.
420    
421     if (obj!=0) {
422     if (FindObjThisEvt(name)==0) {
423     if (obj->IsOnHeap()) {
424     fEventObjs.Add(new TAMEvtObj(obj,name));
425 loizides 1.5 return kTRUE;
426 loizides 1.1 } else {
427     Error("AddObjThisEvt",
428     "Object [%s] does not appear to be on heap. "
429     "Can not add object.",
430     name);
431     }
432     } else {
433     Error("AddObjThisEvt",
434     "Object with name [%s] already added to this event.",
435     name);
436     }
437     } else {
438     Error("AddObjThisEvt",
439     "Can not add null object to event.");
440     }
441     return kFALSE;
442     }
443    
444    
445     //______________________________________________________________________________
446     void TAMSelector::AddNewOutputLists()
447     {
448     // Make a new hierarchy of TAMOutput objects and add it to
449     // the output of this selector. Note that this should only be called
450     // ONCE (by SlaveBegin())!!
451    
452     fAModules->NewOutputList(GetOutputList());
453     }
454    
455    
456     //______________________________________________________________________________
457     void TAMSelector::Begin(TTree */*tree*/)
458     {
459     // The Begin() function is called at the start of the query.
460     // When running with PROOF Begin() is only called on the client.
461     // The tree argument is deprecated (on PROOF 0 is passed).
462     // It checks the tree of modules to be sure each module is set
463     // up to use this selector.
464     // Adds the tree of modules to the input list of this selector.
465     // Builds the hierarchy of TAMOutput objects in the output list.
466     // Then calls Begin() on the TAModule(s). At the end, it adds the
467     // list of loaders to the input list.
468    
469     if (fAnalysisAborted) {
470     return;
471     } else {
472     if (fInput==0) {
473     fOwnInput = new TList;
474     SetInputList(fOwnInput);
475     }
476     CopyModsFromInput();
477     fAModules->SetSelector(this);
478     if (fAModules->CheckSelectors(this)) {
479     if (fModAborted) {
480     fAModules->ResetAllActiveFlags();
481     fModAborted = kFALSE;
482     }
483     fAModules->ExecuteTask(&TAModule::kExecBegin);
484     if (fEventObjs.IsEmpty()==kFALSE) fEventObjs.Delete();
485     } else {
486     Error("Begin",
487     "One or more module is not set up to use this selector.");
488     R__ASSERT(0);
489     }
490     }
491    
492     fInput->AddLast(&fLoaders);
493     }
494    
495 loizides 1.14
496 bendavid 1.10 //______________________________________________________________________________
497     void TAMSelector::CleanObjTable(TProcessID *pid, UInt_t lastKeptUID) const
498     {
499 loizides 1.14 // Clear from the object table of a given process id, all pointers for objects
500     // with uids after lastKeptUid.
501 bendavid 1.11
502 bendavid 1.10 TObjArray *objTable = pid->GetObjects();
503    
504     Int_t lastIdxKept = lastKeptUID & 0xffffff;
505     Int_t last = objTable->GetLast();
506    
507     if (last==-1)
508     return;
509    
510     if (lastIdxKept>=0 && lastIdxKept<=last) {
511     TObject **cont = objTable->GetObjectRef(0);
512     memset(&cont[lastIdxKept+1],0,(last-lastIdxKept)*sizeof(TObject*));
513     objTable->SetLast(lastIdxKept);
514     }
515     else
516     Error("TAMSelector::CleanObjTable",
517     "Out of Bounds trying to clean object table from Process ID.");
518     }
519 loizides 1.1
520 loizides 1.14
521 loizides 1.1 //______________________________________________________________________________
522     void TAMSelector::ClearAllLoaders()
523     {
524     // Go through the list of requested branches and clear the
525     // loaders that were used.
526    
527     TIter nextBr(fBranchTable.MakeIterator());
528 paus 1.4 while (TAMBranchInfo *br = dynamic_cast<TAMBranchInfo*>(nextBr())) {
529 loizides 1.1
530     if (!br->fIsLoaded) continue;
531     // don't bother with branches that weren't loaded
532    
533     TAMVirtualBranchLoader *l = br->GetLoader();
534     if (l) {
535     l->Clear();
536     } else {
537     Error("ClearAllLoaders",
538     "Could not get loader for [%s]. Unable to "
539     "try to clear loader for this branch.",
540     br->GetName());
541     }
542     }
543     }
544    
545    
546     //______________________________________________________________________________
547     void TAMSelector::CopyModsFromInput()
548     {
549     // Find all TAModules in the input list and copy them to fAModules.
550    
551     R__ASSERT(fInput!=0);
552    
553     if (fInput->IsEmpty()) return;
554    
555     TObject *obj=fInput->First(), *tobj=0;
556     while (obj!=0) {
557     tobj = fInput->After(obj);
558     if (obj->InheritsFrom(TAModule::Class())) {
559     AddInput(dynamic_cast<TAModule*>(obj));
560     }
561     obj = tobj;
562     }
563     }
564    
565    
566     //______________________________________________________________________________
567     Bool_t TAMSelector::FindLoader(TAMBranchInfo *brInfo)
568     {
569     // Find the loader responsible for the given branch info. Return kFALSE
570     // if none could be found in the list of loaders.
571    
572     const char *bname = brInfo->GetName();
573     TIter next(&fLoaders);
574    
575     while( TAMVirtualLoader *l = dynamic_cast<TAMVirtualLoader*>(next()) ) {
576     TAMVirtualBranchLoader *bl = l->CreateBranchLoader(fTree, brInfo);
577     if (bl != 0) {
578     Info("FindLoader","branch '%s' use loader %s",
579     bname, bl->IsA()->GetName());
580     brInfo->SetLoader(bl);
581     return kTRUE;
582     }
583     }
584    
585     Error("FindLoader","No loader found for branch '%s'", bname);
586     return kFALSE;
587     }
588    
589    
590     //______________________________________________________________________________
591     TAMOutput *TAMSelector::FindModOutput(const TAModule *mod)
592     {
593     // Find the TAMOutput object corresponding to the given module.
594    
595     TAMOutput *obj=0, *tobj=0;
596     TIter nextObj(GetOutputList());
597     while ( (obj = dynamic_cast<TAMOutput*>(nextObj())) ) {
598     tobj = obj->FindModOutput(mod); // search mod and sub mods
599     if (tobj!=0) return tobj;
600     }
601     return 0;
602     }
603    
604    
605     //______________________________________________________________________________
606     TObject *TAMSelector::FindObjThisEvt(const Char_t *name) const
607     {
608     // Looks for the object with the specified name that was added to
609     // this event. If not found, returns zero.
610    
611     TAMEvtObj *eo = dynamic_cast<TAMEvtObj*>(fEventObjs.FindObject(name));
612     if (eo!=0) {
613     return eo->fObj;
614     }
615     return 0;
616     }
617    
618    
619     //______________________________________________________________________________
620     TObject *TAMSelector::FindPublicObj(const Char_t *name) const
621     {
622     // Looks for the public object with the specified name. If not found,
623     // returns zero.
624     // Note: TAModules are not public objects and will not be found by this
625     // function.
626    
627     if (fInput!=0) {
628     TIter nextObj(fInput);
629     TObject *obj=0;
630     TString nm(name);
631     while ( (obj = nextObj()) ) {
632     if ( (nm.CompareTo(obj->GetName())==0) &&
633     (obj->InheritsFrom(TAModule::Class())==kFALSE) ) {
634     return obj;
635     }
636     }
637     } else {
638     Error("FindPublicObj",
639     "Input list is null. No public objects exist.");
640     }
641     return 0;
642     }
643    
644    
645     //______________________________________________________________________________
646 loizides 1.5 const TFile *TAMSelector::GetCurrentFile() const
647     {
648     // Returns the current file that the tree is in.
649    
650     return (fTree) ? (const_cast<const TFile *>(fTree->GetCurrentFile())) : 0;
651     }
652    
653    
654     //______________________________________________________________________________
655     TFile *TAMSelector::GetCurrentFile()
656 loizides 1.1 {
657     // Returns the current file that the tree is in.
658    
659     return (fTree) ? (fTree->GetCurrentFile()) : 0;
660     }
661    
662    
663     //______________________________________________________________________________
664     const TAMOutput *TAMSelector::GetModOutput() const
665     {
666     // Return the top-most module output object.
667    
668     return fAModules->GetModOutput();
669     }
670    
671    
672     //______________________________________________________________________________
673     TAMOutput *TAMSelector::GetModOutput()
674     {
675     // Return the top-most module output object.
676    
677     return fAModules->GetModOutput();
678     }
679    
680 loizides 1.14
681 bendavid 1.12 //______________________________________________________________________________
682     TObject *TAMSelector::GetObjectWithID(UInt_t uid, TProcessID *pid)
683     {
684 loizides 1.14 // This function should be called by a reference class which is designed to
685     // use TAM. Branch loading via TAM is handled through the BranchProxy class.
686     // Optimized loading is used if the object table cleaning has been
687     // enabled in TAM.
688 bendavid 1.12
689 loizides 1.14 return fProxy.GetObjectWithID(uid,pid);
690     }
691 bendavid 1.12
692 loizides 1.1
693 paus 1.20 //__________________________________________________________________________________________________
694     void TAMSelector::Init(TTree *tree)
695     {
696     // Set the tree for this selector. The Init() function is called when the selector needs to
697     // initialize a new tree (or chain). It will be called many times when running with PROOF. The
698     // list of requested branches must have been constructed already.
699 loizides 1.1
700     if (tree==0) {
701     Error("Init", "Specified tree is null.");
702     AbortAnalysis();
703     }
704    
705     fTree = tree;
706    
707     TIter nextBranch(fBranchTable.MakeIterator());
708     while ( TAMBranchInfo *brInfo =
709 paus 1.20 dynamic_cast<TAMBranchInfo*>(nextBranch()) ) {
710     brInfo->Init();
711 loizides 1.1 }
712     }
713    
714 loizides 1.14
715 loizides 1.1 //______________________________________________________________________________
716 bendavid 1.8 void TAMSelector::LoadBranch(const Char_t *bname)
717 loizides 1.1 {
718 bendavid 1.8 // Loads branch by name, getting the pointer to the corresponding
719 loizides 1.9 // TAMBranchInfo and then use it in the call of the protected LoadBranch
720     // function.
721 bendavid 1.8
722 loizides 1.1 if(fCurEvt==-1) {
723     Error("LoadBranch",
724 loizides 1.2 "Can not load branch with name [%s] at this point (fCurEvt==-1).",
725 loizides 1.1 bname);
726     AbortAnalysis();
727     }
728    
729     TAMBranchInfo *brInfo =
730     dynamic_cast<TAMBranchInfo*>( fBranchTable.FindObject(bname) );
731    
732     if (brInfo==0) {
733     Error("LoadBranch",
734     "Could not find branch with name [%s] in requested branch list.",
735     bname);
736     AbortAnalysis();
737 bendavid 1.8 }
738    
739     LoadBranch(brInfo);
740     }
741    
742 loizides 1.14
743 bendavid 1.8 //______________________________________________________________________________
744     void TAMSelector::LoadBranch(TAMBranchInfo* brInfo)
745     {
746     // Loads the selected branch and get the current entry (number fCurEvt)
747     // if the branch has not yet been loaded for this event. It then makes
748     // sure each module's pointer to the branch object point to the address
749     // of this branch.
750 loizides 1.1
751 loizides 1.6 if (brInfo->IsLoaded())
752 loizides 1.1 return;
753    
754     // find loader for branch info if needed and notify
755     if (brInfo->GetLoader() == 0) {
756    
757     if ( !FindLoader(brInfo) ) {
758     Error("LoadBranch","Branch [%s] FindLoader() failed in file [%s].",
759     brInfo->GetName(),
760     (GetCurrentFile()!=0) ? (GetCurrentFile()->GetName()) : "null");
761    
762     AbortAnalysis();
763     }
764    
765     if ( !brInfo->Notify(fTree) )
766     Error("LoadBranch","Branch [%s] Notify() failed in file [%s].",
767     brInfo->GetName(),
768     (GetCurrentFile()!=0) ? (GetCurrentFile()->GetName()) : "null");
769    
770     // have to reset our object counter to take care of objects
771     // created in the notify of the loader
772 loizides 1.6 //fObjCounter=TProcessID::GetObjectCount();
773     }
774 loizides 1.1
775 paus 1.20 // load the entry (using the cache)
776     if (fTreeCache)
777     GetCurrentFile()->SetCacheRead(fTreeCache);
778 loizides 1.1 Int_t ret = brInfo->GetEntry(fCurEvt);
779 paus 1.20 if (fTreeCache)
780     GetCurrentFile()->SetCacheRead(0);
781    
782    
783 loizides 1.1 if(ret<0) {
784     Error("LoadBranch",
785     "Error in file [%s] when accessing branch with name [%s] in "
786     "requested branch list.",
787     (GetCurrentFile()!=0) ? (GetCurrentFile()->GetName())
788 bendavid 1.8 : "null", brInfo->GetName());
789 loizides 1.1
790     AbortAnalysis();
791     }
792 loizides 1.6
793     // make sure autoload proxy is enabled
794 loizides 1.7 if (fDoProxy)
795 loizides 1.14 fProxy.Enable();
796 loizides 1.1 }
797    
798    
799     //______________________________________________________________________________
800     Bool_t TAMSelector::Notify()
801     {
802     // Sets the branch pointers, checks that the types of the module's
803     // pointers correspond to what is in the tree and calls SetBranchAddress
804     // if necessary. Finally, it also notifies the modules.
805     // The Notify() function is called when a new file is opened. This
806     // can be either for a new TTree in a TChain or when a new TTree
807     // is started when using PROOF.
808    
809 loizides 1.6 // we are just in Notify(), therefore we ignore this call
810 loizides 1.1 if(fActNotify) return kTRUE;
811     fActNotify = kTRUE; //"lock" to protect for recursive calls
812    
813     #if ROOT_VERSION_CODE >= ROOT_VERSION(5,0,0) && \
814     ROOT_VERSION_CODE < ROOT_VERSION(5,11,2)
815    
816 loizides 1.6 // workaround for older ROOT: give file name to
817 loizides 1.1 // TFile that accidentally stripped the protocol
818 loizides 1.6 const TUrl *url = (GetCurrentFile()!=0)
819     ? GetCurrentFile()->GetEndpointUrl() : 0;
820 loizides 1.1
821 loizides 1.6 if (url==0) {
822     Warning("Notify","Could not get URL for file [%s].",
823     (GetCurrentFile()!=0) ? (GetCurrentFile()->GetName())
824     : "null");
825     } else {
826 loizides 1.1
827 loizides 1.6 if(GetCurrentFile()!=0)
828     GetCurrentFile()->SetName((const_cast<TUrl*>(url)->GetUrl()));
829     // the const_cast is for some version of
830     // 5.08.00 that didn't have GetUrl const.
831     }
832 loizides 1.1 #endif
833    
834     if (fVerbosity>0) {
835 bendavid 1.19 Info("Notify","Opening file [%s] at current event [%lld].",
836 loizides 1.1 (GetCurrentFile()!=0) ? (GetCurrentFile()->GetName()) : "null",
837     fCurEvt);
838     }
839    
840 loizides 1.6 // status (return value) of notify
841 loizides 1.1 Bool_t notifyStat = kTRUE;
842    
843 paus 1.20 // Set up the caching process
844     if (fCacheSize > 0) {
845     fTreeCache->SetLearnEntries(1);
846     fTree->SetCacheSize(128*1024*1024);
847     printf(" CurrentFile Name: %s\n",fTree->GetCurrentFile()->GetName());
848     if (fTreeCache)
849     delete fTreeCache;
850     fTreeCache = dynamic_cast<TTreeCache*>(fTree->GetCurrentFile()->GetCacheRead());
851     fTree->GetCurrentFile()->SetCacheRead(0);
852    
853     //// Read all data products (even if we don't use them).
854     //// Remove the below lines to read minimal sets of products.
855     //// Removing would cause many more I/O operations and slow the cluster.
856     //fTreeCache->StartLearningPhase();
857     //fTreeCache->AddBranch("*", kTRUE);
858     //fTreeCache->StopLearningPhase();
859     }
860    
861 loizides 1.1 // no event yet processed eg, no loaders assigned,
862     // so that the notify is being delayed to LoadBranch()
863     if(fCurEvt>=0) {
864    
865     TIter nextBranch(fBranchTable.MakeIterator());
866    
867     while ( TAMBranchInfo *brInfo =
868     dynamic_cast<TAMBranchInfo*>(nextBranch()) ) {
869    
870     // notify the branchinfo and its corresponding loader
871     notifyStat &= brInfo->Notify(fTree);
872    
873     // have to reset our object counter to take care of objects
874     // created in the notify of the loader
875 loizides 1.6 //fObjCounter=TProcessID::GetObjectCount();;
876 loizides 1.1
877     if (notifyStat==kFALSE) {
878     Error("Notify","Branch [%s] Notify() failed in file [%s].",
879     brInfo->GetName(),
880     (GetCurrentFile()!=0) ? (GetCurrentFile()->GetName())
881     : "null");
882    
883     AbortAnalysis();
884     }
885     }
886 paus 1.20
887     // status (return value) of notify
888     notifyStat = kTRUE;
889    
890 loizides 1.1 }
891    
892     if (notifyStat && (fAnalysisAborted==kFALSE)) {
893     notifyStat &= fAModules->NotifyAll();
894     }
895    
896     //release "lock" on notify
897     fActNotify = kFALSE;
898    
899     return notifyStat;
900     }
901    
902    
903     //______________________________________________________________________________
904     Bool_t TAMSelector::Process(Long64_t entry)
905     {
906     // Zero's all TAModule(s)'s branch addresses and calls Process()
907     // on the TAModule(s). Then clear any requested clones arrays and
908     // clear the list of event-stored objects.
909     //
910     // The Process() function is called for each entry in the tree (or possibly
911     // keyed object in the case of PROOF) to be processed.
912    
913     if (fAnalysisAborted) {
914     // There is currently no reliable way to exit cleanly back to the
915     // interpreter, so simply return immediately.
916     return kFALSE;
917 loizides 1.2 }
918 loizides 1.1
919 loizides 1.6 // prepare for event
920 loizides 1.2 fCurEvt = entry;
921     ZeroAllBranches();
922     if (fModAborted) {
923     fAModules->ResetAllActiveFlags();
924     fModAborted = kFALSE;
925     } else if (fEventAborted) {
926     fAModules->ResetAllActiveFlags();
927     fEventAborted = kFALSE;
928     }
929 loizides 1.1
930 loizides 1.6 // store object counter for run boundaries
931 loizides 1.2 if (BeginRun()) {
932     fObjCounterRun=TProcessID::GetObjectCount();;
933     fAModules->ExecuteTask(&TAModule::kExecBeginRun);
934     }
935 loizides 1.1
936 loizides 1.6 // store object counter and process event
937 loizides 1.2 fObjCounter=TProcessID::GetObjectCount();;
938 bendavid 1.10 TObjArray *pids = GetCurrentFile()->GetListOfProcessIDs();
939 loizides 1.14
940     // store uids for later cleaning for object tables from file
941 bendavid 1.10 std::map<TProcessID*, UInt_t> lastUIDs;
942 loizides 1.14 if (fDoObjTabClean) {
943     for (Int_t i=0; i<pids->GetEntriesFast(); ++i) {
944     TProcessID *pid = static_cast<TProcessID*>(pids->At(i));
945     Int_t last = pid->GetObjects()->GetLast();
946     if (last==-1)
947     last = 0;
948     lastUIDs[pid] = last;
949     }
950     // set clean object table state variable to true
951     fObjTabClean = kTRUE;
952 bendavid 1.10 }
953    
954 loizides 1.2 if (fVerbosity>9) {
955     if ((entry % 100)==0) {
956     fprintf(stderr,"Processing entry %lld... \r",
957     entry);
958     }
959     }
960 loizides 1.6
961 loizides 1.2 fAModules->ExecuteTask(&TAModule::kExecProcess);
962    
963 loizides 1.6 // update object counter on run boundaries
964 loizides 1.2 if (EndRun()) {
965     fAModules->ExecuteTask(&TAModule::kExecEndRun);
966     fObjCounter=fObjCounterRun;
967     }
968    
969 loizides 1.6 // cleanup
970 loizides 1.2 ClearAllLoaders();
971     if (fEventObjs.IsEmpty()==kFALSE)
972     fEventObjs.Delete();
973 loizides 1.1
974 loizides 1.14 // set clean object table state variable to false
975     fObjTabClean = kFALSE;
976    
977 loizides 1.2 // restore object counter
978     TProcessID::SetObjectCount(fObjCounter);
979 loizides 1.16
980     // Clean object table for current process id:
981     // This guarantees that objects which are not yet loaded in the
982     // current event have null pointers in the object table.
983 loizides 1.14 if (fDoObjTabClean) {
984 bendavid 1.12 CleanObjTable(TProcessID::GetSessionProcessID(), fObjCounter);
985    
986     //clean object tables for process ids being read from the file
987     for (Int_t i=0; i<pids->GetEntriesFast(); ++i) {
988     TProcessID *pid = static_cast<TProcessID*>(pids->At(i));
989     std::map<TProcessID*, UInt_t>::const_iterator lastUID = lastUIDs.find(pid);
990     if (lastUID != lastUIDs.end()) {
991     CleanObjTable(pid, lastUID->second);
992     }
993     else {
994     CleanObjTable(pid, 0);
995     }
996     }
997 bendavid 1.10 }
998 loizides 1.1
999 loizides 1.16 // Clean object table for current process id.
1000     // This guarantees that objects which are not yet loaded in the
1001     // current event have null pointers in the object table.
1002 loizides 1.14 if (fDoObjTabClean) {
1003     CleanObjTable(TProcessID::GetSessionProcessID(), fObjCounter);
1004    
1005     //clean object tables for process ids being read from the file
1006     for (Int_t i=0; i<pids->GetEntriesFast(); ++i) {
1007     TProcessID *pid = static_cast<TProcessID*>(pids->At(i));
1008 loizides 1.16 std::map<TProcessID*, UInt_t>::const_iterator lastUID =
1009     lastUIDs.find(pid);
1010 loizides 1.14 if (lastUID != lastUIDs.end())
1011     CleanObjTable(pid, lastUID->second);
1012     else
1013     CleanObjTable(pid, 0);
1014     }
1015     }
1016    
1017 loizides 1.1 return kTRUE;
1018     }
1019    
1020    
1021     //______________________________________________________________________________
1022     Bool_t TAMSelector::PublishObj(TObject *obj)
1023     {
1024     // Adds an object to a list of objects which is outside the module
1025     // hierarchy. This can be used to pass objects (for example, calibrations)
1026     // between modules. Objects in this list are available before Begin
1027     // until the end of SlaveTerminate. They are not guaranteed to be available
1028     // during or after Terminate.
1029     // Checks (by name) if the object is already in the list. If it is, returns
1030     // kFALSE and does not publish the object.
1031     // NOTE: These objects are NOT owned by the list! Whatever creates these
1032     // objects must take care to (1) remove the object from the list using
1033     // RetractObj() and (2) delete the object.
1034     // Also NOTE: This will not publish TAModule objects.
1035    
1036     if (obj!=0) {
1037     TObject *ob = FindPublicObj(obj->GetName());
1038     if (ob==0) {
1039     if (obj->InheritsFrom(TAModule::Class())) {
1040     // Disallow publishing of TAModules since it breaks the assumption
1041     // that all TAModules in the input list were added by AddInput and
1042     // are intended to be executed. (Such modules would be found by
1043     // TakeModsFromInput and CopyModsFromInput.)
1044     Warning("PublishObj",
1045     "Can not publish a TAModule. Object named [%s] not "
1046     "published.",
1047     obj->GetName());
1048     } else {
1049     if (fInput!=0) {
1050     fInput->Add(obj);
1051     return kTRUE;
1052     } else {
1053     Error("PublishObj",
1054     "Input list is null. Could not publish object named [%s].",
1055     obj->GetName());
1056     }
1057     }
1058     } else {
1059     Warning("PublishObj",
1060     "Public object named [%s] already exists.",
1061     obj->GetName());
1062     }
1063     } else {
1064     Error("PublishObj",
1065     "Can not publish null object.");
1066     }
1067     return kFALSE;
1068     }
1069    
1070    
1071     //______________________________________________________________________________
1072     TObject *TAMSelector::RemoveObjThisEvt(const Char_t *name)
1073     {
1074     // Finds the object with the specified name and removes it from
1075     // the list of objects added to this event.
1076     // Returns the object that was removed.
1077    
1078     TAMEvtObj *eo = dynamic_cast<TAMEvtObj*>(fEventObjs.FindObject(name));
1079     if (eo!=0) {
1080     fEventObjs.Remove(eo);
1081     TObject *obj = eo->fObj;
1082     eo->fObj = 0; // so it won't delete the object
1083     delete eo;
1084     return obj;
1085     }
1086     return 0;
1087     }
1088    
1089    
1090     //______________________________________________________________________________
1091     TObject *TAMSelector::RetractObj(const Char_t *name)
1092     {
1093     // Finds the public object with the specified name and removes it from
1094     // the list of public objects.
1095     // Returns the object that was retracted.
1096     // Note: TAmodules are not public objects and will not be removed by this
1097     // function.
1098    
1099     TObject *obj = FindPublicObj(name);
1100     if (obj!=0) {
1101     return fInput->Remove(obj);
1102     }
1103     return 0;
1104     }
1105    
1106    
1107     //______________________________________________________________________________
1108     void TAMSelector::SlaveBegin(TTree *tree)
1109     {
1110     // Inits the tree and calls SlaveBegin() on the TAModule(s)
1111     //
1112     // The SlaveBegin() function is called after the Begin() function.
1113     // When running with PROOF SlaveBegin() is called on each slave server.
1114     // The tree argument is deprecated (on PROOF 0 is passed).
1115    
1116     if (fAnalysisAborted) {
1117     return;
1118     } else {
1119     // call SlaveBegin first so the modules can call ReqBranch and
1120     // build the fBranchTable
1121     if (fModAborted) {
1122     fAModules->ResetAllActiveFlags();
1123     fModAborted = kFALSE;
1124     }
1125     // remove the modules from the input list and put them in the top-most
1126     // module of this selector. make sure they use this selector
1127     TakeLoadersFromInput();
1128     TakeModsFromInput();
1129     fAModules->SetSelector(this);
1130     if (fAModules->CheckSelectors(this)) {
1131     // build the module output hierarchy
1132     AddNewOutputLists();
1133     // set up the module's output members
1134     if (fAModules->fOutput != 0) {
1135     fAModules->fOutput->SetAllOutputMembers(kFALSE);
1136     }
1137    
1138     fAModules->ExecuteTask(&TAModule::kExecSlaveBegin);
1139     } else {
1140     Error("SlaveBegin",
1141     "One or more module is not set up to use this selector.");
1142     R__ASSERT(0);
1143     }
1144    
1145     // init requires fBranchTable to be built already
1146     // in Proof, this is called automatically when a new tree is loaded
1147 loizides 1.6 if (tree)
1148     Init(tree);
1149     if (fEventObjs.IsEmpty()==kFALSE)
1150     fEventObjs.Delete();
1151 loizides 1.1 }
1152     }
1153    
1154    
1155     //______________________________________________________________________________
1156     void TAMSelector::SlaveTerminate()
1157     {
1158     // Calls SlaveTerminate() on the TAModule(s)
1159     // The SlaveTerminate() function is called after all entries or objects
1160     // have been processed. When running with PROOF SlaveTerminate() is called
1161     // on each slave server.
1162    
1163     if (fAnalysisAborted) {
1164     return;
1165     } else {
1166     if (fModAborted) {
1167     fAModules->ResetAllActiveFlags();
1168     fModAborted = kFALSE;
1169     }
1170     fAModules->ExecuteTask(&TAModule::kExecSlaveTerminate);
1171     if (fEventObjs.IsEmpty()==kFALSE) fEventObjs.Delete();
1172     }
1173     }
1174    
1175    
1176     //______________________________________________________________________________
1177     void TAMSelector::TakeModsFromInput()
1178     {
1179     // Find all TAModules in the input list and move them to fAModules.
1180    
1181     R__ASSERT(fInput!=0);
1182    
1183     if (fInput->IsEmpty()==kFALSE) {
1184     TObject *obj=fInput->First(), *tobj=0;
1185     while (obj!=0) {
1186     tobj = fInput->After(obj);
1187     if (obj->InheritsFrom(TAModule::Class())) {
1188     AddInput(dynamic_cast<TAModule*>(obj));
1189     fInput->Remove(obj);
1190     }
1191     obj = tobj;
1192     }
1193     }
1194     }
1195    
1196    
1197     //______________________________________________________________________________
1198     void TAMSelector::TakeLoadersFromInput()
1199     {
1200     // Find all TAModules in the input list and copy them to fAModules.
1201    
1202     R__ASSERT(fInput!=0);
1203    
1204     TList *loaders = dynamic_cast<TList*>(fInput->FindObject("TAM_LOADERS"));
1205     if (loaders != 0) {
1206     TIter next(loaders);
1207     while ( TAMVirtualLoader *l =
1208     dynamic_cast<TAMVirtualLoader*>(next()) ) {
1209     if (loaders != &fLoaders) fLoaders.AddLast(l);
1210     }
1211     }
1212    
1213     fLoaders.AddLast(new TAMTreeLoader());
1214     }
1215    
1216    
1217     //______________________________________________________________________________
1218     void TAMSelector::Terminate()
1219     {
1220     // Calls Terminate() on the TAModule(s).
1221     // When running under Proof, will copy the TAMOutput object
1222     // from the fOutput list to the top most TAModule object. Assumes
1223     // that the only TAMOutput object in fOutput is the one belonging
1224     // to the top most module.
1225     // The Terminate() function is the last function to be called during
1226     // a query. It always runs on the client, it can be used to present
1227     // the results graphically or save the results to file.
1228    
1229     if (fAnalysisAborted) {
1230     return;
1231     }
1232    
1233     if (fModAborted) {
1234     fAModules->ResetAllActiveFlags();
1235     fModAborted = kFALSE;
1236     }
1237    
1238     if (fAModules->GetModOutput()==0) {
1239     // When running under Proof, copy the module output hierarchy to
1240     // this selector's top most module.
1241     TIter nextObj(GetOutputList());
1242     TAMOutput *tout=0;
1243     TObject *obj=0;
1244     Bool_t alreadySet=kFALSE;
1245     while ( (obj = nextObj()) ) {
1246     if (obj->InheritsFrom(TAMOutput::Class())) {
1247     tout = dynamic_cast<TAMOutput*>(obj);
1248     if (alreadySet) {
1249     Warning("Terminate",
1250     "Found more than one TAMOutput object in the "
1251     "output list. Assuming the first one is from the "
1252     "top-most module. The output list contains:");
1253     GetOutputList()->ls();
1254     } else {
1255     // copy module output hierarchy
1256     fAModules->SetAllModOutput(tout);
1257     // try to set module's pointers to their output objects
1258     tout->SetAllOutputMembers(kTRUE);
1259     alreadySet=kTRUE;
1260     }
1261     }
1262     }
1263     if (alreadySet==kFALSE) {
1264     Error("Terminate",
1265     "Did not find TAMOutput object in the output list! "
1266     "The output list contains:");
1267     GetOutputList()->ls();
1268     }
1269     }
1270    
1271     // move output objs from current objects to stored objects
1272     // incase the module wants to write output to a file
1273     if (fAModules->GetModOutput()!=0) {
1274     fAModules->GetModOutput()->StoreAllOutputObjs();
1275     } else {
1276     Error("Terminate",
1277     "Could not store output objects from this process run.");
1278     }
1279    
1280     fAModules->ExecuteTask(&TAModule::kExecTerminate);
1281     if (fEventObjs.IsEmpty()==kFALSE) fEventObjs.Delete();
1282    
1283     // store output objects again, in case a module added any output
1284     // objects during terminate
1285     if (fAModules->GetModOutput()!=0) {
1286     fAModules->GetModOutput()->StoreAllOutputObjs();
1287     } else {
1288     Error("Terminate",
1289     "Could not store output objects after terminate.");
1290     }
1291 paus 1.20
1292     // Delete the tree cache
1293     if (fTreeCache)
1294     delete fTreeCache;
1295 loizides 1.1 }
1296    
1297    
1298     //______________________________________________________________________________
1299     void TAMSelector::ZeroAllBranches()
1300     {
1301     // Loops through all branches in fBranchTable and sets all user addresses
1302     // for each branch to zero.
1303    
1304     TIter nextBranch(fBranchTable.MakeIterator());
1305 paus 1.4 while (TAMBranchInfo *brInfo =
1306     dynamic_cast<TAMBranchInfo*>(nextBranch())) {
1307 loizides 1.1 brInfo->ZeroUsrAddrs();
1308     }
1309     }