ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/MitAna/TAM/src/TAMSelector.cxx
Revision: 1.6
Committed: Mon Oct 6 16:50:07 2008 UTC (16 years, 7 months ago) by loizides
Content type: text/plain
Branch: MAIN
Changes since 1.5: +186 -28 lines
Log Message:
Introduce auto branch loading (of TRefs) via TAM (via LoadBranch)

File Contents

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