ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/MitAna/TAM/src/TAMSelector.cxx
Revision: 1.14
Committed: Mon Mar 2 12:34:51 2009 UTC (16 years, 2 months ago) by loizides
Content type: text/plain
Branch: MAIN
CVS Tags: Mit_009, Mit_008, Mit_008pre2, Mit_008pre1
Changes since 1.13: +244 -144 lines
Log Message:
BranchProxy that combines previous AutoLoad class with recent changes in TAMSelector.

File Contents

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