ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/MitAna/TreeMod/src/AnaFwkMod.cc
Revision: 1.7
Committed: Thu Mar 12 18:25:35 2009 UTC (16 years, 2 months ago) by loizides
Content type: text/plain
Branch: MAIN
CVS Tags: Mit_008pre2
Changes since 1.6: +6 -7 lines
Log Message:
Read skimmed headers properly.

File Contents

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