ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/MitAna/TAM/src/TAMSelector.cxx
Revision: 1.18
Committed: Sun Sep 19 18:27:29 2010 UTC (14 years, 7 months ago) by bendavid
Content type: text/plain
Branch: MAIN
CVS Tags: Mit_018, Mit_017, Mit_017pre3, Mit_017pre2, Mit_017pre1, Mit_016, Mit_015b, Mit_015a, Mit_015, Mit_014e
Changes since 1.17: +6 -4 lines
Log Message:
Fix bug in auto-branch loading through TAM when reading multiple files

File Contents

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