ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/MitAna/TreeMod/src/OutputMod.cc
Revision: 1.9
Committed: Thu Mar 12 18:24:33 2009 UTC (16 years, 1 month ago) by loizides
Content type: text/plain
Branch: MAIN
CVS Tags: Mit_008pre2
Changes since 1.8: +37 -19 lines
Log Message:
Reskimming now properly done

File Contents

# Content
1 // $Id: OutputMod.cc,v 1.8 2009/03/11 10:07:12 loizides Exp $
2
3 #include "MitAna/TreeMod/interface/OutputMod.h"
4 #include "MitAna/TreeMod/interface/HLTFwkMod.h"
5 #include "MitAna/DataUtil/interface/Debug.h"
6 #include "MitAna/DataTree/interface/Names.h"
7 #include "MitAna/DataUtil/interface/TreeWriter.h"
8 #include "MitAna/TreeMod/interface/TreeBranchLoader.h"
9
10 using namespace mithep;
11 using namespace std;
12
13 ClassImp(mithep::OutputMod)
14
15 //--------------------------------------------------------------------------------------------------
16 OutputMod::OutputMod(const char *name, const char *title) :
17 BaseMod(name,title),
18 fTreeName(Names::gkEvtTreeName),
19 fPrefix("skimtest"),
20 fPathName("."),
21 fMaxSize(1024),
22 fCompLevel(9),
23 fSplitLevel(99),
24 fBranchSize(64*1024),
25 fDoReset(kFALSE),
26 fCheckTamBr(kTRUE),
27 fKeepTamBr(kTRUE),
28 fTreeWriter(0),
29 fEventHeader(0),
30 fAllEventHeader(0),
31 fRunInfo(0),
32 fLaHeader(0),
33 fBranches(0),
34 fNBranchesMax(1024),
35 fRunTree(0),
36 fLATree(0),
37 fAllTree(0),
38 fSkimmedIn(0),
39 fL1Tree(0),
40 fHltTree(0),
41 fRunEntries(0),
42 fOrigL1Entry(-1),
43 fL1Entries(0),
44 fOrigHltEntry(-1),
45 fHltEntries(0),
46 fFileNum(0),
47 fLastWrittenEvt(-1),
48 fLastSeenEvt(-1),
49 fCounter(0)
50 {
51 // Constructor.
52 }
53
54 //--------------------------------------------------------------------------------------------------
55 void OutputMod::BeginRun()
56 {
57 // Create HLT tree if HLTFwkMod is being run.
58
59 if (!HasHLTInfo())
60 return;
61
62 if (!fHltTree) {
63 HLTFwkMod *hm = const_cast<HLTFwkMod*>(GetHltFwkMod());
64 fTreeWriter->AddBranchToTree(hm->HLTTreeName(), hm->HLTTabName(),
65 TClass::GetClass(typeid(*hm->fHLTTab))->GetName(),
66 &(hm->fHLTTab), 32000, 0);
67 fTreeWriter->AddBranchToTree(hm->HLTTreeName(), hm->HLTLabName(),
68 TClass::GetClass(typeid(*hm->fHLTLab))->GetName(),
69 &(hm->fHLTLab), 32000, 0);
70 fTreeWriter->SetAutoFill(hm->HLTTreeName(), 0);
71 fHltTree = fTreeWriter->GetTree(hm->HLTTreeName());
72 }
73 }
74
75 //--------------------------------------------------------------------------------------------------
76 void OutputMod::CheckAndAddBranch(const char *bname, const char *cname)
77 {
78 // Check if the given branch should be kept or dropped.
79
80 if (IsAcceptedBranch(bname))
81 return;
82
83 // populate regular expression list if this was not yet done
84 if (fCmdReList.size() != fCmdList.size()) {
85 for (UInt_t i=0; i<fCmdList.size(); ++i) {
86 const char *ptr = fCmdList.at(i).c_str();
87 fCmdReList.push_back(TRegexp(ptr+5,kTRUE));
88 if (ptr[0]=='k')
89 fCmdDeList.push_back(kTRUE);
90 else
91 fCmdDeList.push_back(kFALSE);
92 }
93 }
94
95 // decide whether given branch name should be kept or dropped
96 TString brname(bname);
97 Bool_t decision = kFALSE;
98 Bool_t decision_found = kFALSE;
99
100 for (UInt_t i=0; i<fCmdList.size(); ++i) {
101 TRegexp &re(fCmdReList.at(i));
102 if (brname.Index(re) == kNPOS)
103 continue;
104 decision = fCmdDeList.at(i);
105 decision_found = kTRUE;
106 }
107
108 if (!decision_found) { // no decision found: still drop branch
109 Warning("CheckAndAddBranch",
110 "No decision found for branch '%s' and class '%s'. Branch therefore dropped!",
111 bname, cname);
112 return;
113 }
114
115 if (!decision) { // drop branch according to request
116 SendError(kWarning, "CheckAndAddBranch",
117 "Dropped branch '%s' and class '%s'.", bname, cname);
118 return;
119 }
120
121 // add branch to accepted branch list
122 SendError(kWarning, "CheckAndAddBranch",
123 "Kept branch '%s' and class '%s'.", bname, cname);
124
125 fBrNameList.push_back(string(bname));
126 fBrClassList.push_back(string(cname));
127
128 // request branch
129 RequestBranch(bname);
130 }
131
132 //--------------------------------------------------------------------------------------------------
133 void OutputMod::CheckAndResolveDep(Bool_t solve)
134 {
135 // Check if TAM has loaded additional branches. If requested try to solve the the dependency
136 // by adding the branch to the list of branches.
137
138 const THashTable &ht = GetSel()->GetBranchTable();
139
140 TIter iter(ht.MakeIterator());
141 const TAMBranchInfo *next = dynamic_cast<const TAMBranchInfo*>(iter.Next());
142
143 while (next) {
144 const TAMBranchInfo *cur = next;
145 next = dynamic_cast<const TAMBranchInfo*>(iter.Next());
146 Bool_t isloaded = cur->IsLoaded();
147 if (!isloaded)
148 continue;
149
150 const char *bname = cur->GetName();
151 if (IsAcceptedBranch(bname))
152 continue;
153
154 TreeBranchLoader *loader = dynamic_cast<TreeBranchLoader*>(cur->GetLoader());
155 if (!loader)
156 continue;
157
158 TBranch *br = loader->GetBranch();
159 if (!br)
160 continue;
161
162 const char *cname = br->GetClassName();
163
164 if (solve) {
165 SendError(kWarning, "CheckAndResolveDep",
166 "Resolving dependency for loaded branch %s and class %s", bname,cname);
167
168 fBrNameList.push_back(string(bname));
169 fBrClassList.push_back(string(cname));
170 fBranches[GetNBranches()-1] = reinterpret_cast<TObject*>(loader->GetAddress());
171
172 } else {
173 Warning("CheckAndResolveDep", "Unresolved dependency for loaded branch %s and class %s",
174 bname,cname);
175 }
176 }
177 }
178
179 //--------------------------------------------------------------------------------------------------
180 void OutputMod::EndRun()
181 {
182 // Nothing to be done at this point.
183 }
184
185 //--------------------------------------------------------------------------------------------------
186 void OutputMod::FillAllEventHeader(Bool_t isremoved)
187 {
188 // Fill event header into the all-event-header tree.
189
190 if (!fTreeWriter->BeginEvent(kFALSE)) {
191 SendError(kAbortAnalysis, "FillAllEventHeader", "Begin event failed!");
192 return;
193 }
194
195 if (fSkimmedIn) { // copy alread skimmed headers if any there
196 for(UInt_t i=0; i<fSkimmedIn->Entries(); ++i) {
197 const EventHeader *eh = fSkimmedIn->At(i);
198 fAllEventHeader->SetEvtNum(eh->EvtNum());
199 fAllEventHeader->SetLumiSec(eh->LumiSec());
200 fAllEventHeader->SetRunNum(eh->RunNum());
201 fAllEventHeader->SetRunEntry(eh->RunEntry());
202 fAllEventHeader->SetSkimmed(eh->Skimmed()+1);
203 fAllTree->Fill();
204 }
205 }
206
207 const EventHeader *eh = GetEventHeader();
208 fAllEventHeader->SetEvtNum(eh->EvtNum());
209 fAllEventHeader->SetLumiSec(eh->LumiSec());
210 fAllEventHeader->SetRunNum(eh->RunNum());
211 if (isremoved) {
212 fAllEventHeader->SetRunEntry(-1);
213 fAllEventHeader->SetSkimmed(eh->Skimmed()+1);
214 } else {
215 fAllEventHeader->SetRunEntry(eh->RunEntry());
216 fAllEventHeader->SetSkimmed(eh->Skimmed());
217 }
218
219 fAllTree->Fill();
220 }
221
222 //--------------------------------------------------------------------------------------------------
223 void OutputMod::FillL1Info()
224 {
225 // Not doing anything here until the production writes out L1 information.
226
227 if (!fL1Tree)
228 return;
229 }
230
231 //--------------------------------------------------------------------------------------------------
232 void OutputMod::FillHltInfo()
233 {
234 // Write HLT trigger table if needed.
235
236 if (!fHltTree)
237 return;
238
239 if (fOrigHltEntry == GetHltFwkMod()->fCurEnt)
240 return;
241
242 fHltTree->Fill();
243 fOrigHltEntry = GetHltFwkMod()->fCurEnt;
244 ++fHltEntries;
245 }
246
247 //--------------------------------------------------------------------------------------------------
248 Bool_t OutputMod::IsAcceptedBranch(const char *bname)
249 {
250 // Return true if given branch is already in branch list. Also return true if a special
251 // branch like the "EventHeader" branch is reqested.
252
253 // search in branch list
254 for (UInt_t i=0; i<GetNBranches(); ++i) {
255 if (fBrNameList.at(i).compare(bname) == 0)
256 return kTRUE;
257 }
258
259 // check if special branch that we take care of ourselves
260 string name(bname);
261 if (name.compare("EventHeader") == 0) {
262 return kTRUE;
263 }
264
265 return kFALSE;
266 }
267
268 //--------------------------------------------------------------------------------------------------
269 Bool_t OutputMod::Notify()
270 {
271 // On first notify, loop over list of branches to determine the list of kept branches.
272
273 if (GetNEventsProcessed() != 0)
274 return kTRUE;
275
276 TTree *tree=const_cast<TTree*>(GetSel()->GetTree());
277 if (!tree)
278 return kFALSE;
279
280 TObjArray *arr = tree->GetTree()->GetListOfBranches();
281 if (!arr)
282 return kFALSE;
283
284 for (Int_t i=0; i<arr->GetEntries(); ++i) {
285 TBranch *br = dynamic_cast<TBranch*>(arr->At(i));
286 if (!br && !br->GetMother())
287 continue;
288 br = br->GetMother();
289 TClass *cls = TClass::GetClass(br->GetClassName());
290 if (!cls)
291 continue;
292
293 if (!cls->InheritsFrom("TObject")) {
294 Warning("Notify", "Found branch %s where class %s does not derive from TObject.",
295 br->GetName(), br->GetClassName());
296 continue;
297 }
298
299 CheckAndAddBranch(br->GetName(), br->GetClassName());
300 }
301
302 return kTRUE;
303 }
304
305 //--------------------------------------------------------------------------------------------------
306 void OutputMod::LoadBranches()
307 {
308 // Loop over requested branches and load them.
309
310 for (UInt_t i=0; i<GetNBranches(); ++i) {
311 LoadBranch(fBrNameList.at(i).c_str());
312 }
313 }
314
315 //--------------------------------------------------------------------------------------------------
316 void OutputMod::Process()
317 {
318 // Write out the kept branches of the current event. Make sure the meta information is
319 // correctly updated.
320
321 if (GetSel()->GetCurEvt() == fLastSeenEvt) {
322 Warning("Process", "Event with %ul already seen", fLastSeenEvt);
323 return;
324 }
325 fLastSeenEvt = GetSel()->GetCurEvt();
326
327 if (GetSel()->GetCurEvt() == fLastWrittenEvt) {
328 Warning("Process", "Event with %ul already written", fLastWrittenEvt);
329 return;
330 }
331 fLastWrittenEvt = GetSel()->GetCurEvt();
332 ++fCounter;
333
334 // prepare for tree filling
335 if (!fTreeWriter->BeginEvent(fDoReset)) {
336 SendError(kAbortAnalysis, "Process", "Begin event failed!");
337 return;
338 }
339
340 if (GetNEventsProcessed() == 0 && fCheckTamBr) {
341 CheckAndResolveDep(fKeepTamBr);
342 }
343
344 // load all our branches
345 LoadBranches();
346
347 // pass our branches to tree writer if on first event
348 if (GetNEventsProcessed() == 0) {
349 SetupBranches();
350 }
351
352 // reset per file quantities if a new file was opened
353 if (fTreeWriter->GetFileNumber()!=fFileNum) {
354 fRunmap.clear();
355 fRunEntries = 0;
356 fL1Entries = -1;
357 fHltEntries = -1;
358 fFileNum = fTreeWriter->GetFileNumber();
359 }
360
361 UInt_t runnum = GetEventHeader()->RunNum();
362
363 // store look ahead information
364 if (fRunEntries>0) {
365 fLaHeader->SetRunNum(runnum);
366 fLATree->Fill();
367 }
368
369 // fill event header
370 fEventHeader->SetEvtNum(GetEventHeader()->EvtNum());
371 fEventHeader->SetLumiSec(GetEventHeader()->LumiSec());
372 fEventHeader->SetRunNum(runnum);
373
374 // fill all event header
375 FillAllEventHeader(kFALSE);
376
377 // look-up if entry is in map
378 map<UInt_t,Int_t>::iterator riter = fRunmap.find(runnum);
379 if (riter != fRunmap.end()) { // found existing run info
380 Int_t runentry = riter->second;
381 fEventHeader->SetRunEntry(runentry);
382
383 IncNEventsProcessed();
384 fTreeWriter->EndEvent(fDoReset);
385 return;
386 }
387
388 // fill new run info
389 Int_t runentry = fRunEntries;
390 ++fRunEntries;
391 fEventHeader->SetRunEntry(runentry);
392 fRunmap.insert(pair<UInt_t,Int_t>(runnum,runentry));
393 fRunInfo->SetRunNum(runnum);
394
395 Int_t l1entry = fL1Entries;
396 FillL1Info();
397 fRunInfo->SetL1Entry(l1entry);
398
399 Int_t hltentry = fHltEntries;
400 FillHltInfo();
401 fRunInfo->SetHltEntry(hltentry);
402
403 fRunTree->Fill();
404
405 IncNEventsProcessed();
406
407 if (!fTreeWriter->EndEvent(fDoReset)) {
408 SendError(kAbortAnalysis, "Process", "End event failed!");
409 return;
410 }
411 }
412
413 //--------------------------------------------------------------------------------------------------
414 void OutputMod::ProcessAll()
415 {
416 // Called by the Selector class for events that were skipped.
417
418 if (GetSel()->GetCurEvt() == fLastSeenEvt)
419 return;
420
421 fLastSeenEvt = GetSel()->GetCurEvt();
422 ++fCounter;
423
424 // prepare for tree filling
425 FillAllEventHeader(kTRUE);
426 }
427
428 //--------------------------------------------------------------------------------------------------
429 void OutputMod::RequestBranch(const char *bname)
430 {
431 // Request given branch from TAM.
432
433 if (GetNBranches()>=fNBranchesMax) {
434 SendError(kAbortAnalysis, "RequestBranch", "Can not request branch for %bname"
435 "since maximum number of branches [%d] is reached", bname, fNBranchesMax);
436 return;
437 }
438
439 fBranches[GetNBranches()-1] = 0;
440 TAModule::ReqBranch(bname, fBranches[GetNBranches()-1]);
441 }
442
443 //--------------------------------------------------------------------------------------------------
444 void OutputMod::SetupBranches()
445 {
446 // Setup branches in tree writer.
447
448 for (UInt_t i=0; i<GetNBranches(); ++i) {
449 const char *bname = fBrNameList.at(i).c_str();
450 const char *cname = fBrClassList.at(i).c_str();
451 if (!fBranches[i]) {
452 SendError(kWarning, "SetupBranches",
453 "Pointer for branch '%s' and class '%s' is NULL.", bname, cname);
454 continue;
455 }
456 fTreeWriter->AddBranch(bname, cname, &fBranches[i]);
457 }
458 }
459
460 //--------------------------------------------------------------------------------------------------
461 void OutputMod::SlaveBegin()
462 {
463 // Setup the tree writer and create branches that can already be created at this point.
464
465 // setup tree writer
466 fTreeWriter = new TreeWriter(fTreeName, kFALSE);
467 fTreeWriter->SetBaseURL(fPathName);
468 fTreeWriter->SetPrefix(fPrefix);
469 fTreeWriter->SetMaxSize(fMaxSize*1024*1024);
470 fTreeWriter->SetCompressLevel(fCompLevel);
471 fTreeWriter->SetDefaultSL(fSplitLevel);
472 fTreeWriter->SetDefaultBrSize(fBranchSize);
473 fTreeWriter->AddTree(fTreeName);
474 fTreeWriter->DoBranchRef(fTreeName);
475
476 // deal with my own tree objects
477 fEventHeader = new EventHeader;
478 fTreeWriter->AddBranch(GetSel()->GetEvtHdrName(), &fEventHeader);
479
480 // deal with other trees
481 const char *tname = 0;
482 fRunInfo = new RunInfo;
483 tname = GetSel()->GetRunTreeName();
484 fTreeWriter->AddBranchToTree(tname, GetSel()->GetRunInfoName(), &fRunInfo);
485 fTreeWriter->SetAutoFill(tname, 0);
486 fRunTree = fTreeWriter->GetTree(tname);
487 fLaHeader = new LAHeader;
488 tname = GetSel()->GetLATreeName();
489 fTreeWriter->AddBranchToTree(tname, GetSel()->GetLAHdrName(), &fLaHeader);
490 fTreeWriter->SetAutoFill(tname, 0);
491 fLATree = fTreeWriter->GetTree(tname);
492 fAllEventHeader = new EventHeader;
493 tname = GetSel()->GetAllEvtTreeName();
494 fTreeWriter->AddBranchToTree(tname, GetSel()->GetAllEvtHdrBrn(), &fAllEventHeader);
495 fAllTree = fTreeWriter->GetTree(tname);
496 fTreeWriter->SetAutoFill(tname, 0);
497
498 // get pointer to all event headers
499 fSkimmedIn = GetPublicObj<EventHeaderCol>(Names::gkSkimmedHeaders);
500
501 // deal here with published objects
502 // todo
503
504 // create TObject space for TAM
505 fBranches = new TObject*[fNBranchesMax];
506
507 // adjust checks for TAM branches
508 if (fKeepTamBr)
509 fCheckTamBr = kTRUE;
510 }
511
512 //--------------------------------------------------------------------------------------------------
513 void OutputMod::SlaveTerminate()
514 {
515 // Terminate tree writing and do cleanup.
516
517 RetractObj(Names::gkSkimmedHeaders);
518
519 delete fTreeWriter;
520 fTreeWriter = 0;
521
522 delete fEventHeader;
523 delete fRunInfo;
524 delete fLaHeader;
525 delete fAllEventHeader;
526
527 delete[] fBranches;
528
529 Double_t frac = 100.*GetNEventsProcessed()/fCounter;
530 SendError(kWarning, "SlaveTerminate", "Stored %.2g%% events (%ld out of %ld)",
531 frac, GetNEventsProcessed(), fCounter);
532 }