ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/MitAna/TreeMod/src/OutputMod.cc
Revision: 1.20
Committed: Wed Mar 28 12:15:38 2012 UTC (13 years, 1 month ago) by paus
Content type: text/plain
Branch: MAIN
CVS Tags: Mit_032, Mit_031, Mit_030, Mit_029c, Mit_029b, Mit_030_pre1, Mit_029a, Mit_029, Mit_029_pre1, Mit_028a, Mit_028, Mit_027a, Mit_027, Mit_026, Mit_025e, HEAD
Changes since 1.19: +64 -60 lines
Log Message:
Enable skimming.

File Contents

# Content
1 // $Id: OutputMod.cc,v 1.19 2011/03/11 04:03:54 bendavid 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/BranchTable.h"
7 #include "MitAna/DataTree/interface/EventHeaderCol.h"
8 #include "MitAna/DataTree/interface/Names.h"
9 #include "MitAna/DataUtil/interface/TreeWriter.h"
10 #include "MitAna/TreeMod/interface/TreeBranchLoader.h"
11
12 #include "MitAna/DataTree/interface/PhotonCol.h"
13
14 using namespace mithep;
15 using namespace std;
16
17 ClassImp(mithep::OutputMod)
18
19 //--------------------------------------------------------------------------------------------------
20 OutputMod::OutputMod(const char *name, const char *title) :
21 BaseMod (name,title),
22 fTreeName (Names::gkEvtTreeName),
23 fPrefix ("skimtest"),
24 fPathName ("."),
25 fMaxSize (1024),
26 fCompLevel (9),
27 fSplitLevel (99),
28 fBranchSize (16*1024),
29 fDoReset (kFALSE),
30 fCheckBrDep (kTRUE),
31 fUseBrDep (kTRUE),
32 fCheckTamBr (kTRUE),
33 fKeepTamBr (kFALSE),
34 fTreeWriter (0),
35 fEventHeader (0),
36 fAllEventHeader(0),
37 fRunInfo (0),
38 fLaHeader (0),
39 fBranchTable (0),
40 fBranches (0),
41 fNBranchesMax (1024),
42 fRunTree (0),
43 fLATree (0),
44 fAllTree (0),
45 fSkimmedIn (0),
46 fHltTree (0),
47 fHLTTab (new vector<string>),
48 fHLTLab (new vector<string>),
49 fRunEntries (0),
50 fHltEntries (0),
51 fFileNum (-1),
52 fLastWrittenEvt(-1),
53 fLastSeenEvt (-1),
54 fCounter (0),
55 fGoodPhotons (0)
56 {
57 // Constructor.
58 }
59
60 //--------------------------------------------------------------------------------------------------
61 void OutputMod::BeginRun()
62 {
63 // Create HLT tree if HLTFwkMod is being run.
64
65 if (!HasHLTInfo())
66 return;
67
68 if (!fHltTree) {
69 HLTFwkMod *hm = const_cast<HLTFwkMod*>(GetHltFwkMod());
70 fTreeWriter->AddBranchToTree(hm->HLTTreeName(), hm->HLTTabName(),
71 TClass::GetClass(typeid(*fHLTTab))->GetName(),
72 &fHLTTab, 32000, 0);
73 fTreeWriter->AddBranchToTree(hm->HLTTreeName(), hm->HLTLabName(),
74 TClass::GetClass(typeid(*fHLTLab))->GetName(),
75 &fHLTLab, 32000, 0);
76 fTreeWriter->SetAutoFill(hm->HLTTreeName(), 0);
77 fHltTree = fTreeWriter->GetTree(hm->HLTTreeName());
78 }
79 }
80
81 //--------------------------------------------------------------------------------------------------
82 void OutputMod::CheckAndAddBranch(const char *bname, const char *cname)
83 {
84 // Check if the given branch should be kept or dropped.
85
86 if (IsAcceptedBranch(bname))
87 return;
88
89 // populate regular expression list if this was not yet done
90 if (fCmdReList.size() != fCmdList.size()) {
91 for (UInt_t i=0; i<fCmdList.size(); ++i) {
92 const char *ptr = fCmdList.at(i).c_str();
93 fCmdReList.push_back(TRegexp(ptr+5,kTRUE));
94 if (ptr[0]=='k')
95 fCmdDeList.push_back(kTRUE);
96 else
97 fCmdDeList.push_back(kFALSE);
98 }
99 }
100
101 // decide whether given branch name should be kept or dropped
102 TString brname(bname);
103 Bool_t decision = kFALSE;
104 Bool_t decision_found = kFALSE;
105
106 for (UInt_t i=0; i<fCmdList.size(); ++i) {
107 TRegexp &re(fCmdReList.at(i));
108 if (brname.Index(re) == kNPOS)
109 continue;
110 decision = fCmdDeList.at(i);
111 decision_found = kTRUE;
112 }
113
114 if (!decision_found) { // no decision found: still drop branch
115 Warning("CheckAndAddBranch",
116 "No decision found for branch '%s' and class '%s'. Branch therefore dropped!",
117 bname, cname);
118 return;
119 }
120
121 if (!decision) { // drop branch according to request
122 Info("CheckAndAddBranch",
123 "Dropped branch '%s' and class '%s'", bname, cname);
124 return;
125 }
126
127 // add branch to accepted branch list
128 Info("CheckAndAddBranch",
129 "Kept branch '%s' and class '%s'", bname, cname);
130
131 fBrNameList .push_back(string(bname));
132 fBrClassList.push_back(string(cname));
133 }
134
135 //--------------------------------------------------------------------------------------------------
136 Bool_t OutputMod::CheckAndResolveBranchDep()
137 {
138 // Checks dependency in BranchTable. Resolve dependency automatically if fUserBrDep is kTRUE.
139
140 TFile *cfile = const_cast<TFile*>(GetSel()->GetCurrentFile());
141 if (!cfile) {
142 SendError(kAbortAnalysis, "CheckAndResolveBranchDep", "Could not get pointer to current file!");
143 return kFALSE;
144 }
145
146 BranchTable *br = dynamic_cast<BranchTable*>(cfile->Get(Names::gkBranchTable));
147 if (!br) {
148 SendError(kAbortAnalysis, "CheckAndResolveBranchDep", "Could not get pointer to branch table!");
149 return kFALSE;
150 }
151
152 TList *blist = br->GetBranches();
153 if (!blist) {
154 SendError(kAbortAnalysis, "CheckAndResolveBranchDep", "Could not get list of branches!");
155 return kFALSE;
156 }
157
158 fBranchTable = new BranchTable;
159 fBranchTable->SetName(Names::gkBranchTable);
160 fBranchTable->SetOwner();
161
162 TList sht;
163 sht.SetOwner(kTRUE);
164 for (UInt_t i=0; i<fBrNameList.size(); ++i) {
165 sht.Add(new TObjString(fBrNameList.at(i).c_str()));
166 }
167
168 for (UInt_t i=0; i<fBrNameList.size(); ++i) {
169 TString brname(fBrNameList.at(i));
170 if (!blist->FindObject(brname))
171 continue;
172 TList *bdeps = br->GetDepBranches(brname);
173 if (!bdeps)
174 continue;
175
176 // check dependency
177 TIter iter(bdeps->MakeIterator());
178 const TObjString *n = dynamic_cast<const TObjString*>(iter.Next());
179 while (n) {
180 if (sht.FindObject(n->GetName())) {
181 // dependent branch is already accepted
182 fBranchTable->Add(new BranchName(brname,n->GetName()));
183 } else {
184 if (fUseBrDep) {
185 const TObjArray *arr = GetSel()->GetTree()->GetTree()->GetListOfBranches();
186 TBranch *br = dynamic_cast<TBranch*>(arr->FindObject(n->GetName()));
187 if (!br) {
188 Error("CheckAndResolveBranchDep",
189 "Could not get branch '%s' to resolve dependency for branch '%s'",
190 n->GetName(), brname.Data());
191 } else {
192 Info("CheckAndResolveBranchDep",
193 "Adding branch '%s' to resolve dependency for branch '%s'",
194 n->GetName(), brname.Data());
195 fBrNameList .push_back(string(n->GetName()));
196 fBrClassList.push_back(br->GetClassName());
197 sht.Add(new TObjString(n->GetName()));
198 fBranchTable->Add(new BranchName(brname,n->GetName()));
199 }
200 } else {
201 Warning("CheckAndResolveBranchDep",
202 "Unresolved dependency of branch '%s' and '%s' ",
203 n->GetName(), brname.Data());
204 }
205 }
206 n = dynamic_cast<const TObjString*>(iter.Next());
207 }
208 delete bdeps;
209 }
210 delete blist;
211 return kTRUE;
212 }
213
214 //--------------------------------------------------------------------------------------------------
215 void OutputMod::CheckAndResolveTAMDep(Bool_t solve)
216 {
217 // Check if TAM has loaded additional branches. If requested try to solve the the dependency
218 // by adding the branch to the list of branches.
219
220 const THashTable &ht = GetSel()->GetBranchTable();
221
222 TIter iter(ht.MakeIterator());
223 const TAMBranchInfo *next = dynamic_cast<const TAMBranchInfo*>(iter.Next());
224
225 while (next) {
226 const TAMBranchInfo *cur = next;
227 next = dynamic_cast<const TAMBranchInfo*>(iter.Next());
228 Bool_t isloaded = cur->IsLoaded();
229 if (!isloaded)
230 continue;
231
232 const char *bname = cur->GetName();
233 if (IsAcceptedBranch(bname))
234 continue;
235
236 TreeBranchLoader *loader = dynamic_cast<TreeBranchLoader*>(cur->GetLoader());
237 if (!loader)
238 continue;
239
240 TBranch *br = loader->GetBranch();
241 if (!br)
242 continue;
243
244 const char *cname = br->GetClassName();
245
246 if (solve) {
247 Info("CheckAndResolveTAMDep",
248 "Resolving dependency for loaded branch '%s' and class '%s'", bname,cname);
249
250 fBrNameList. push_back(string(bname));
251 fBrClassList.push_back(string(cname));
252 fBranches[GetNBranches()-1] = reinterpret_cast<TObject*>(loader->GetAddress());
253
254 } else {
255 Warning("CheckAndResolveTAMDep",
256 "Unresolved dependency for loaded branch '%s' and class '%s'",
257 bname,cname);
258 }
259 }
260 }
261
262 //--------------------------------------------------------------------------------------------------
263 void OutputMod::EndRun()
264 {
265 // Nothing to be done at this point.
266 }
267
268 //--------------------------------------------------------------------------------------------------
269 void OutputMod::FillAllEventHeader(Bool_t isremoved)
270 {
271 // Fill event header into the all-event-header tree.
272
273 if (!fTreeWriter->BeginEvent(kFALSE)) {
274 SendError(kAbortAnalysis, "FillAllEventHeader", "Begin event failed!");
275 return;
276 }
277
278 if (fSkimmedIn) { // copy alread skimmed headers if any there
279 const UInt_t n = fSkimmedIn->GetEntries();
280 for (UInt_t i=0; i<n; ++i) {
281 const EventHeader *eh = fSkimmedIn->At(i);
282 *fAllEventHeader = *eh;
283 fAllEventHeader->SetSkimmed(eh->Skimmed()+1);
284 fAllTree->Fill();
285 }
286 }
287
288 const EventHeader *eh = GetEventHeader();
289 *fAllEventHeader = *eh;
290 if (isremoved) {
291 fAllEventHeader->SetRunEntry(-1);
292 fAllEventHeader->SetSkimmed(eh->Skimmed()+1);
293 } else {
294 fAllEventHeader->SetRunEntry(eh->RunEntry());
295 fAllEventHeader->SetSkimmed(eh->Skimmed());
296 }
297
298 fAllTree->Fill();
299 }
300
301 //--------------------------------------------------------------------------------------------------
302 void OutputMod::FillHltInfo()
303 {
304 // Write HLT trigger table if needed.
305
306 if (!fHltTree)
307 return;
308
309 HLTFwkMod *hm = const_cast<HLTFwkMod*>(GetHltFwkMod());
310 vector<string> *trigtable = hm->fHLTTab;
311 vector<string> *labels = hm->fHLTLab;
312
313 Bool_t doCopy = kFALSE;
314 if (fHLTTab->size()==0) {
315 doCopy = kTRUE;
316 } else {
317 // check if existing table contains all necessary paths:
318 // if so keep it, otherwise store the new one
319
320 if ((fHLTTab->size() != trigtable->size()) ||
321 (fHLTLab->size() != labels->size())) {
322 doCopy = kTRUE;
323 } else {
324 // need to check more thoroughly
325
326 for (UInt_t i=0; i<trigtable->size(); ++i) {
327 if (trigtable->at(i) != fHLTTab->at(i)) {
328 doCopy = kTRUE;
329 break;
330 }
331 }
332 if (!doCopy) {
333 for (UInt_t i=0; i<labels->size(); ++i) {
334 if (labels->at(i) != fHLTLab->at(i)) {
335 doCopy = kTRUE;
336 break;
337 }
338 }
339 }
340 }
341 }
342
343 if (!doCopy)
344 return;
345
346 fHLTTab->resize(trigtable->size());
347 copy(trigtable->begin(),trigtable->end(), fHLTTab->begin());
348 fHLTLab->resize(labels->size());
349 copy(labels->begin(),labels->end(), fHLTLab->begin());
350
351 ++fHltEntries;
352 fHltTree->Fill();
353 }
354
355 //--------------------------------------------------------------------------------------------------
356 Bool_t OutputMod::IsAcceptedBranch(const char *bname)
357 {
358 // Return true if given branch is already in branch list. Also return true if a special
359 // branch like the "EventHeader" branch is reqested.
360
361 // search in branch list
362 for (UInt_t i=0; i<GetNBranches(); ++i) {
363 if (fBrNameList.at(i).compare(bname) == 0)
364 return kTRUE;
365 }
366
367 // check if special branch that we take care of ourselves
368 string name(bname);
369 if (name.compare("EventHeader") == 0) {
370 return kTRUE;
371 }
372
373 return kFALSE;
374 }
375
376 //--------------------------------------------------------------------------------------------------
377 Bool_t OutputMod::Notify()
378 {
379 // On first notify, loop over list of branches to determine the list of kept branches.
380
381 if (GetNEventsProcessed() != 0)
382 return kTRUE;
383
384 const TTree *tree=GetSel()->GetTree();
385 if (!tree)
386 return kFALSE;
387
388 const TObjArray *arr = tree->GetTree()->GetListOfBranches();
389 if (!arr)
390 return kFALSE;
391
392 for (Int_t i=0; i<arr->GetEntries(); ++i) {
393 TBranch *br = dynamic_cast<TBranch*>(arr->At(i));
394 if (!br && !br->GetMother())
395 continue;
396 br = br->GetMother();
397 TClass *cls = TClass::GetClass(br->GetClassName());
398 if (!cls)
399 continue;
400
401 if (!cls->InheritsFrom("TObject")) {
402 Warning("Notify", "Found branch '%s' where class '%s' does not derive from TObject",
403 br->GetName(), br->GetClassName());
404 continue;
405 }
406
407 CheckAndAddBranch(br->GetName(), br->GetClassName());
408 }
409
410 if (fCheckBrDep && !CheckAndResolveBranchDep())
411 return kFALSE;
412
413 RequestBranches();
414 return kTRUE;
415 }
416
417 //--------------------------------------------------------------------------------------------------
418 void OutputMod::LoadBranches()
419 {
420 // Loop over requested branches and load them.
421
422 for (UInt_t i=0; i<GetNBranches(); ++i)
423 LoadBranch(fBrNameList.at(i).c_str());
424 }
425
426 //--------------------------------------------------------------------------------------------------
427 void OutputMod::Process()
428 {
429 // Write out the kept branches of the current event. Make sure the meta information is
430 // correctly updated.
431
432 if (GetSel()->GetCurEvt() == fLastSeenEvt) {
433 Warning("Process", "Event with %lld already seen", fLastSeenEvt);
434 return;
435 }
436 fLastSeenEvt = GetSel()->GetCurEvt();
437
438 if (GetSel()->GetCurEvt() == fLastWrittenEvt) {
439 Warning("Process", "Event with %lld already written", fLastWrittenEvt);
440 return;
441 }
442 fLastWrittenEvt = GetSel()->GetCurEvt();
443 ++fCounter;
444
445 // prepare for tree filling
446 if (!fTreeWriter->BeginEvent(fDoReset)) {
447 SendError(kAbortAnalysis, "Process", "Begin event failed!");
448 return;
449 }
450
451 if (GetNEventsProcessed() == 0 && fCheckTamBr) {
452 CheckAndResolveTAMDep(fKeepTamBr);
453 }
454
455 // load all our branches
456 LoadBranches();
457
458 // pass our branches to tree writer if on first event
459 if (GetNEventsProcessed() == 0) {
460 SetupBranches();
461 }
462
463 // reset per file quantities if a new file was opened
464 if (fTreeWriter->GetFileNumber()!=fFileNum) {
465 fRunmap.clear();
466 fHLTTab->clear();
467 fHLTLab->clear();
468 fRunEntries = 0;
469 fHltEntries = 0;
470 fFileNum = fTreeWriter->GetFileNumber();
471 if (fBranchTable)
472 fTreeWriter->StoreObject(fBranchTable);
473 }
474
475 UInt_t runnum = GetEventHeader()->RunNum();
476
477 // store look ahead information
478 if (fRunEntries>0) {
479 fLaHeader->SetRunNum(runnum);
480 fLATree->Fill();
481 }
482
483 // fill event header
484 *fEventHeader = *GetEventHeader();
485
486 // fill all event header
487 FillAllEventHeader(kFALSE);
488
489 // look-up if entry is in map
490 map<UInt_t,Int_t>::iterator riter = fRunmap.find(runnum);
491 if (riter != fRunmap.end()) { // found existing run info
492 Int_t runEntry = riter->second;
493 fEventHeader->SetRunEntry(runEntry);
494
495 IncNEventsProcessed();
496 fTreeWriter->EndEvent(fDoReset);
497 return;
498 }
499
500 // fill new run info
501 Int_t runEntry = fRunEntries;
502 ++fRunEntries;
503 fEventHeader->SetRunEntry(runEntry);
504 fRunmap.insert(pair<UInt_t,Int_t>(runnum,runEntry));
505 fRunInfo->SetRunNum(runnum);
506
507 Int_t hltEntry = fHltEntries;
508 FillHltInfo();
509 if (hltEntry < fHltEntries)
510 fRunInfo->SetHltEntry(hltEntry);
511 else
512 fRunInfo->SetHltEntry(hltEntry-1);
513
514 fRunTree->Fill();
515
516 IncNEventsProcessed();
517
518 if (!fTreeWriter->EndEvent(fDoReset)) {
519 SendError(kAbortAnalysis, "Process", "End event failed!");
520 return;
521 }
522 }
523
524 //--------------------------------------------------------------------------------------------------
525 void OutputMod::ProcessAll()
526 {
527 // Called by the Selector class for events that were skipped.
528
529 if (GetSel()->GetCurEvt() == fLastSeenEvt)
530 return;
531
532 fLastSeenEvt = GetSel()->GetCurEvt();
533 ++fCounter;
534
535 // prepare for tree filling
536 FillAllEventHeader(kTRUE);
537 }
538
539 //--------------------------------------------------------------------------------------------------
540 void OutputMod::RequestBranches()
541 {
542 // Loop over requested branches and request them.
543
544 for (UInt_t i=0; i<GetNBranches(); ++i) {
545 if (i>=fNBranchesMax) {
546 SendError(kAbortAnalysis, "RequestBranches", "Cannot request branch '%s' "
547 "since maximum number of branches [%d] is reached",
548 fBrNameList.at(i).c_str(), fNBranchesMax);
549 return;
550 }
551 fBranches[i] = 0;
552 TAModule::ReqBranch(fBrNameList.at(i).c_str(), fBranches[i]);
553 }
554 }
555
556 //--------------------------------------------------------------------------------------------------
557 void OutputMod::SetupBranches()
558 {
559 // Setup branches in tree writer.
560
561 for (UInt_t i=0; i<GetNBranches(); ++i) {
562 const char *bname = fBrNameList.at(i).c_str();
563 const char *cname = fBrClassList.at(i).c_str();
564 if (!fBranches[i]) {
565 SendError(kWarning, "SetupBranches",
566 "Pointer for branch '%s' and class '%s' is NULL", bname, cname);
567 continue;
568 }
569 Int_t bsize = fBranchSize;
570 TString cnamestr(cname);
571 if ((bsize<128*1024) && (cnamestr.Contains("mithep::MCParticle"))) {
572 bsize=128*1024;
573 } else if ((bsize<32*1024) && (cnamestr.Contains("mithep::CaloTower"))) {
574 bsize=32*1024;
575 }
576
577 fTreeWriter->AddBranch(bname, cname, &fBranches[i], bsize);
578 }
579 }
580
581 //--------------------------------------------------------------------------------------------------
582 void OutputMod::SlaveBegin()
583 {
584 // Setup the tree writer and create branches that can already be created at this point.
585
586 // setup tree writer
587 fTreeWriter = new TreeWriter (fTreeName, kFALSE);
588 fTreeWriter->SetBaseURL (fPathName);
589 fTreeWriter->SetPrefix (fPrefix);
590 fTreeWriter->SetMaxSize (fMaxSize*1024*1024);
591 fTreeWriter->SetCompressLevel(fCompLevel);
592 fTreeWriter->SetDefaultSL (fSplitLevel);
593 fTreeWriter->SetDefaultBrSize(fBranchSize);
594 fTreeWriter->AddTree (fTreeName);
595 fTreeWriter->DoBranchRef (fTreeName);
596
597 // deal with my own tree objects
598 fEventHeader = new EventHeader;
599 fTreeWriter->AddBranch(GetSel()->GetEvtHdrName(), &fEventHeader);
600
601 // deal with other trees
602 const char *tname = 0;
603 fRunInfo = new RunInfo;
604 tname = GetSel()->GetRunTreeName();
605 fTreeWriter->AddBranchToTree(tname, GetSel()->GetRunInfoName(), &fRunInfo);
606 fTreeWriter->SetAutoFill(tname, 0);
607 // the run tree first
608 fRunTree = fTreeWriter->GetTree(tname);
609 fLaHeader = new LAHeader;
610 tname = GetSel()->GetLATreeName();
611 fTreeWriter->AddBranchToTree(tname, GetSel()->GetLAHdrName(), &fLaHeader);
612 fTreeWriter->SetAutoFill(tname, 0);
613 // the Look Ahead tree next
614 fLATree = fTreeWriter->GetTree(tname);
615 fAllEventHeader = new EventHeader;
616 tname = GetSel()->GetAllEvtTreeName();
617 fTreeWriter->AddBranchToTree(tname, GetSel()->GetAllEvtHdrBrn(), &fAllEventHeader);
618 fAllTree = fTreeWriter->GetTree(tname);
619 fTreeWriter->SetAutoFill(tname, 0);
620
621 // get pointer to all event headers
622 fSkimmedIn = GetPublicObj<EventHeaderCol>(Names::gkSkimmedHeaders);
623
624 // create TObject space for TAM
625 fBranches = new TObject*[fNBranchesMax + fAddList.size()];
626
627 // deal here with additional published objects
628 for (UInt_t i=0; i<fAddList.size(); ++i) {
629 TString objname(fAddList.at(i));
630 TObject *obj = FindPublicObj(objname);
631 if (obj) {
632 fBranches[fNBranchesMax+i] = obj;
633 fTreeWriter->AddBranch(objname, &fBranches[fNBranchesMax+i], 64*1024, 0);
634 Info("SlaveBegin", "Adding additional branch named '%s' as requested", objname.Data());
635 } else {
636 SendError(kAbortAnalysis, "SlaveBegin",
637 "Object named '%s' for additional branch is NULL", objname.Data());
638 }
639 }
640
641 // adjust checks for TAM branches
642 if (fKeepTamBr)
643 fCheckTamBr = kTRUE;
644 }
645
646 //--------------------------------------------------------------------------------------------------
647 void OutputMod::SlaveTerminate()
648 {
649 // Terminate tree writing and do cleanup.
650
651 RetractObj(Names::gkSkimmedHeaders);
652
653 delete fTreeWriter;
654 fTreeWriter = 0;
655
656 delete fEventHeader;
657 delete fRunInfo;
658 delete fLaHeader;
659 delete fAllEventHeader;
660
661 delete[] fBranches;
662
663 Double_t frac = 100.*GetNEventsProcessed()/fCounter;
664 Info("SlaveTerminate", "Stored %.2f%% events (%d out of %lld)",
665 frac, GetNEventsProcessed(), fCounter);
666 }