ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/MitAna/TreeMod/src/HLTFwkMod.cc
Revision: 1.13
Committed: Mon Oct 26 11:04:56 2009 UTC (15 years, 6 months ago) by loizides
Content type: text/plain
Branch: MAIN
CVS Tags: Mit_012b, Mit_012a, Mit_012
Changes since 1.12: +3 -2 lines
Log Message:
Use GetEntries

File Contents

# User Rev Content
1 loizides 1.13 // $Id: HLTFwkMod.cc,v 1.12 2009/08/11 15:24:39 loizides Exp $
2 loizides 1.1
3     #include "MitAna/TreeMod/interface/HLTFwkMod.h"
4     #include <TFile.h>
5     #include <TTree.h>
6 loizides 1.6 #include "MitAna/DataUtil/interface/Debug.h"
7 loizides 1.1 #include "MitAna/DataTree/interface/Names.h"
8     #include "MitAna/DataTree/interface/TriggerName.h"
9 loizides 1.10 #include "MitAna/DataTree/interface/TriggerTable.h"
10     #include "MitAna/DataTree/interface/TriggerObjectBaseCol.h"
11     #include "MitAna/DataTree/interface/TriggerObjectRelCol.h"
12 loizides 1.9 #include "MitAna/DataTree/interface/TriggerObjectCol.h"
13 loizides 1.10 #include "MitAna/DataTree/interface/TriggerObjectsTable.h"
14 loizides 1.1
15     using namespace mithep;
16    
17     ClassImp(mithep::HLTFwkMod)
18    
19    
20     //--------------------------------------------------------------------------------------------------
21     HLTFwkMod::HLTFwkMod(const char *name, const char *title) :
22     BaseMod(name,title),
23     fHLTTreeName(Names::gkHltTreeName),
24     fHLTTabName(Names::gkHltTableBrn),
25     fHLTLabName(Names::gkHltLabelBrn),
26     fObjsName(Names::gkHltObjBrn),
27     fRelsName(Form("%sRelation",fObjsName.Data())),
28     fHLTTabNamePub(Form("%sFwk",fHLTTabName.Data())),
29     fHLTLabNamePub(Form("%sFwk",fHLTLabName.Data())),
30     fObjsNamePub(Form("%sFwk",fObjsName.Data())),
31     fNMaxTriggers(256),
32     fObjs(0),
33     fRels(0),
34     fReload(0),
35     fHLTTree(0),
36     fHLTTab(0),
37     fHLTLab(0),
38     fCurEnt(-2),
39     fTriggers(new TriggerTable(fNMaxTriggers)),
40     fLabels(new TriggerTable(fNMaxTriggers*16)),
41     fTrigObjs(new TriggerObjectsTable(fTriggers,fNMaxTriggers))
42     {
43     // Constructor.
44    
45     fTriggers->SetName(fHLTTabNamePub);
46     fTriggers->SetOwner();
47     fLabels->SetName(fHLTLabNamePub);
48     fLabels->SetOwner();
49     fTrigObjs->SetName(fObjsNamePub);
50     fTrigObjs->SetOwner();
51     }
52    
53     //--------------------------------------------------------------------------------------------------
54     HLTFwkMod::~HLTFwkMod()
55     {
56     // Destructor.
57    
58     fReload = 0;
59     fHLTTree = 0;
60     fHLTTab = 0;
61     fHLTLab = 0;
62     fCurEnt = -2;
63     delete fTriggers;
64     fTriggers = 0;
65     delete fLabels;
66     fLabels = 0;
67     delete fTrigObjs;
68     fTrigObjs = 0;
69     }
70    
71     //--------------------------------------------------------------------------------------------------
72     void HLTFwkMod::BeginRun()
73     {
74     // Get HLT tree and set branches if new file was opened. Read next entry in HLT key
75     // depending on entry in RunInfo.
76    
77     if (fReload) {
78    
79     // reset to be (re-)loaded variables
80     fReload = 0;
81     fHLTTree = 0;
82     fHLTTab = 0;
83     fHLTLab = 0;
84     fCurEnt = -2;
85    
86     // get current file
87     TFile *file = GetCurrentFile();
88     if (!file)
89     return;
90    
91     // get HLT tree
92     fHLTTree = dynamic_cast<TTree*>(file->Get(fHLTTreeName));
93     if (!fHLTTree) {
94     SendError(kAbortAnalysis, "BeginRun",
95 loizides 1.6 "Cannot find HLT tree with name %s.", fHLTTreeName.Data());
96 loizides 1.1 }
97    
98     // get HLT trigger name branch
99     if (fHLTTree->GetBranch(fHLTTabName)) {
100     fHLTTree->SetBranchAddress(fHLTTabName, &fHLTTab);
101     } else {
102     SendError(kAbortAnalysis, "BeginRun",
103 loizides 1.6 "Cannot find HLT tree branch with name %s.", fHLTTabName.Data());
104 loizides 1.1 }
105    
106     // get HLT module labels branch
107     if (fHLTTree->GetBranch(fHLTLabName)) {
108     fHLTTree->SetBranchAddress(fHLTLabName, &fHLTLab);
109     } else {
110     SendError(kAbortAnalysis, "BeginRun",
111 loizides 1.6 "Cannot find HLT tree branch with name %s.", fHLTLabName.Data());
112 loizides 1.1 }
113     }
114    
115     // get entry for HLT info
116     const RunInfo *runinfo = GetRunInfo();
117     if (!runinfo) {
118     SendError(kAbortAnalysis, "BeginRun",
119 loizides 1.6 "Cannot obtain run info object from selector.");
120 loizides 1.1 return;
121     }
122    
123     // load trigger table
124     if (runinfo->HltEntry()!=fCurEnt) {
125 loizides 1.7 MDB(kAnalysis, 1)
126 loizides 1.6 Info("BeginRun", "Loading trigger table for run %ld", runinfo->RunNum());
127    
128 loizides 1.1 fCurEnt = runinfo->HltEntry();
129     Bool_t load = LoadTriggerTable();
130     if (!load) {
131     SendError(kAbortAnalysis, "BeginRun",
132 loizides 1.6 "Cannot load trigger table for next entry (%ld).", fCurEnt);
133 loizides 1.1 return;
134     }
135 loizides 1.6
136 loizides 1.7 MDB(kAnalysis, 2) {
137 loizides 1.6 Info("BeginRun", "Printing tables for run %ld", runinfo->RunNum());
138     cout << " --- Trigger table ---" << endl;
139     fTriggers->Print();
140     cout << " --- Module lables ---" << endl;
141     fLabels->Print();
142     }
143 loizides 1.1 }
144     }
145    
146     //--------------------------------------------------------------------------------------------------
147     Bool_t HLTFwkMod::LoadTriggerTable()
148     {
149     // Load next trigger (and module) table from HLT tree.
150    
151     if (fCurEnt<0)
152     return kFALSE;
153    
154     // delete old tables
155     fTriggers->Delete();
156     fLabels->Delete();
157    
158     // load next event in HLT tree
159     fHLTTab = 0;
160     fHLTLab = 0;
161     Int_t ret = fHLTTree->GetEvent(fCurEnt);
162     if (ret<0 || fHLTTab==0 || fHLTTab==0 ) {
163 loizides 1.8 SendError(kAbortAnalysis, "LoadTriggerTable",
164     "Could not get trigger data for next entry (%ld).", fCurEnt);
165 loizides 1.1 return kFALSE;
166     }
167    
168     // check size of trigger table
169 loizides 1.5 if (fHLTTab->size()>fNMaxTriggers) {
170 loizides 1.1 SendError(kAbortAnalysis, "LoadTriggerTable",
171 loizides 1.6 "Size of trigger table (%ld) larger than maximum (%ld).",
172 loizides 1.7 fHLTTab->size(), fNMaxTriggers);
173 loizides 1.1 return kFALSE;
174     }
175    
176     // add trigger names
177 loizides 1.5 for (UInt_t i=0; i<fHLTTab->size(); ++i) {
178     TriggerName *tname = new TriggerName(fHLTTab->at(i),i);
179 loizides 1.1 fTriggers->Add(tname);
180     }
181    
182     // add module labels
183 loizides 1.5 for (UInt_t i=0; i<fHLTLab->size(); ++i) {
184     TriggerName *tname = new TriggerName(fHLTLab->at(i),i);
185 loizides 1.1 fLabels->Add(tname);
186     }
187    
188     return kTRUE;
189     }
190    
191     //--------------------------------------------------------------------------------------------------
192     Bool_t HLTFwkMod::Notify()
193     {
194     // Save notification of a new file, which will trigger the reloading of the tables and bitmasks.
195    
196     fReload = kTRUE;
197     return kTRUE;
198     }
199    
200     //--------------------------------------------------------------------------------------------------
201     void HLTFwkMod::Process()
202     {
203     // Read trigger objects and relation branch and fill our object table.
204    
205     fTrigObjs->Delete();
206    
207 loizides 1.2 LoadBranch(fObjsName);
208     LoadBranch(fRelsName);
209 loizides 1.1
210 loizides 1.13 const UInt_t n = fRels->GetEntries();
211     for (UInt_t i=0; i<n; ++i) {
212 loizides 1.1 const TriggerObjectRel *rel = fRels->At(i);
213 loizides 1.3 if (!rel) continue;
214 loizides 1.1
215     const TriggerObjectBase *ob = fObjs->At(rel->ObjInd());
216 loizides 1.3 if (!ob) continue;
217 loizides 1.1
218     TriggerObject *obj = new TriggerObject(rel->TrgId(), rel->Type(), ob->Id(),
219     ob->Pt(), ob->Eta(), ob->Phi(), ob->Mass());
220    
221 loizides 1.5 obj->SetTrigName(fHLTTab->at(rel->TrgId()).c_str());
222     obj->SetModuleName(fHLTLab->at(rel->ModInd()).c_str());
223     obj->SetFilterName(fHLTLab->at(rel->FilterInd()).c_str());
224 loizides 1.11 if (obj->TagInd()>=0)
225     obj->SetTagName(fHLTLab->at(obj->TagInd()).c_str());
226     else
227     obj->SetTagName("Unknown");
228    
229 loizides 1.1 fTrigObjs->Add(obj);
230     }
231     }
232    
233     //--------------------------------------------------------------------------------------------------
234     void HLTFwkMod::SlaveBegin()
235     {
236     // Request branches for trigger objects and relation, and publish our tables.
237    
238 loizides 1.12 if (fObjsName != Names::gkHltObjBrn)
239     fRelsName = Form("%sRelation",fObjsName.Data());
240    
241 loizides 1.2 ReqBranch(fObjsName, fObjs);
242     ReqBranch(fRelsName, fRels);
243 loizides 1.1
244     if (!PublishObj(fTriggers)) {
245     SendError(kAbortAnalysis, "SlaveBegin",
246     "Could not publish HLT trigger table with name %s.", fTriggers->GetName());
247     return;
248     }
249     if (!PublishObj(fTrigObjs)) {
250     SendError(kAbortAnalysis, "SlaveBegin",
251     "Could not publish HLT trigger objects table with name %s.", fTrigObjs->GetName());
252     return;
253     }
254     if (!PublishObj(fLabels)) {
255     SendError(kAbortAnalysis, "SlaveBegin",
256     "Could not publish HLT labels with name %s.", fLabels->GetName());
257     return;
258     }
259     }
260    
261     //--------------------------------------------------------------------------------------------------
262     void HLTFwkMod::SlaveTerminate()
263     {
264     // Retract our published objects.
265    
266     RetractObj(fTriggers->GetName());
267     RetractObj(fLabels->GetName());
268     RetractObj(fTrigObjs->GetName());
269     }