ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/MitAna/TAM/src/TAMSelector.cxx
Revision: 1.2
Committed: Mon Jun 23 19:39:14 2008 UTC (16 years, 10 months ago) by loizides
Content type: text/plain
Branch: MAIN
Changes since 1.1: +38 -27 lines
Log Message:
Added Begin and EndRun to core TAM framework.

File Contents

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