ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/MitAna/TAM/src/TAMSelector.cxx
Revision: 1.4
Committed: Mon Jul 7 16:41:53 2008 UTC (16 years, 10 months ago) by paus
Content type: text/plain
Branch: MAIN
CVS Tags: MITHEP_2_0_x
Changes since 1.3: +6 -7 lines
Log Message:
Implement simple Ascii file driven catalog. Easily extensible, though.

File Contents

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