Hall A ROOT/C++ Analyzer (podd)
Loading...
Searching...
No Matches
THaScalerEvtHandler.cxx
Go to the documentation of this file.
1
2//
3// THaScalerEvtHandler
4//
5// Event handler for Hall A scalers.
6// R. Michaels, Sept, 2014
7//
8// This class does the following
9// For a particular set of event types (here, event type 140)
10// decode the scalers and put some variables into global variables.
11// The global variables can then appear in the Podd output tree T.
12// In addition, a tree "TS" is created by this class; it contains
13// just the scaler data by itself. Note, the "fName" is concatenated
14// with "TS" to ensure the tree is unique; further, "fName" is
15// concatenated with the name of the global variables, for uniqueness.
16// The list of global variables and how they are tied to the
17// scaler module and channels is in the scaler.map file, or could
18// be hardcoded here.
19// NOTE: if you don't have the scaler map file (e.g. db_LeftScalevt.dat)
20// there will be no variable output to the Trees.
21//
22// To use in the analyzer, your setup script needs something like this
23// gHaEvtHandlers->Add (new THaScalerEvtHandler("Left","HA scaler event type 140"));
24//
25// To enable debugging you may try this in the setup script
26//
27// THaScalerEvtHandler *lscaler = new THaScalerEvtHandler("Left","HA scaler event type 140");
28// lscaler->SetDebugFile("LeftScaler.txt");
29// gHaEvtHandlers->Add (lscaler);
30//
32
33#include "THaScalerEvtHandler.h"
34#include "Scaler3800.h"
35#include "Scaler3801.h"
36#include "Scaler1151.h"
37#include "Scaler560.h"
38#include "THaEvData.h"
39#include "TString.h"
40#include <cstdio>
41#include <cstdlib>
42#include <iostream>
43#include <string>
44#include <cctype> // isspace
45#include "THaVarList.h"
46#include "VarDef.h"
47#include "THaString.h"
48#include "Textvars.h" // Podd::vsplit
49#include "Helper.h"
50#include "TTree.h"
51
52using namespace std;
53using namespace Decoder;
55using Podd::vsplit;
56
57THaScalerEvtHandler::THaScalerEvtHandler(const char *name, const char* description)
58 : THaEvtTypeHandler(name, description)
59 , evcount(0)
60 , fNormIdx(kMaxUInt)
61 , fNormSlot(kMaxUInt)
62 , dvars(nullptr)
63 , fScalerTree(nullptr)
64{}
65
67{
68 Podd::DeleteContainer(scalerloc);
69 Podd::DeleteContainer(scalers);
70 delete [] dvars;
71 delete fScalerTree;
72}
73
75{
76 if( fScalerTree )
78 return 0;
79}
80
82{
83 if( !IsMyEvent(evdata->GetEvType()) )
84 return -1;
85
86 if( fDebugFile ) {
87 *fDebugFile << endl << "---------------------------------- " << endl << endl;
88 *fDebugFile << "\nEnter THaScalerEvtHandler for fName = " << fName << endl;
89 EvDump(evdata);
90 }
91
92 if( !fScalerTree ) {
93
94 TString sname1 = "TS";
95 TString sname2 = sname1 + fName;
96 TString sname3 = fName + " Scaler Data";
97
98 if( fDebugFile ) {
99 *fDebugFile << "\nAnalyze 1st time for fName = " << fName << endl;
100 *fDebugFile << sname2 << " " << sname3 << endl;
101 }
102
103 fScalerTree = new TTree(sname2.Data(),sname3.Data());
104 fScalerTree->SetAutoSave(200000000);
105
106 TString name = "evcount";
107 TString tinfo = name + "/D";
108 fScalerTree->Branch(name.Data(), &evcount, tinfo.Data(), 4000);
109
110 for( size_t i = 0; i < scalerloc.size(); i++) {
111 name = scalerloc[i]->name;
112 tinfo = name + "/D";
113 fScalerTree->Branch(name.Data(), &dvars[i], tinfo.Data(), 4000);
114 }
115
116 } // if (!fScalerTree)
117
118
119 // Parse the data, load local data arrays.
120
121 UInt_t ndata = evdata->GetEvLength();
122 if( ndata >= MAXTEVT ) {
123 cout << "THaScalerEvtHandler:: ERROR: Event length crazy " << endl;
124 ndata = MAXTEVT - 1;
125 }
126
127 if (fDebugFile)
128 *fDebugFile << endl << endl << "THaScalerEvtHandler:: "
129 << "Debugging event type " << dec << evdata->GetEvType() << endl << endl;
130
131 const UInt_t *p = evdata->GetRawDataBuffer();
132 const UInt_t *pstop = p+ndata;
133 Bool_t ifound = false;
134
135 while( p < pstop ) {
136 if( fDebugFile ) {
137 *fDebugFile << "p and pstop " << p << " " << pstop
138 << " " << hex << *p << " " << dec << endl;
139 }
140 Int_t nskip = 1;
141 for( size_t j = 0; j < scalers.size(); j++ ) {
142 nskip = scalers[j]->Decode(p);
143 if( fDebugFile && nskip > 1 ) {
144 *fDebugFile << "\n===== Scaler # " << j << " fName = " << fName
145 << " nskip = " << nskip << endl;
146 scalers[j]->DebugPrint(fDebugFile);
147 }
148 if( nskip > 1 ) {
149 ifound = true;
150 break;
151 }
152 }
153 p = p + nskip;
154 }
155
156 if( fDebugFile ) {
157 *fDebugFile << "Finished with decoding. " << endl;
158 *fDebugFile << " Found flag = " << ifound << endl;
159 }
160
161 // L-HRS has headers which are different from R-HRS, but both are
162 // event type 140 and come here. If you found no headers, it was
163 // the other arms event type. (The arm is fName).
164
165 if( !ifound )
166 return 0;
167
168 // The correspondence between dvars and the scaler and the channel
169 // will be driven by a scaler.map file, or could be hard-coded.
170
171 for( size_t i = 0; i < scalerloc.size(); i++ ) {
172 UInt_t idx = scalerloc[i]->index;
173 UInt_t ichan = scalerloc[i]->ichan;
174 if( fDebugFile )
175 *fDebugFile << "Debug dvars " << i << " " << idx << " " << ichan << endl;
176 if( idx < scalers.size() && ichan < MAXCHAN ) {
177 if( scalerloc[i]->ikind == ICOUNT )
178 dvars[i] = scalers[idx]->GetData(ichan);
179 if( scalerloc[i]->ikind == IRATE )
180 dvars[i] = scalers[idx]->GetRate(ichan);
181 if( fDebugFile )
182 *fDebugFile << " dvars " << scalerloc[i]->ikind
183 << " " << dvars[i] << endl;
184 } else {
185 cout << "THaScalerEvtHandler:: ERROR:: incorrect index " << i
186 << " " << idx << " " << ichan << endl;
187 }
188 }
189
190 evcount += 1.0;
191
192 for( auto* s: scalers )
193 s->Clear();
194
195 if( fDebugFile )
196 *fDebugFile << "scaler tree ptr " << fScalerTree << endl;
197
198 if( fScalerTree )
199 fScalerTree->Fill();
200
201 return 1;
202}
203
204// Helper functions for Init()
205void THaScalerEvtHandler::ParseVariable( const vector<string>& dbline )
206{
207 string sdesc;
208 for( size_t j = 5; j < dbline.size(); j++ ) {
209 sdesc += " ";
210 sdesc += dbline[j];
211 }
212 Int_t islot = atoi(dbline[1].c_str());
213 Int_t ichan = atoi(dbline[2].c_str());
214 Int_t ikind = atoi(dbline[3].c_str());
215 if( fDebugFile )
216 *fDebugFile << "add var " << dbline[1] << " desc = " << sdesc
217 << " islot= " << islot << " " << ichan
218 << " " << ikind << endl;
219 TString tsname(dbline[4].c_str());
220 TString tsdesc(sdesc.c_str());
221 AddVars(tsname, tsdesc, islot, ichan, ikind);
222}
223
224void THaScalerEvtHandler::ParseMap( const char* cbuf, const vector<string>& dbline )
225{
226 using ssiz_t = decltype(scalers)::size_type;
227 Int_t imodel = 0;
228 UInt_t icrate = 0, islot = 0, inorm = 0, header = 0, mask = 0;
229 char cdum[21];
230 sscanf(cbuf, "%20s %d %u %u %x %x %u \n",
231 cdum, &imodel, &icrate, &islot, &header, &mask, &inorm);
232 if( fNormSlot != inorm )
233 cout << "THaScalerEvtHandler:: WARN: contradictory norm slot " << inorm << endl;
234 fNormSlot = inorm; // slot number used for normalization. This variable is not used but is checked.
235 UInt_t clkchan = kMaxUInt;
236 Double_t clkfreq = 1;
237 if( dbline.size() > 8 ) {
238 clkchan = atoi(dbline[7].c_str());
239 clkfreq = static_cast<Double_t>( atoi(dbline[8].c_str()) );
240 }
241 if( fDebugFile ) {
242 *fDebugFile << "map line " << dec << imodel << " " << icrate
243 << " " << islot << endl;
244 *fDebugFile << " header 0x" << hex << header << " 0x" << mask
245 << dec << " " << inorm << " " << clkchan
246 << " " << clkfreq << endl;
247 }
248 switch( imodel ) {
249 case 560:
250 scalers.push_back(new Scaler560(icrate, islot));
251 break;
252 case 1151:
253 scalers.push_back(new Scaler1151(icrate, islot));
254 break;
255 case 3800:
256 scalers.push_back(new Scaler3800(icrate, islot));
257 break;
258 case 3801:
259 scalers.push_back(new Scaler3801(icrate, islot));
260 break;
261 default:
262 break;
263 }
264 if( !scalers.empty() ) {
265 ssiz_t idx = scalers.size() - 1;
266 scalers[idx]->SetHeader(header, mask);
267// The normalization slot has the clock in it, so we automatically recognize it.
268// fNormIdx is the index in scaler[] and
269// fNormSlot is the slot#, checked for consistency
270 if( clkchan != kMaxUInt ) {
271 scalers[idx]->SetClock(defaultDT, clkchan, clkfreq);
272 fNormIdx = idx;
273 if( islot != fNormSlot )
274 cout << "THaScalerEvtHandler:: WARN: contradictory norm slot ! "
275 << islot << endl;
276 if( fDebugFile )
277 *fDebugFile << "Setting scaler clock ... channel = "
278 << clkchan << " ... freq = " << clkfreq
279 << " fNormIdx = " << fNormIdx
280 << " fNormSlot = " << fNormSlot
281 << " slot = " << islot << endl;
282 }
283 }
284}
285
287{
288 if( fNormIdx < scalers.size() ) {
289 UInt_t i = 0;
290 for( auto* scaler: scalers ) {
291 if( i++ == fNormIdx )
292 continue;
293 scaler->LoadNormScaler(scalers[fNormIdx]);
294 }
295 }
296}
297
299{
300 for( size_t i1 = 0; i1 < scalers.size() - 1; i1++ ) {
301 for( size_t i2 = i1 + 1; i2 < scalers.size(); i2++ ) {
302 if( scalers[i1]->GetSlot() == scalers[i2]->GetSlot() )
303 cout << "THaScalerEvtHandler:: WARN: same slot defined twice" << endl;
304 }
305 }
306}
307
309{
310 for( size_t i = 0; i < scalers.size(); i++ ) {
311 for( auto& loc: scalerloc ) {
312 if( loc->islot == scalers[i]->GetSlot() )
313 loc->index = i;
314 }
315 }
316}
317
319{
320 const int LEN = 200;
321 char cbuf[LEN];
322
323 fStatus = kOK;
324 fNormIdx = -1;
325
326 // cout << "Howdy ! We are initializing THaScalerEvtHandler !! name = "
327 // << fName << endl;
328
329 AddEvtType(140); // what events to look for
330
331// Parse the map file which defines what scalers exist and the global variables.
332
333 TString sname0 = "Scalevt";
334 TString sname = fName+sname0;
335
336 FILE *fi = Podd::OpenDBFile(sname.Data(), date);
337 if ( !fi ) {
338 cout << "Cannot find db file for " << fName << " scaler event handler" << endl;
339 return kFileError;
340 }
341
342 const char comment = '#';
343 const string svariable = "variable";
344 const string smap = "map";
345
346 while( fgets(cbuf, LEN, fi) ) {
347 if (fDebugFile) *fDebugFile << "string input "<<cbuf<<endl;
348 const vector<string> dbline = vsplit(cbuf);
349 if( dbline.empty() )
350 continue;
351 assert(!dbline.front().empty() && !isspace(dbline.front()[0])); // else bug in vsplit
352 if( dbline.front()[0] == comment )
353 continue;
354
355 if( dbline.size() > 4 && CmpNoCase(dbline.front(), svariable) == 0 ) {
356 ParseVariable(dbline);
357 }
358 else if( dbline.size() > 6 && CmpNoCase(dbline.front(), smap) == 0 ) {
359 ParseMap(cbuf, dbline);
360 }
361 }
362 // need to do LoadNormScaler after scalers created and if fNormIdx found.
364
365#ifdef HARDCODED
366 // This code is superseded by the parsing of a map file above. It's another way ...
367 if (fName == "Left") {
368 AddVars("TSbcmu1", "BCM x1 counts", 1, 4, ICOUNT);
369 AddVars("TSbcmu1r","BCM x1 rate", 1, 4, IRATE);
370 AddVars("TSbcmu3", "BCM u3 counts", 1, 5, ICOUNT);
371 AddVars("TSbcmu3r", "BCM u3 rate", 1, 5, IRATE);
372 } else {
373 AddVars("TSbcmu1", "BCM x1 counts", 0, 4, ICOUNT);
374 AddVars("TSbcmu1r","BCM x1 rate", 0, 4, IRATE);
375 AddVars("TSbcmu3", "BCM u3 counts", 0, 5, ICOUNT);
376 AddVars("TSbcmu3r", "BCM u3 rate", 0, 5, IRATE);
377 }
378#endif
379
380
381 DefVars();
382
383#ifdef HARDCODED
384 // This code is superseded by the parsing of a map file above. It's another way ...
385 if (fName == "Left") {
386 scalers.push_back(new Scaler1151(1,0));
387 scalers.push_back(new Scaler3800(1,1));
388 scalers.push_back(new Scaler3800(1,2));
389 scalers.push_back(new Scaler3800(1,3));
390 scalers[0]->SetHeader(0xabc00000, 0xffff0000);
391 scalers[1]->SetHeader(0xabc10000, 0xffff0000);
392 scalers[2]->SetHeader(0xabc20000, 0xffff0000);
393 scalers[3]->SetHeader(0xabc30000, 0xffff0000);
394 scalers[0]->LoadNormScaler(scalers[1]);
395 scalers[1]->SetClock(4, 7, 1024);
396 scalers[2]->LoadNormScaler(scalers[1]);
397 scalers[3]->LoadNormScaler(scalers[1]);
398 } else {
399 scalers.push_back(new Scaler3800(2,0));
400 scalers.push_back(new Scaler3800(2,0));
401 scalers.push_back(new Scaler1151(2,1));
402 scalers.push_back(new Scaler1151(2,2));
403 scalers[0]->SetHeader(0xceb00000, 0xffff0000);
404 scalers[1]->SetHeader(0xceb10000, 0xffff0000);
405 scalers[2]->SetHeader(0xceb20000, 0xffff0000);
406 scalers[3]->SetHeader(0xceb30000, 0xffff0000);
407 scalers[0]->SetClock(4, 7, 1024);
408 scalers[1]->LoadNormScaler(scalers[0]);
409 scalers[2]->LoadNormScaler(scalers[0]);
410 scalers[3]->LoadNormScaler(scalers[0]);
411 }
412#endif
413
414 // Verify that the slots are not defined twice
415 VerifySlots();
416
417 // Identify indices of scalers[] vector to variables.
418 SetIndices();
419
420 if(fDebugFile) {
421 *fDebugFile << "THaScalerEvtHandler:: Name of scaler bank "<<fName<<endl;
422 for (size_t i=0; i<scalers.size(); i++) {
423 *fDebugFile << "Scaler # "<<i<<endl;
424 scalers[i]->SetDebugFile(fDebugFile);
425 scalers[i]->DebugPrint(fDebugFile);
426 }
427 }
428
429 fclose(fi); // NOLINT(cert-err33-c) // read-only file, should always succeed
430 return kOK;
431}
432
433void THaScalerEvtHandler::AddVars( const TString& name, const TString& desc,
434 UInt_t iscal, UInt_t ichan, UInt_t ikind)
435{
436 // need to add fName here to make it a unique variable. (Left vs Right HRS, for example)
437 // We don't yet know the correspondence between index of scalers[] and slots.
438 // Will put that in later.
439 scalerloc.push_back(new ScalerLoc(fName + name, fName + desc, 0, iscal, ichan, ikind));
440}
441
443{
444 // called after AddVars has finished being called.
445 size_t Nvars = scalerloc.size();
446 if( Nvars == 0 )
447 return;
448 delete [] dvars;
449 dvars = new Double_t[Nvars]; // dvars is a member of this class
450 memset(dvars, 0, Nvars * sizeof(Double_t));
451 if( gHaVars ) {
452 if( fDebugFile )
453 *fDebugFile << "THaScalerEVtHandler:: Have gHaVars " << gHaVars << endl;
454 } else {
455 cout << "No gHaVars ?! Well, that's a problem !!" << endl;
456 return;
457 }
458 if( fDebugFile )
459 *fDebugFile << "THaScalerEvtHandler:: scalerloc size " << scalerloc.size() << endl;
460 const Int_t* count = nullptr;
461 for( size_t i = 0; i < scalerloc.size(); i++ ) {
463 scalerloc[i]->description.Data(),
464 &dvars[i], kDouble, count);
465 }
466}
467
int Int_t
unsigned int UInt_t
bool Bool_t
const UInt_t kMaxUInt
double Double_t
winID h TVirtualViewer3D TVirtualGLPainter p
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t mask
char name[80]
R__EXTERN class THaVarList * gHaVars
Definition THaGlobals.h:11
static const UInt_t IRATE
static const UInt_t defaultDT
static const UInt_t ICOUNT
static const UInt_t MAXCHAN
static const UInt_t MAXTEVT
UInt_t GetEvType() const
Definition THaEvData.h:53
const UInt_t * GetRawDataBuffer() const
Definition THaEvData.h:41
UInt_t GetEvLength() const
Definition THaEvData.h:54
virtual void AddEvtType(UInt_t evtype)
virtual Bool_t IsMyEvent(UInt_t type) const
virtual void EvDump(THaEvData *evdata) const
std::ofstream * fDebugFile
std::vector< Decoder::GenScaler * > scalers
void AddVars(const TString &name, const TString &desc, UInt_t iscal, UInt_t ichan, UInt_t ikind)
virtual Int_t Analyze(THaEvData *evdata)
THaScalerEvtHandler(const char *name, const char *description)
void ParseMap(const char *cbuf, const std::vector< std::string > &dbline)
virtual Int_t End(THaRunBase *r=nullptr)
void ParseVariable(const std::vector< std::string > &dbline)
std::vector< ScalerLoc * > scalerloc
virtual THaVar * DefineByType(const char *name, const char *desc, const void *loc, VarType type, const Int_t *count, const char *errloc="DefineByType")
TString fName
virtual void Clear(Option_t *="")
const char * Data() const
virtual Int_t Fill()
virtual void SetAutoSave(Long64_t autos=-300000000)
Int_t Write(const char *name=nullptr, Int_t option=0, Int_t bufsize=0) const override
virtual Int_t Branch(const char *folder, Int_t bufsize=32000, Int_t splitlevel=99)
int CmpNoCase(const string &r, const string &s)
Definition THaString.cxx:19
STL namespace.
ClassImp(TPyArg)