ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/MitAna/TreeMod/src/HLTFwkMod.cc
Revision: 1.8
Committed: Mon Mar 23 22:15:15 2009 UTC (16 years, 1 month ago) by loizides
Content type: text/plain
Branch: MAIN
CVS Tags: Mit_009a, Mit_009, Mit_008
Changes since 1.7: +3 -9 lines
Log Message:
Cosmetics

File Contents

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