1 |
// $Id: AnaFwkMod.cc,v 1.13 2009/06/22 15:36:47 loizides Exp $
|
2 |
|
3 |
#include "MitAna/TreeMod/interface/AnaFwkMod.h"
|
4 |
#include "MitAna/DataUtil/interface/Debug.h"
|
5 |
#include "MitAna/DataTree/interface/Names.h"
|
6 |
#include <TFile.h>
|
7 |
#include <TH1D.h>
|
8 |
#include <TStopwatch.h>
|
9 |
#include <TTree.h>
|
10 |
|
11 |
using namespace mithep;
|
12 |
|
13 |
ClassImp(mithep::AnaFwkMod)
|
14 |
|
15 |
//--------------------------------------------------------------------------------------------------
|
16 |
AnaFwkMod::AnaFwkMod(const char *name, const char *title) :
|
17 |
BaseMod(name,title),
|
18 |
fAllHeadTreeName(Names::gkAllEvtTreeName),
|
19 |
fAllHeadBrName(Names::gkAllEvtHeaderBrn),
|
20 |
fSkipNEvents(0),
|
21 |
fPrintScale(1),
|
22 |
fSWtotal(0),
|
23 |
fSWevent(0),
|
24 |
fAllHeaders(0,Names::gkSkimmedHeaders),
|
25 |
fAllHeadTree(0),
|
26 |
fAllEventHeader(0),
|
27 |
fReload(kFALSE),
|
28 |
fCurEnt(-2),
|
29 |
fNEventsSkimmed(0),
|
30 |
fNEventsSkipped(0)
|
31 |
{
|
32 |
// Constructor.
|
33 |
}
|
34 |
|
35 |
//--------------------------------------------------------------------------------------------------
|
36 |
void AnaFwkMod::BeginRun()
|
37 |
{
|
38 |
// Get HLT tree and set branches if new file was opened. Read next entry in HLT key
|
39 |
// depending on entry in RunInfo.
|
40 |
|
41 |
if (fReload) {
|
42 |
// reset to be (re-)loaded variables
|
43 |
fReload = 0;
|
44 |
fAllHeadTree = 0;
|
45 |
fAllEventHeader = 0;
|
46 |
fCurEnt = 0;
|
47 |
|
48 |
// get current file
|
49 |
TFile *file = GetCurrentFile();
|
50 |
if (!file)
|
51 |
return;
|
52 |
|
53 |
// get all event header tree
|
54 |
fAllHeadTree = dynamic_cast<TTree*>(file->Get(fAllHeadTreeName));
|
55 |
if (!fAllHeadTree) {
|
56 |
SendError(kWarning, "BeginRun",
|
57 |
"Cannot find tree '%s' in file '%s'",
|
58 |
fAllHeadTreeName.Data(),file->GetName());
|
59 |
return;
|
60 |
}
|
61 |
|
62 |
// get all event header branch
|
63 |
if (fAllHeadTree->GetBranch(fAllHeadBrName)) {
|
64 |
fAllHeadTree->SetBranchAddress(fAllHeadBrName, &fAllEventHeader);
|
65 |
} else {
|
66 |
SendError(kWarning, "BeginRun",
|
67 |
"Cannot find branch '%s' in tree '%s'",
|
68 |
fAllHeadBrName.Data(), fAllHeadTreeName.Data());
|
69 |
return;
|
70 |
}
|
71 |
}
|
72 |
}
|
73 |
|
74 |
//--------------------------------------------------------------------------------------------------
|
75 |
void AnaFwkMod::CopyAllEventHeaders()
|
76 |
{
|
77 |
// Deal with the headers from all events if needed.
|
78 |
|
79 |
fAllHeaders.Reset();
|
80 |
|
81 |
const EventHeader *curev = GetEventHeader();
|
82 |
if (!curev) {
|
83 |
SendError(kWarning, "CopyAllEventHeaders", "Cannot obtain current event");
|
84 |
return;
|
85 |
}
|
86 |
|
87 |
if (fAllHeadTree) {
|
88 |
const Int_t nemax = fAllHeadTree->GetEntries();
|
89 |
if (fCurEnt == nemax) {
|
90 |
SendError(kAbortEvent, "CopyAllEventHeaders",
|
91 |
"End of all events tree reached (%d=%d)", fCurEnt, nemax);
|
92 |
return;
|
93 |
}
|
94 |
fAllHeadTree->GetEntry(fCurEnt++);
|
95 |
while(fCurEnt<=nemax && fAllEventHeader->Skimmed()) {
|
96 |
EventHeader *eh = fAllHeaders.AddNew();
|
97 |
eh->SetRunNum(fAllEventHeader->RunNum());
|
98 |
eh->SetEvtNum(fAllEventHeader->EvtNum());
|
99 |
eh->SetLumiSec(fAllEventHeader->LumiSec());
|
100 |
eh->SetRunEntry(fAllEventHeader->RunEntry());
|
101 |
eh->SetSkimmed(fAllEventHeader->Skimmed());
|
102 |
fAllHeadTree->GetEntry(fCurEnt++);
|
103 |
}
|
104 |
if ((fAllEventHeader->RunNum()!=curev->RunNum()) ||
|
105 |
(fAllEventHeader->EvtNum()!=curev->EvtNum()) ||
|
106 |
(fAllEventHeader->LumiSec()!=curev->LumiSec()) ||
|
107 |
(fAllEventHeader->RunEntry()!=curev->RunEntry())) {
|
108 |
SendError(kWarning, "CopyAllEventHeaders",
|
109 |
"Event header information for entry %d inconsistent: "
|
110 |
"%d==%d, %d==%d, %d==%d, %d==%d",
|
111 |
fCurEnt,
|
112 |
fAllEventHeader->RunNum(), curev->RunNum(),
|
113 |
fAllEventHeader->EvtNum(), curev->EvtNum(),
|
114 |
fAllEventHeader->LumiSec(), curev->LumiSec(),
|
115 |
fAllEventHeader->RunEntry(), curev->RunEntry());
|
116 |
return;
|
117 |
}
|
118 |
|
119 |
// read-ahead to check if more events are coming
|
120 |
if (fCurEnt<nemax) {
|
121 |
Int_t testEnt = fCurEnt;
|
122 |
fAllHeadTree->GetEntry(testEnt++);
|
123 |
while(testEnt<=nemax && fAllEventHeader->Skimmed())
|
124 |
fAllHeadTree->GetEntry(testEnt++);
|
125 |
if (testEnt==nemax+1) { // need to add remaining skimmed events
|
126 |
fAllHeadTree->GetEntry(fCurEnt++);
|
127 |
while(fCurEnt<=nemax) {
|
128 |
EventHeader *eh = fAllHeaders.AddNew();
|
129 |
eh->SetRunNum(fAllEventHeader->RunNum());
|
130 |
eh->SetEvtNum(fAllEventHeader->EvtNum());
|
131 |
eh->SetLumiSec(fAllEventHeader->LumiSec());
|
132 |
eh->SetRunEntry(fAllEventHeader->RunEntry());
|
133 |
eh->SetSkimmed(fAllEventHeader->Skimmed());
|
134 |
fAllHeadTree->GetEntry(fCurEnt++);
|
135 |
}
|
136 |
if (fCurEnt != nemax+1) {
|
137 |
SendError(kAbortEvent, "CopyAllEventHeaders",
|
138 |
"End of all events tree unexpectedly not reached (%d!=%d)", fCurEnt, nemax);
|
139 |
return;
|
140 |
}
|
141 |
}
|
142 |
}
|
143 |
}
|
144 |
}
|
145 |
|
146 |
//--------------------------------------------------------------------------------------------------
|
147 |
Bool_t AnaFwkMod::Notify()
|
148 |
{
|
149 |
// Make sure to get the new "AllEvents" tree when the file changes.
|
150 |
|
151 |
fReload = kTRUE;
|
152 |
return kTRUE;
|
153 |
}
|
154 |
|
155 |
//--------------------------------------------------------------------------------------------------
|
156 |
void AnaFwkMod::Process()
|
157 |
{
|
158 |
// Do event skipping and counting and print out timing information.
|
159 |
|
160 |
// counting processed events
|
161 |
IncNEventsProcessed();
|
162 |
|
163 |
// get skimmed event headers
|
164 |
CopyAllEventHeaders();
|
165 |
fNEventsSkimmed += fAllHeaders.GetEntries();
|
166 |
|
167 |
// check if events should be skipped
|
168 |
if (fNEventsSkipped<fSkipNEvents) {
|
169 |
++fNEventsSkipped;
|
170 |
MDB(kAnalysis, 3) {
|
171 |
Info("Process", "Skipping (aborting) %d out of %d first events.",
|
172 |
fNEventsSkipped, fSkipNEvents);
|
173 |
}
|
174 |
AbortEvent();
|
175 |
return;
|
176 |
}
|
177 |
|
178 |
// check if printout should be done
|
179 |
Bool_t doPrint = 0;
|
180 |
UInt_t nProcessed = fPrintScale;
|
181 |
|
182 |
MDB(kAnalysis, 4) {
|
183 |
if (GetNEventsProcessed() % (fPrintScale) == 0)
|
184 |
doPrint = 1;
|
185 |
} else {
|
186 |
MDB(kAnalysis, 3) {
|
187 |
if (GetNEventsProcessed() % (fPrintScale*10) == 0)
|
188 |
doPrint = 1;
|
189 |
nProcessed = fPrintScale*10;
|
190 |
} else {
|
191 |
MDB(kAnalysis, 2) {
|
192 |
if (GetNEventsProcessed() % (fPrintScale*100) == 0)
|
193 |
doPrint = 1;
|
194 |
nProcessed = fPrintScale*100;
|
195 |
} else {
|
196 |
MDB(kAnalysis, 1) {
|
197 |
if (GetNEventsProcessed() % (fPrintScale*1000) == 0)
|
198 |
doPrint = 1;
|
199 |
nProcessed = fPrintScale*1000;
|
200 |
}
|
201 |
}
|
202 |
}
|
203 |
}
|
204 |
|
205 |
if (doPrint) {
|
206 |
fSWevent->Stop();
|
207 |
Info("Process",
|
208 |
"Events %d -> %.2gs real, %.2gs cpu (%.2g real, %.2g cpu per event)",
|
209 |
GetNEventsProcessed(), fSWevent->RealTime(), fSWevent->CpuTime(),
|
210 |
fSWevent->RealTime()/nProcessed,
|
211 |
fSWevent->CpuTime()/nProcessed);
|
212 |
fSWevent->Start();
|
213 |
}
|
214 |
}
|
215 |
|
216 |
//--------------------------------------------------------------------------------------------------
|
217 |
void AnaFwkMod::SlaveBegin()
|
218 |
{
|
219 |
// Book our histogram and start the stop watches.
|
220 |
|
221 |
fSWtotal = new TStopwatch;
|
222 |
fSWevent = new TStopwatch;
|
223 |
|
224 |
if (!PublishObj(&fAllHeaders)) {
|
225 |
SendError(kAbortAnalysis, "SlaveBegin",
|
226 |
"Could not publish all event headers with name %s.", fAllHeaders.GetName());
|
227 |
return;
|
228 |
}
|
229 |
}
|
230 |
|
231 |
//--------------------------------------------------------------------------------------------------
|
232 |
void AnaFwkMod::SlaveTerminate()
|
233 |
{
|
234 |
// Fill event histogram and printout timing information.
|
235 |
|
236 |
RetractObj(fAllHeaders.GetName());
|
237 |
|
238 |
SaveNEventsProcessed();
|
239 |
TH1D *hDAllEvents = new TH1D("hDAllEvents","Sum of processed and skimmed events",1,-0.5,0.5);
|
240 |
hDAllEvents->Fill(0.0,fNEventsSkimmed+GetNEventsProcessed());
|
241 |
hDAllEvents->SetEntries(fNEventsSkimmed+GetNEventsProcessed());
|
242 |
AddOutput(hDAllEvents);
|
243 |
|
244 |
TH1D *hDSkippedEvents = new TH1D("hDSkippedEvents","Number of skipped events",1,-0.5,0.5);
|
245 |
hDSkippedEvents->Fill(0.0,fNEventsSkipped);
|
246 |
hDSkippedEvents->SetEntries(fNEventsSkipped);
|
247 |
AddOutput(hDSkippedEvents);
|
248 |
|
249 |
fSWtotal->Stop();
|
250 |
fSWevent->Stop();
|
251 |
|
252 |
MDB(kAnalysis, 1)
|
253 |
Info("SlaveTerminate",
|
254 |
"Events %d -> %.2gs real, %.2gs cpu (%.2gs real, %.2gs cpu per event)",
|
255 |
GetNEventsProcessed(), fSWtotal->RealTime(), fSWtotal->CpuTime(),
|
256 |
fSWtotal->RealTime()/GetNEventsProcessed(),
|
257 |
fSWtotal->CpuTime()/GetNEventsProcessed());
|
258 |
|
259 |
delete fSWtotal;
|
260 |
delete fSWevent;
|
261 |
}
|