CbmRoot
Loading...
Searching...
No Matches
CbmStar2019MonitorTof.cxx
Go to the documentation of this file.
1/* Copyright (C) 2018-2021 Facility for Antiproton and Ion Research in Europe, Darmstadt
2 SPDX-License-Identifier: GPL-3.0-only
3 Authors: Pierre-Alain Loizeau [committer] */
4
5// -----------------------------------------------------------------------------
6// ----- -----
7// ----- CbmStar2019MonitorTof -----
8// ----- Created 10.07.2018 by P.-A. Loizeau -----
9// ----- -----
10// -----------------------------------------------------------------------------
11
13
14#include "CbmHistManager.h"
15#include "CbmStar2019TofPar.h"
16#include "FairRootManager.h"
17#include "FairRun.h"
18#include "FairRunOnline.h"
19#include "FairRuntimeDb.h"
20#include "Rtypes.h"
21#include "TCanvas.h"
22#include "TClonesArray.h"
23#include "TF1.h"
24#include "TH1.h"
25#include "TH2.h"
26#include "THStack.h"
27#include "THttpServer.h"
28#include "TMath.h"
29#include "TPaveStats.h"
30#include "TProfile.h"
31#include "TProfile2D.h"
32#include "TROOT.h"
33#include "TString.h"
34#include "TStyle.h"
35
36#include <Logger.h>
37
38#include <algorithm>
39#include <cstdint>
40#include <ctime>
41#include <iomanip>
42#include <iostream>
43
50
56 , fbIgnoreOverlapMs(kFALSE)
58 , fuTotalMsNb(0)
59 , fuOverlapMsNb(0)
60 , fuCoreMs(0)
61 , fdMsSizeInNs(0.0)
62 , fdTsCoreSizeInNs(0.0)
63 , fuMinNbGdpb(0)
64 , fuCurrNbGdpb(0)
65 , fUnpackPar()
66 , fuNrOfGdpbs(0)
71 , fuNrOfGet4(0)
74 , fuRawDataPrintMsgNb(100000)
76 , fbPrintAllHitsEnable(kFALSE)
78 , fbPulserModeEnable(kFALSE)
79 , fbCoincMapsEnable(kFALSE)
81 , fuCurrentMs(0)
82 , fdMsIndex(0)
83 , fuGdpbId(0)
84 , fuGdpbNr(0)
85 , fuGet4Id(0)
86 , fuGet4Nr(0)
87 , fiEquipmentId(0)
88 , fviMsgCounter(11, 0)
89 , // length of enum MessageTypes initialized with 0
110 , dMinDt(-1. * (kuNbBinsDt * gdpbv100::kdBinSize / 2.) - gdpbv100::kdBinSize / 2.)
111 , dMaxDt(1. * (kuNbBinsDt * gdpbv100::kdBinSize / 2.) + gdpbv100::kdBinSize / 2.)
112 , fuNbFeePlot(2)
114 , fdStartTime(-1.)
115 , fdStartTimeLong(-1.)
116 , fdStartTimeMsSz(-1.)
117 , fuHistoryHistoSize(1800)
120 , fdFitZoomWidthPs(0.0)
121 , fcMsSizeAll(NULL)
122 , fvhMsSzPerLink(12, NULL)
123 , fvhMsSzTimePerLink(12, NULL)
124 , fhMessType(NULL)
125 , fhSysMessType(NULL)
126 , fhGet4MessType(NULL)
127 , fhGet4ChanScm(NULL)
128 , fhGet4ChanErrors(NULL)
129 , fhGet4EpochFlags(NULL)
130 , fhGdpbMessType(NULL)
131 , fhGdpbSysMessType(NULL)
133 , fhGdpbEpochFlags(NULL)
134 , fhGdpbEpochSyncEvo(NULL)
135 , fhGdpbEpochMissEvo(NULL)
139 , fhScmScalerCounters(NULL)
141 , fhScmSeuCounters(NULL)
142 , fhScmSeuCountersEvo(NULL)
143 , fhPatternMissmatch(NULL)
144 , fhPatternEnable(NULL)
145 , fhPatternResync(NULL)
149 , fvhRawFt_gDPB()
164 , fvhModRate()
174 , fhTimeMeanPulser(NULL)
175 , fhTimeRmsPulser(NULL)
177 , fhTimeResFitPuls(NULL)
180 , fvuPadiToGet4()
181 , fvuGet4ToPadi()
185{
186}
187
189
191{
192 LOG(info) << "Initializing Get4 monitor";
193
194 FairRootManager* ioman = FairRootManager::Instance();
195 if (ioman == NULL) { LOG(fatal) << "No FairRootManager instance"; } // if( ioman == NULL )
196
197 return kTRUE;
198}
199
201{
202 LOG(info) << "Setting parameter containers for " << GetName();
203 fUnpackPar = (CbmStar2019TofPar*) (FairRun::Instance()->GetRuntimeDb()->getContainer("CbmStar2019TofPar"));
204}
205
207{
208 LOG(info) << "Init parameter containers for " << GetName();
209 Bool_t initOK = ReInitContainers();
210
212
217 for (UInt_t i = 0; i < fuNrOfGdpbs; ++i) {
218 for (UInt_t j = 0; j < fuNrOfGet4PerGdpb; ++j) {
222 } // for( UInt_t j = 0; j < fuNrOfGet4PerGdpb; ++j )
223 } // for( UInt_t i = 0; i < fuNrOfGdpbs; ++i )
224
225 return initOK;
226}
227
229{
230 LOG(info) << "ReInit parameter containers for " << GetName();
231
232 fuNrOfGdpbs = fUnpackPar->GetNrOfGdpbs();
233 LOG(info) << "Nr. of Tof GDPBs: " << fuNrOfGdpbs;
235
236 fuNrOfFeePerGdpb = fUnpackPar->GetNrOfFeesPerGdpb();
237 LOG(info) << "Nr. of FEEs per Tof GDPB: " << fuNrOfFeePerGdpb;
238
239 fuNrOfGet4PerFee = fUnpackPar->GetNrOfGet4PerFee();
240 LOG(info) << "Nr. of GET4 per Tof FEE: " << fuNrOfGet4PerFee;
241
242 fuNrOfChannelsPerGet4 = fUnpackPar->GetNrOfChannelsPerGet4();
243 LOG(info) << "Nr. of channels per GET4: " << fuNrOfChannelsPerGet4;
244
246 LOG(info) << "Nr. of channels per FEE: " << fuNrOfChannelsPerFee;
247
249 LOG(info) << "Nr. of GET4s: " << fuNrOfGet4;
250
252 LOG(info) << "Nr. of GET4s per GDPB: " << fuNrOfGet4PerGdpb;
253
255 LOG(info) << "Nr. of channels per GDPB: " << fuNrOfChannelsPerGdpb;
256
257 fGdpbIdIndexMap.clear();
258 for (UInt_t i = 0; i < fuNrOfGdpbs; ++i) {
259 fGdpbIdIndexMap[fUnpackPar->GetGdpbId(i)] = i;
260 LOG(info) << "GDPB Id of TOF " << i << " : " << std::hex << fUnpackPar->GetGdpbId(i) << std::dec;
261 } // for( UInt_t i = 0; i < fuNrOfGdpbs; ++i )
262
263 fuNrOfGbtx = fUnpackPar->GetNrOfGbtx();
264 LOG(info) << "Nr. of GBTx: " << fuNrOfGbtx;
265
266 fuNrOfModules = fUnpackPar->GetNrOfModules();
267 LOG(info) << "Nr. of GBTx: " << fuNrOfModules;
268
269 fviRpcType.resize(fuNrOfGbtx);
270 fviModuleId.resize(fuNrOfGbtx);
271 fviNrOfRpc.resize(fuNrOfGbtx);
272 fviRpcSide.resize(fuNrOfGbtx);
273 for (UInt_t uGbtx = 0; uGbtx < fuNrOfGbtx; ++uGbtx) {
274 fviNrOfRpc[uGbtx] = fUnpackPar->GetNrOfRpc(uGbtx);
275 fviRpcType[uGbtx] = fUnpackPar->GetRpcType(uGbtx);
276 fviRpcSide[uGbtx] = fUnpackPar->GetRpcSide(uGbtx);
277 fviModuleId[uGbtx] = fUnpackPar->GetModuleId(uGbtx);
278 } // for( UInt_t uGbtx = 0; uGbtx < uNrOfGbtx; ++uGbtx)
279
280 LOG(info) << "Nr. of RPCs per GBTx: ";
281 std::stringstream ss;
282 for (UInt_t uGbtx = 0; uGbtx < fuNrOfGbtx; ++uGbtx)
283 ss << Form(" %2d", fviNrOfRpc[uGbtx]);
284 LOG(info) << ss.str();
285
286 LOG(info) << "RPC type per GBTx: ";
287 std::stringstream ss;
288 for (UInt_t uGbtx = 0; uGbtx < fuNrOfGbtx; ++uGbtx)
289 ss << Form(" %2d", fviRpcType[uGbtx]);
290 LOG(info) << ss.str();
291
292 LOG(info) << "RPC side per GBTx: ";
293 std::stringstream ss;
294 for (UInt_t uGbtx = 0; uGbtx < fuNrOfGbtx; ++uGbtx)
295 ss << Form(" %2d", fviRpcSide[uGbtx]);
296 LOG(info) << ss.str();
297
298 LOG(info) << "Module ID per GBTx: ";
299 std::stringstream ss;
300 for (UInt_t uGbtx = 0; uGbtx < fuNrOfGbtx; ++uGbtx)
301 ss << Form(" %2d", fviModuleId[uGbtx]);
302 LOG(info) << ss.str();
303
304 fuTotalMsNb = fUnpackPar->GetNbMsTot();
305 fuOverlapMsNb = fUnpackPar->GetNbMsOverlap();
307 fdMsSizeInNs = fUnpackPar->GetSizeMsInNs();
309 LOG(info) << "Timeslice parameters: " << fuTotalMsNb << " MS per link, of which " << fuOverlapMsNb
310 << " overlap MS, each MS is " << fdMsSizeInNs << " ns";
311
322 for (UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb) {
323 fvulGdpbTsMsb[uGdpb] = 0;
324 fvulGdpbTsLsb[uGdpb] = 0;
325 fvulStarTsMsb[uGdpb] = 0;
326 fvulStarTsMid[uGdpb] = 0;
327 fvulGdpbTsFullLast[uGdpb] = 0;
328 fvulStarTsFullLast[uGdpb] = 0;
329 fvuStarTokenLast[uGdpb] = 0;
330 fvuStarDaqCmdLast[uGdpb] = 0;
331 fvuStarTrigCmdLast[uGdpb] = 0;
332 } // for (Int_t iGdpb = 0; iGdpb < fuNrOfGdpbs; ++iGdpb)
333
335
337 if (kTRUE == fbPulserModeEnable) {
340 } // if( kTRUE == fbPulserModeEnable )
341
343 if (kTRUE == fbCoincMapsEnable) {
346 for (UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb) {
348 fvdCoincTsLastHit[uGdpb].resize(fuNrOfChannelsPerGdpb, 0.0);
349 } // for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
350 } // if( kTRUE == fbCoincMapsEnable )
351
355 /*
356 UInt_t uGet4topadi[32] = {
357 4, 3, 2, 1, // provided by Jochen
358 24, 23, 22, 21,
359 8, 7, 6, 5,
360 28, 27, 26, 25,
361 12, 11, 10, 9,
362 32, 31, 30, 29,
363 16, 15, 14, 13,
364 20, 19, 18, 17 };
365*/
367 UInt_t uGet4topadi[32] = {4, 3, 2, 1, // provided by Jochen
368 8, 7, 6, 5, 12, 11, 10, 9, 16, 15, 14, 13, 20, 19,
369 18, 17, 24, 23, 22, 21, 28, 27, 26, 25, 32, 31, 30, 29};
370
371 UInt_t uPaditoget4[32] = {4, 3, 2, 1, // provided by Jochen
372 12, 11, 10, 9, 20, 19, 18, 17, 28, 27, 26, 25, 32, 31,
373 30, 29, 8, 7, 6, 5, 16, 15, 14, 13, 24, 23, 22, 21};
374
375 for (UInt_t uChan = 0; uChan < fuNrOfChannelsPerFee; ++uChan) {
376 fvuPadiToGet4[uChan] = uPaditoget4[uChan] - 1;
377 fvuGet4ToPadi[uChan] = uGet4topadi[uChan] - 1;
378 } // for( UInt_t uChan = 0; uChan < fuNrOfChannelsPerFee; ++uChan )
379
380
384 UInt_t kuElinkToGet4[kuNbGet4PerGbtx] = {27, 2, 7, 3, 31, 26, 30, 1, 33, 37, 32, 13, 9, 14,
385 10, 15, 17, 21, 16, 35, 34, 38, 25, 24, 0, 6, 20, 23,
386 18, 22, 28, 4, 29, 5, 19, 36, 39, 8, 12, 11};
387 UInt_t kuGet4ToElink[kuNbGet4PerGbtx] = {24, 7, 1, 3, 31, 33, 25, 2, 37, 12, 14, 39, 38, 11,
388 13, 15, 18, 16, 28, 34, 26, 17, 29, 27, 23, 22, 5, 0,
389 30, 32, 6, 4, 10, 8, 20, 19, 35, 9, 21, 36};
390
391 for (UInt_t uLinkAsic = 0; uLinkAsic < kuNbGet4PerGbtx; ++uLinkAsic) {
392 fvuElinkToGet4[uLinkAsic] = kuElinkToGet4[uLinkAsic];
393 fvuGet4ToElink[uLinkAsic] = kuGet4ToElink[uLinkAsic];
394 } // for( UInt_t uChan = 0; uChan < fuNrOfChannelsPerFee; ++uChan )
395
396 return kTRUE;
397}
398
399
400void CbmStar2019MonitorTof::AddMsComponentToList(size_t component, UShort_t usDetectorId)
401{
403 for (UInt_t uCompIdx = 0; uCompIdx < fvMsComponentsList.size(); ++uCompIdx)
404 if (component == fvMsComponentsList[uCompIdx]) return;
405
407 fvMsComponentsList.push_back(component);
408
410 if (NULL == fvhMsSzPerLink[component]) {
411 TString sMsSzName = Form("MsSz_link_%02lu", component);
412 TString sMsSzTitle = Form("Size of MS from link %02lu; Ms Size [bytes]", component);
413 fvhMsSzPerLink[component] = new TH1F(sMsSzName.Data(), sMsSzTitle.Data(), 160000, 0., 20000.);
414
415 sMsSzName = Form("MsSzTime_link_%02lu", component);
416 sMsSzTitle = Form("Size of MS vs time for gDPB of link %02lu; Time[s] ; Ms Size [bytes]", component);
417 fvhMsSzTimePerLink[component] =
418 new TProfile(sMsSzName.Data(), sMsSzTitle.Data(), 100 * fuHistoryHistoSize, 0., 2 * fuHistoryHistoSize);
419 THttpServer* server = FairRunOnline::Instance()->GetHttpServer();
420 if (server) {
421 server->Register("/FlibRaw", fvhMsSzPerLink[component]);
422 server->Register("/FlibRaw", fvhMsSzTimePerLink[component]);
423 } // if( server )
424 if (NULL != fcMsSizeAll) {
425 fcMsSizeAll->cd(1 + component);
426 gPad->SetLogy();
427 fvhMsSzTimePerLink[component]->Draw("hist le0");
428 } // if( NULL != fcMsSizeAll )
429 LOG(info) << "Added MS size histo for component (link): " << component;
430 } // if( NULL == fvhMsSzPerLink[ component ] )
431}
432void CbmStar2019MonitorTof::SetNbMsInTs(size_t uCoreMsNb, size_t uOverlapMsNb)
433{
434 fuNbCoreMsPerTs = uCoreMsNb;
435 fuNbOverMsPerTs = uOverlapMsNb;
436
437 // UInt_t uNbMsTotal = fuNbCoreMsPerTs + fuNbOverMsPerTs;
438}
439
441{
442 LOG(info) << "create Histos for " << fuNrOfGdpbs << " gDPBs ";
443
444 THttpServer* server = FairRunOnline::Instance()->GetHttpServer();
445
446 TString name {""};
447 TString title {""};
448
449 // Full Fee time difference test
450 UInt_t uNbBinsDt = kuNbBinsDt + 1; // To account for extra bin due to shift by 1/2 bin of both ranges
451
453 Double_t dBinSzG4v2 = (6250. / 112.);
454 dMinDt = -1. * (kuNbBinsDt * dBinSzG4v2 / 2.) - dBinSzG4v2 / 2.;
455 dMaxDt = 1. * (kuNbBinsDt * dBinSzG4v2 / 2.) + dBinSzG4v2 / 2.;
456
457 /*******************************************************************/
458 name = "hMessageType";
459 title = "Nb of message for each type; Type";
460 // Test Big Data readout with plotting
461 fhMessType = new TH1I(name, title, 1 + gdpbv100::MSG_STAR_TRI_D, 0., 1 + gdpbv100::MSG_STAR_TRI_D);
462 fhMessType->GetXaxis()->SetBinLabel(1 + gdpbv100::MSG_HIT, "HIT");
463 fhMessType->GetXaxis()->SetBinLabel(1 + gdpbv100::MSG_EPOCH, "EPOCH");
464 fhMessType->GetXaxis()->SetBinLabel(1 + gdpbv100::MSG_SLOWC, "SLOWC");
465 fhMessType->GetXaxis()->SetBinLabel(1 + gdpbv100::MSG_SYST, "SYST");
466 fhMessType->GetXaxis()->SetBinLabel(1 + gdpbv100::MSG_STAR_TRI_A, "TRI_A");
467 fhMessType->GetXaxis()->SetBinLabel(1 + gdpbv100::MSG_STAR_TRI_B, "TRI_B");
468 fhMessType->GetXaxis()->SetBinLabel(1 + gdpbv100::MSG_STAR_TRI_C, "TRI_C");
469 fhMessType->GetXaxis()->SetBinLabel(1 + gdpbv100::MSG_STAR_TRI_D, "TRI_D");
470
471 /*******************************************************************/
472 name = "hSysMessType";
473 title = "Nb of system message for each type; System Type";
474 fhSysMessType = new TH1I(name, title, 1 + gdpbv100::SYS_PATTERN, 0., 1 + gdpbv100::SYS_PATTERN);
475 fhSysMessType->GetXaxis()->SetBinLabel(1 + gdpbv100::SYS_GET4_ERROR, "GET4 ERROR");
476 fhSysMessType->GetXaxis()->SetBinLabel(1 + gdpbv100::SYS_GDPB_UNKWN, "UNKW GET4 MSG");
477 fhSysMessType->GetXaxis()->SetBinLabel(1 + gdpbv100::SYS_GET4_SYNC_MISS, "SYS_GET4_SYNC_MISS");
478 fhSysMessType->GetXaxis()->SetBinLabel(1 + gdpbv100::SYS_PATTERN, "SYS_PATTERN");
479
480 /*******************************************************************/
481 name = "hGet4MessType";
482 title = "Nb of message for each type per GET4; GET4 chip # ; Type";
483 fhGet4MessType = new TH2I(name, title, fuNrOfGet4, 0., fuNrOfGet4, 4, 0., 4.);
484 fhGet4MessType->GetYaxis()->SetBinLabel(1, "DATA 32b");
485 fhGet4MessType->GetYaxis()->SetBinLabel(2, "EPOCH");
486 fhGet4MessType->GetYaxis()->SetBinLabel(3, "S.C. M");
487 fhGet4MessType->GetYaxis()->SetBinLabel(4, "ERROR");
488 // fhGet4MessType->GetYaxis()->SetBinLabel( 5, "DATA 24b");
489 // fhGet4MessType->GetYaxis()->SetBinLabel( 6, "STAR Trigger");
490
491 /*******************************************************************/
492 name = "hGet4ChanScm";
493 title = "SC messages per GET4 channel; GET4 channel # ; SC type";
495 new TH2I(name, title, 2 * fuNrOfGet4 * fuNrOfChannelsPerGet4, 0., fuNrOfGet4 * fuNrOfChannelsPerGet4, 5, 0., 5.);
496 fhGet4ChanScm->GetYaxis()->SetBinLabel(1, "Hit Scal");
497 fhGet4ChanScm->GetYaxis()->SetBinLabel(2, "Deadtime");
498 fhGet4ChanScm->GetYaxis()->SetBinLabel(3, "SPI");
499 fhGet4ChanScm->GetYaxis()->SetBinLabel(4, "SEU Scal");
500 fhGet4ChanScm->GetYaxis()->SetBinLabel(5, "START");
501
502 /*******************************************************************/
503 name = "hGet4ChanErrors";
504 title = "Error messages per GET4 channel; GET4 channel # ; Error";
506 new TH2I(name, title, fuNrOfGet4 * fuNrOfChannelsPerGet4, 0., fuNrOfGet4 * fuNrOfChannelsPerGet4, 21, 0., 21.);
507 fhGet4ChanErrors->GetYaxis()->SetBinLabel(1, "0x00: Readout Init ");
508 fhGet4ChanErrors->GetYaxis()->SetBinLabel(2, "0x01: Sync ");
509 fhGet4ChanErrors->GetYaxis()->SetBinLabel(3, "0x02: Epoch count sync");
510 fhGet4ChanErrors->GetYaxis()->SetBinLabel(4, "0x03: Epoch ");
511 fhGet4ChanErrors->GetYaxis()->SetBinLabel(5, "0x04: FIFO Write ");
512 fhGet4ChanErrors->GetYaxis()->SetBinLabel(6, "0x05: Lost event ");
513 fhGet4ChanErrors->GetYaxis()->SetBinLabel(7, "0x06: Channel state ");
514 fhGet4ChanErrors->GetYaxis()->SetBinLabel(8, "0x07: Token Ring state");
515 fhGet4ChanErrors->GetYaxis()->SetBinLabel(9, "0x08: Token ");
516 fhGet4ChanErrors->GetYaxis()->SetBinLabel(10, "0x09: Error Readout ");
517 fhGet4ChanErrors->GetYaxis()->SetBinLabel(11, "0x0a: SPI ");
518 fhGet4ChanErrors->GetYaxis()->SetBinLabel(12, "0x0b: DLL Lock error "); // <- From GET4 v1.2
519 fhGet4ChanErrors->GetYaxis()->SetBinLabel(13, "0x0c: DLL Reset invoc."); // <- From GET4 v1.2
520 fhGet4ChanErrors->GetYaxis()->SetBinLabel(14, "0x11: Overwrite ");
521 fhGet4ChanErrors->GetYaxis()->SetBinLabel(15, "0x12: ToT out of range");
522 fhGet4ChanErrors->GetYaxis()->SetBinLabel(16, "0x13: Event Discarded ");
523 fhGet4ChanErrors->GetYaxis()->SetBinLabel(17, "0x14: Add. Rising edge"); // <- From GET4 v1.3
524 fhGet4ChanErrors->GetYaxis()->SetBinLabel(18, "0x15: Unpaired Falling"); // <- From GET4 v1.3
525 fhGet4ChanErrors->GetYaxis()->SetBinLabel(19, "0x16: Sequence error "); // <- From GET4 v1.3
526 fhGet4ChanErrors->GetYaxis()->SetBinLabel(20, "0x7f: Unknown ");
527 fhGet4ChanErrors->GetYaxis()->SetBinLabel(21, "Corrupt/unsuprtd error");
528
529 /*******************************************************************/
530 name = "hGet4EpochFlags";
531 title = "Epoch flags per GET4; GET4 chip # ; Type";
532 fhGet4EpochFlags = new TH2I(name, title, fuNrOfGet4, 0., fuNrOfGet4, 4, 0., 4.);
533 fhGet4EpochFlags->GetYaxis()->SetBinLabel(1, "SYNC");
534 fhGet4EpochFlags->GetYaxis()->SetBinLabel(2, "Ep LOSS");
535 fhGet4EpochFlags->GetYaxis()->SetBinLabel(3, "Da LOSS");
536 fhGet4EpochFlags->GetYaxis()->SetBinLabel(4, "MISSMAT");
537
538 /*******************************************************************/
539 name = "hGdpbMessageType";
540 title = "Nb of message for each type per Gdpb; Type; Gdpb Idx []";
541 // Test Big Data readout with plotting
543 -0.5, fuNrOfGdpbs - 0.5);
544 fhGdpbMessType->GetXaxis()->SetBinLabel(1 + gdpbv100::MSG_HIT, "HIT");
545 fhGdpbMessType->GetXaxis()->SetBinLabel(1 + gdpbv100::MSG_EPOCH, "EPOCH");
546 fhGdpbMessType->GetXaxis()->SetBinLabel(1 + gdpbv100::MSG_SLOWC, "SLOWC");
547 fhGdpbMessType->GetXaxis()->SetBinLabel(1 + gdpbv100::MSG_SYST, "SYST");
548 fhGdpbMessType->GetXaxis()->SetBinLabel(1 + gdpbv100::MSG_STAR_TRI_A, "TRI_A");
549 fhGdpbMessType->GetXaxis()->SetBinLabel(1 + gdpbv100::MSG_STAR_TRI_B, "TRI_B");
550 fhGdpbMessType->GetXaxis()->SetBinLabel(1 + gdpbv100::MSG_STAR_TRI_C, "TRI_C");
551 fhGdpbMessType->GetXaxis()->SetBinLabel(1 + gdpbv100::MSG_STAR_TRI_D, "TRI_D");
552
553 /*******************************************************************/
554 name = "hGdpbSysMessType";
555 title = "Nb of system message for each type per Gdpb; System Type; Gdpb Idx []";
556 fhGdpbSysMessType = new TH2I(name, title, 1 + gdpbv100::SYS_PATTERN, 0., 1 + gdpbv100::SYS_PATTERN, fuNrOfGdpbs, -0.5,
557 fuNrOfGdpbs - 0.5);
558 fhGdpbSysMessType->GetXaxis()->SetBinLabel(1 + gdpbv100::SYS_GET4_ERROR, "GET4 ERROR");
559 fhGdpbSysMessType->GetXaxis()->SetBinLabel(1 + gdpbv100::SYS_GDPB_UNKWN, "UNKW GET4 MSG");
560 fhGdpbSysMessType->GetXaxis()->SetBinLabel(1 + gdpbv100::SYS_GET4_SYNC_MISS, "SYS_GET4_SYNC_MISS");
561 fhGdpbSysMessType->GetXaxis()->SetBinLabel(1 + gdpbv100::SYS_PATTERN, "SYS_PATTERN");
562
563 /*******************************************************************/
564 name = "hGdpbSysMessPattType";
565 title = "Nb of pattern message for each type per Gdpb; Pattern Type; Gdpb Idx []";
567 -0.5, fuNrOfGdpbs - 0.5);
568 fhGdpbSysMessPattType->GetXaxis()->SetBinLabel(1 + gdpbv100::PATT_MISSMATCH, "PATT_MISSMATCH");
569 fhGdpbSysMessPattType->GetXaxis()->SetBinLabel(1 + gdpbv100::PATT_ENABLE, "PATT_ENABLE");
570 fhGdpbSysMessPattType->GetXaxis()->SetBinLabel(1 + gdpbv100::PATT_RESYNC, "PATT_RESYNC");
571
572 /*******************************************************************/
573 name = "hGdpbEpochFlags";
574 title = "Epoch flags per gDPB; gDPB # ; Type";
575 fhGdpbEpochFlags = new TH2I(name, title, fuNrOfGdpbs, 0., fuNrOfGdpbs, 4, 0., 4.);
576 fhGdpbEpochFlags->GetYaxis()->SetBinLabel(1, "SYNC");
577 fhGdpbEpochFlags->GetYaxis()->SetBinLabel(2, "Ep LOSS");
578 fhGdpbEpochFlags->GetYaxis()->SetBinLabel(3, "Da LOSS");
579 fhGdpbEpochFlags->GetYaxis()->SetBinLabel(4, "MISSMAT");
580
581 /*******************************************************************/
582 name = Form("hGdpbEpochSyncEvo");
583 title = Form("Epoch SYNC per second and gDPB; Time[s]; gDPB #; SYNC Nb");
585 new TH2D(name.Data(), title.Data(), fuHistoryHistoSize, 0, fuHistoryHistoSize, fuNrOfGdpbs, 0., fuNrOfGdpbs);
586
587 /*******************************************************************/
588 name = Form("hGdpbEpochMissEvo");
589 title = Form("Epoch Missmatch per second and gDPB; Time[s]; gDPB #; Missmatch Nb");
591 new TH2D(name.Data(), title.Data(), fuHistoryHistoSize, 0, fuHistoryHistoSize, fuNrOfGdpbs, 0., fuNrOfGdpbs);
592
593 /*******************************************************************/
594 // Slow control messages analysis
595 name = "hScmScalerCounters";
596 title = "Content of Scaler counter SC messages; Scaler counter [hit]; Channel []";
597 fhScmScalerCounters = new TH2I(name, title, fuNrOfGet4 * fuNrOfChannelsPerGet4 * 2, 0.,
598 fuNrOfGet4 * fuNrOfChannelsPerGet4, 8192, 0., 8192.);
599
600 name = "hScmDeadtimeCounters";
601 title = "Content of Deadtime counter SC messages; Deadtime [Clk Cycles]; "
602 "Channel []";
603 fhScmDeadtimeCounters = new TH2I(name, title, fuNrOfGet4 * fuNrOfChannelsPerGet4 * 2, 0.,
604 fuNrOfGet4 * fuNrOfChannelsPerGet4, 8192, 0., 8192.);
605
606 name = "hScmSeuCounters";
607 title = "Content of SEU counter SC messages; SEU []; Channel []";
608 fhScmSeuCounters = new TH2I(name, title, fuNrOfGet4 * fuNrOfChannelsPerGet4 * 2, 0.,
609 fuNrOfGet4 * fuNrOfChannelsPerGet4, 8192, 0., 8192.);
610
611 name = "hScmSeuCountersEvo";
612 title = "SEU counter rate from SC messages; Time in Run [s]; Channel []; SEU []";
613 fhScmSeuCountersEvo = new TH2I(name, title, fuNrOfGet4 * fuNrOfChannelsPerGet4 * 2, 0.,
615
616 /*******************************************************************/
617 name = "hPatternMissmatch";
618 title = "Missmatch pattern integral per Gdpb; ASIC Pattern []; Gdpb Idx []";
620 new TH2I(name, title, fuNrOfGet4PerGdpb, 0., fuNrOfGet4PerGdpb, fuNrOfGdpbs, -0.5, fuNrOfGdpbs - 0.5);
621 name = "hPatternEnable";
622 title = "Enable pattern integral per Gdpb; ASIC Pattern []; Gdpb Idx []";
624 new TH2I(name, title, fuNrOfGet4PerGdpb, 0., fuNrOfGet4PerGdpb, fuNrOfGdpbs, -0.5, fuNrOfGdpbs - 0.5);
625 name = "hPatternResync";
626 title = "Resync pattern integral per Gdpb; ASIC Pattern []; Gdpb Idx []";
628 new TH2I(name, title, fuNrOfGet4PerGdpb, 0., fuNrOfGet4PerGdpb, fuNrOfGdpbs, -0.5, fuNrOfGdpbs - 0.5);
629
630
631 /*******************************************************************/
632 for (UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb) {
633
634 /*******************************************************************/
635 name = Form("hGdpbGet4MessType_%02u", uGdpb);
636 title = Form("Nb of message for each type per GET4 in gDPB %02u; GET4 chip # ; Type", uGdpb);
637 fvhGdpbGet4MessType.push_back(new TH2I(name, title, fuNrOfGet4PerGdpb, 0., fuNrOfGet4PerGdpb, 4, 0., 4.));
638 fvhGdpbGet4MessType[uGdpb]->GetYaxis()->SetBinLabel(1, "DATA 32b");
639 fvhGdpbGet4MessType[uGdpb]->GetYaxis()->SetBinLabel(2, "EPOCH");
640 fvhGdpbGet4MessType[uGdpb]->GetYaxis()->SetBinLabel(3, "S.C. M");
641 fvhGdpbGet4MessType[uGdpb]->GetYaxis()->SetBinLabel(4, "ERROR");
642
643 /*******************************************************************/
644 name = Form("hGdpbGet4ChanScm_%02u", uGdpb);
645 title = Form("SC messages per GET4 channel in gDPB %02u; GET4 channel # ; SC type", uGdpb);
646 fvhGdpbGet4ChanScm.push_back(
647 new TH2I(name, title, 2 * fuNrOfChannelsPerGdpb, 0., fuNrOfChannelsPerGdpb, 5, 0., 5.));
648 fvhGdpbGet4ChanScm[uGdpb]->GetYaxis()->SetBinLabel(1, "Hit Scal");
649 fvhGdpbGet4ChanScm[uGdpb]->GetYaxis()->SetBinLabel(2, "Deadtime");
650 fvhGdpbGet4ChanScm[uGdpb]->GetYaxis()->SetBinLabel(3, "SPI");
651 fvhGdpbGet4ChanScm[uGdpb]->GetYaxis()->SetBinLabel(4, "SEU Scal");
652 fvhGdpbGet4ChanScm[uGdpb]->GetYaxis()->SetBinLabel(5, "START");
653
654 /*******************************************************************/
655 name = Form("hGdpbGet4ChanErrors_%02u", uGdpb);
656 title = Form("Error messages per GET4 channel in gDPB %02u; GET4 channel # ; Error", uGdpb);
657 fvhGdpbGet4ChanErrors.push_back(
658 new TH2I(name, title, fuNrOfChannelsPerGdpb, 0., fuNrOfChannelsPerGdpb, 22, 0., 22.));
659 fvhGdpbGet4ChanErrors[uGdpb]->GetYaxis()->SetBinLabel(1, "0x00: Readout Init ");
660 fvhGdpbGet4ChanErrors[uGdpb]->GetYaxis()->SetBinLabel(2, "0x01: Sync ");
661 fvhGdpbGet4ChanErrors[uGdpb]->GetYaxis()->SetBinLabel(3, "0x02: Epoch count sync");
662 fvhGdpbGet4ChanErrors[uGdpb]->GetYaxis()->SetBinLabel(4, "0x03: Epoch ");
663 fvhGdpbGet4ChanErrors[uGdpb]->GetYaxis()->SetBinLabel(5, "0x04: FIFO Write ");
664 fvhGdpbGet4ChanErrors[uGdpb]->GetYaxis()->SetBinLabel(6, "0x05: Lost event ");
665 fvhGdpbGet4ChanErrors[uGdpb]->GetYaxis()->SetBinLabel(7, "0x06: Channel state ");
666 fvhGdpbGet4ChanErrors[uGdpb]->GetYaxis()->SetBinLabel(8, "0x07: Token Ring state");
667 fvhGdpbGet4ChanErrors[uGdpb]->GetYaxis()->SetBinLabel(9, "0x08: Token ");
668 fvhGdpbGet4ChanErrors[uGdpb]->GetYaxis()->SetBinLabel(10, "0x09: Error Readout ");
669 fvhGdpbGet4ChanErrors[uGdpb]->GetYaxis()->SetBinLabel(11, "0x0a: SPI ");
670 fvhGdpbGet4ChanErrors[uGdpb]->GetYaxis()->SetBinLabel(12, "0x0b: DLL Lock error "); // <- From GET4 v1.2
671 fvhGdpbGet4ChanErrors[uGdpb]->GetYaxis()->SetBinLabel(13, "0x0c: DLL Reset invoc."); // <- From GET4 v1.2
672 fvhGdpbGet4ChanErrors[uGdpb]->GetYaxis()->SetBinLabel(14, "0x11: Overwrite "); // <- From GET4 v1.0 to 1.3
673 fvhGdpbGet4ChanErrors[uGdpb]->GetYaxis()->SetBinLabel(15, "0x12: ToT out of range");
674 fvhGdpbGet4ChanErrors[uGdpb]->GetYaxis()->SetBinLabel(16, "0x13: Event Discarded ");
675 fvhGdpbGet4ChanErrors[uGdpb]->GetYaxis()->SetBinLabel(17, "0x14: Add. Rising edge"); // <- From GET4 v1.3
676 fvhGdpbGet4ChanErrors[uGdpb]->GetYaxis()->SetBinLabel(18, "0x15: Unpaired Falling"); // <- From GET4 v1.3
677 fvhGdpbGet4ChanErrors[uGdpb]->GetYaxis()->SetBinLabel(19, "0x16: Sequence error "); // <- From GET4 v1.3
678 fvhGdpbGet4ChanErrors[uGdpb]->GetYaxis()->SetBinLabel(20, "0x17: Epoch Overflow "); // <- From GET4 v2.0
679 fvhGdpbGet4ChanErrors[uGdpb]->GetYaxis()->SetBinLabel(21, "0x7f: Unknown ");
680 fvhGdpbGet4ChanErrors[uGdpb]->GetYaxis()->SetBinLabel(22, "Corrupt/unsuprtd error");
681
682 /*******************************************************************/
683 name = Form("hGdpbPatternMissmatchEvo_%02u", uGdpb);
684 title = Form("Missmatch pattern vs TS index in gDPB %02u; TS # ; ASIC Pattern []", uGdpb);
686 new TH2I(name, title, 10000, 0., 100000, fuNrOfGet4PerGdpb, 0., fuNrOfGet4PerGdpb));
687
688 name = Form("hGdpbPatternEnableEvo_%02u", uGdpb);
689 title = Form("Enable pattern vs TS index in gDPB %02u; TS # ; ASIC Pattern []", uGdpb);
690 fvhGdpbPatternEnableEvo.push_back(
691 new TH2I(name, title, 10000, 0., 100000, fuNrOfGet4PerGdpb, 0., fuNrOfGet4PerGdpb));
692
693 name = Form("hGdpbPatternResyncEvo%02u", uGdpb);
694 title = Form("Resync pattern vs TS index in gDPB %02u; TS # ; ASIC Pattern []", uGdpb);
695 fvhGdpbPatternResyncEvo.push_back(
696 new TH2I(name, title, 10000, 0., 100000, fuNrOfGet4PerGdpb, 0., fuNrOfGet4PerGdpb));
697
699 name = Form("RawFt_gDPB_%02u", uGdpb);
700 title = Form("Raw FineTime gDPB %02u Plot 0; channel; FineTime [bin]", uGdpb);
701 fvhRawFt_gDPB.push_back(
702 new TH2F(name.Data(), title.Data(), fuNrOfChannelsPerGdpb, 0, fuNrOfChannelsPerGdpb, 128, 0, 128));
703
705 name = Form("RawTot_gDPB_%02u", uGdpb);
706 title = Form("Raw TOT gDPB %02u; channel; TOT [bin]", uGdpb);
707 fvhRawTot_gDPB.push_back(
708 new TH2F(name.Data(), title.Data(), fuNrOfChannelsPerGdpb, 0, fuNrOfChannelsPerGdpb, 256, 0, 256));
709
711 name = Form("ChCount_gDPB_%02u", uGdpb);
712 title = Form("Channel counts gDPB %02u; channel; Hits", uGdpb);
713 fvhChCount_gDPB.push_back(new TH1I(name.Data(), title.Data(), fuNrOfFeePerGdpb * fuNrOfChannelsPerFee, 0,
715
717 name = Form("ChRate_gDPB_%02u", uGdpb);
718 title = Form("Channel rate gDPB %02u; Time in run [s]; channel; Rate [1/s]", uGdpb);
719 fvhChannelRate_gDPB.push_back(new TH2D(name.Data(), title.Data(), fuHistoryHistoSize, 0, fuHistoryHistoSize,
722
724 name = Form("RemapTot_gDPB_%02u", uGdpb);
725 title = Form("Raw TOT gDPB %02u remapped; PADI channel; TOT [bin]", uGdpb);
726 fvhRemapTot_gDPB.push_back(
727 new TH2F(name.Data(), title.Data(), fuNrOfChannelsPerGdpb, 0, fuNrOfChannelsPerGdpb, 256, 0, 256));
728
730 name = Form("RemapChCount_gDPB_%02u", uGdpb);
731 title = Form("Channel counts gDPB %02u remapped; PADI channel; Hits", uGdpb);
732 fvhRemapChCount_gDPB.push_back(new TH1I(name.Data(), title.Data(), fuNrOfFeePerGdpb * fuNrOfChannelsPerFee, 0,
734
736 name = Form("RemapChRate_gDPB_%02u", uGdpb);
737 title = Form("PADI channel rate gDPB %02u; Time in run [s]; PADI channel; Rate [1/s]", uGdpb);
738 fvhRemapChRate_gDPB.push_back(new TH2D(name.Data(), title.Data(), fuHistoryHistoSize, 0, fuHistoryHistoSize,
741
743 for (UInt_t uFee = 0; uFee < fuNrOfFeePerGdpb; ++uFee) {
744 name = Form("FeeRate_gDPB_g%02u_f%1u", uGdpb, uFee);
745 title = Form("Counts per second in Fee %1u of gDPB %02u; Time[s] ; Counts", uFee, uGdpb);
746 fvhFeeRate_gDPB.push_back(new TH1D(name.Data(), title.Data(), fuHistoryHistoSize, 0, fuHistoryHistoSize));
747
748 name = Form("FeeErrorRate_gDPB_g%02u_f%1u", uGdpb, uFee);
749 title = Form("Error Counts per second in Fee %1u of gDPB %02u; Time[s] ; "
750 "Error Counts",
751 uFee, uGdpb);
752 fvhFeeErrorRate_gDPB.push_back(new TH1D(name.Data(), title.Data(), fuHistoryHistoSize, 0, fuHistoryHistoSize));
753
754 name = Form("FeeErrorRatio_gDPB_g%02u_f%1u", uGdpb, uFee);
755 title = Form("Error to data ratio per second in Fee %1u of gDPB %02u; "
756 "Time[s] ; Error ratio[]",
757 uFee, uGdpb);
758 fvhFeeErrorRatio_gDPB.push_back(
759 new TProfile(name.Data(), title.Data(), fuHistoryHistoSize, 0, fuHistoryHistoSize));
760
761 name = Form("FeeRateLong_gDPB_g%02u_f%1u", uGdpb, uFee);
762 title = Form("Counts per minutes in Fee %1u of gDPB %02u; Time[min] ; Counts", uFee, uGdpb);
763 fvhFeeRateLong_gDPB.push_back(
764 new TH1D(name.Data(), title.Data(), fuHistoryHistoSizeLong, 0, fuHistoryHistoSizeLong));
765
766 name = Form("FeeErrorRateLong_gDPB_g%02u_f%1u", uGdpb, uFee);
767 title = Form("Error Counts per minutes in Fee %1u of gDPB %02u; "
768 "Time[min] ; Error Counts",
769 uFee, uGdpb);
770 fvhFeeErrorRateLong_gDPB.push_back(
771 new TH1D(name.Data(), title.Data(), fuHistoryHistoSizeLong, 0, fuHistoryHistoSizeLong));
772
773 name = Form("FeeErrorRatioLong_gDPB_g%02u_f%1u", uGdpb, uFee);
774 title = Form("Error to data ratio per minutes in Fee %1u of gDPB %02u; "
775 "Time[min] ; Error ratio[]",
776 uFee, uGdpb);
778 new TProfile(name.Data(), title.Data(), fuHistoryHistoSizeLong, 0, fuHistoryHistoSizeLong));
779 } // for (UInt_t uFee = 0; uFee < fuNrOfFeePerGdpb; uFee++)
780
783 name = Form("hTokenMsgType_gDPB_%02u", uGdpb);
784 title = Form("STAR trigger Messages type gDPB %02u; Type ; Counts", uGdpb);
785 fvhTokenMsgType.push_back(new TH1F(name, title, 4, 0, 4));
786 fvhTokenMsgType[uGdpb]->GetXaxis()->SetBinLabel(1, "A"); // gDPB TS high
787 fvhTokenMsgType[uGdpb]->GetXaxis()->SetBinLabel(2, "B"); // gDPB TS low, STAR TS high
788 fvhTokenMsgType[uGdpb]->GetXaxis()->SetBinLabel(3, "C"); // STAR TS mid
789 fvhTokenMsgType[uGdpb]->GetXaxis()->SetBinLabel(4, "D"); // STAR TS low, token, CMDs
790 if (server) server->Register("/StarRaw", fvhTokenMsgType[uGdpb]);
791
792 name = Form("hTriggerRate_gDPB_%02u", uGdpb);
793 title = Form("STAR trigger signals per second gDPB %02u; Time[s] ; Counts", uGdpb);
794 fvhTriggerRate.push_back(new TH1F(name, title, fuHistoryHistoSize, 0, fuHistoryHistoSize));
795 if (server) server->Register("/StarRaw", fvhTriggerRate[uGdpb]);
796
797 name = Form("hCmdDaqVsTrig_gDPB_%02u", uGdpb);
798 title = Form("STAR daq command VS STAR trigger command gDPB %02u; DAQ ; TRIGGER", uGdpb);
799 fvhCmdDaqVsTrig.push_back(new TH2I(name, title, 16, 0, 16, 16, 0, 16));
800 fvhCmdDaqVsTrig[uGdpb]->GetXaxis()->SetBinLabel(1, "0x0: no-trig "); // idle link
801 fvhCmdDaqVsTrig[uGdpb]->GetXaxis()->SetBinLabel(
802 2, "0x1: clear "); // clears redundancy counters on the readout boards
803 fvhCmdDaqVsTrig[uGdpb]->GetXaxis()->SetBinLabel(3, "0x2: mast-rst"); // general reset of the whole front-end logic
804 fvhCmdDaqVsTrig[uGdpb]->GetXaxis()->SetBinLabel(4, "0x3: spare "); // reserved
805 fvhCmdDaqVsTrig[uGdpb]->GetXaxis()->SetBinLabel(
806 5, "0x4: trigg. 0"); // Default physics readout, all det support required
807 fvhCmdDaqVsTrig[uGdpb]->GetXaxis()->SetBinLabel(6, "0x5: trigg. 1"); //
808 fvhCmdDaqVsTrig[uGdpb]->GetXaxis()->SetBinLabel(7, "0x6: trigg. 2"); //
809 fvhCmdDaqVsTrig[uGdpb]->GetXaxis()->SetBinLabel(8, "0x7: trigg. 3"); //
810 fvhCmdDaqVsTrig[uGdpb]->GetXaxis()->SetBinLabel(9, "0x8: puls. 0"); //
811 fvhCmdDaqVsTrig[uGdpb]->GetXaxis()->SetBinLabel(10, "0x9: puls. 1"); //
812 fvhCmdDaqVsTrig[uGdpb]->GetXaxis()->SetBinLabel(11, "0xA: puls. 2"); //
813 fvhCmdDaqVsTrig[uGdpb]->GetXaxis()->SetBinLabel(12, "0xB: puls. 3"); //
814 fvhCmdDaqVsTrig[uGdpb]->GetXaxis()->SetBinLabel(
815 13,
816 "0xC: config "); // housekeeping trigger: return geographic info of FE
817 fvhCmdDaqVsTrig[uGdpb]->GetXaxis()->SetBinLabel(14, "0xD: abort "); // aborts and clears an active event
818 fvhCmdDaqVsTrig[uGdpb]->GetXaxis()->SetBinLabel(15, "0xE: L1accept"); //
819 fvhCmdDaqVsTrig[uGdpb]->GetXaxis()->SetBinLabel(16, "0xF: L2accept"); //
820 fvhCmdDaqVsTrig[uGdpb]->GetYaxis()->SetBinLabel(1, "0x0: 0"); // To be filled at STAR
821 fvhCmdDaqVsTrig[uGdpb]->GetYaxis()->SetBinLabel(2, "0x1: 1"); // To be filled at STAR
822 fvhCmdDaqVsTrig[uGdpb]->GetYaxis()->SetBinLabel(3, "0x2: 2"); // To be filled at STAR
823 fvhCmdDaqVsTrig[uGdpb]->GetYaxis()->SetBinLabel(4, "0x3: 3"); // To be filled at STAR
824 fvhCmdDaqVsTrig[uGdpb]->GetYaxis()->SetBinLabel(5, "0x4: 4"); // To be filled at STAR
825 fvhCmdDaqVsTrig[uGdpb]->GetYaxis()->SetBinLabel(6, "0x5: 5"); // To be filled at STAR
826 fvhCmdDaqVsTrig[uGdpb]->GetYaxis()->SetBinLabel(7, "0x6: 6"); // To be filled at STAR
827 fvhCmdDaqVsTrig[uGdpb]->GetYaxis()->SetBinLabel(8, "0x7: 7"); // To be filled at STAR
828 fvhCmdDaqVsTrig[uGdpb]->GetYaxis()->SetBinLabel(9, "0x8: 8"); // To be filled at STAR
829 fvhCmdDaqVsTrig[uGdpb]->GetYaxis()->SetBinLabel(10, "0x9: 9"); // To be filled at STAR
830 fvhCmdDaqVsTrig[uGdpb]->GetYaxis()->SetBinLabel(11, "0xA: 10"); // To be filled at STAR
831 fvhCmdDaqVsTrig[uGdpb]->GetYaxis()->SetBinLabel(12, "0xB: 11"); // To be filled at STAR
832 fvhCmdDaqVsTrig[uGdpb]->GetYaxis()->SetBinLabel(13, "0xC: 12"); // To be filled at STAR
833 fvhCmdDaqVsTrig[uGdpb]->GetYaxis()->SetBinLabel(14, "0xD: 13"); // To be filled at STAR
834 fvhCmdDaqVsTrig[uGdpb]->GetYaxis()->SetBinLabel(15, "0xE: 14"); // To be filled at STAR
835 fvhCmdDaqVsTrig[uGdpb]->GetYaxis()->SetBinLabel(16, "0xF: 15"); // To be filled at STAR
836 if (server) server->Register("/StarRaw", fvhCmdDaqVsTrig[uGdpb]);
837
838 name = Form("hStarTokenEvo_gDPB_%02u", uGdpb);
839 title = Form("STAR token value VS time gDPB %02u; Time in Run [s] ; STAR "
840 "Token; Counts",
841 uGdpb);
842 fvhStarTokenEvo.push_back(new TH2I(name, title, fuHistoryHistoSize, 0, fuHistoryHistoSize, 410, 0, 4100));
843
844
845 name = Form("hStarTrigGdpbTsEvo_gDPB_%02u", uGdpb);
846 title = Form("gDPB TS in STAR triger tokens for gDPB %02u; Time in Run [s] ; gDPB TS;", uGdpb);
847 fvhStarTrigGdpbTsEvo.push_back(new TProfile(name, title, fuHistoryHistoSize, 0, fuHistoryHistoSize));
848
849 name = Form("hStarTrigStarTsEvo_gDPB_%02u", uGdpb);
850 title = Form("STAR TS in STAR triger tokens for gDPB %02u; Time in Run [s] ; STAR TS;", uGdpb);
851 fvhStarTrigStarTsEvo.push_back(new TProfile(name, title, fuHistoryHistoSize, 0, fuHistoryHistoSize));
852 } // for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
853
854 /*******************************************************************/
856 for (UInt_t uMod = 0; uMod < fuNrOfModules; uMod++) {
857 name = Form("RemapTotSideA_mod_%02u", uMod);
858 title = Form("Raw TOT module %02u Side A; PADI channel; TOT [bin]", uMod);
859 fvhRemapTotSideA_mod.push_back(new TH2F(name.Data(), title.Data(), kuNbFeeSide * fuNrOfChannelsPerFee, 0,
860 kuNbFeeSide * fuNrOfChannelsPerFee, 256, 0, 256));
861 name = Form("RemapTotSideB_mod_%02u", uMod);
862 title = Form("Raw TOT module %02u Side B; PADI channel; TOT [bin]", uMod);
863 fvhRemapTotSideB_mod.push_back(new TH2F(name.Data(), title.Data(), kuNbFeeSide * fuNrOfChannelsPerFee, 0,
864 kuNbFeeSide * fuNrOfChannelsPerFee, 256, 0, 256));
865
866 name = Form("ModRate_gDPB_m%02u", uMod);
867 title = Form("Counts per second in Module %02u; Time[s] ; Counts", uMod);
868 fvhModRate.push_back(new TH1D(name.Data(), title.Data(), fuHistoryHistoSize, 0, fuHistoryHistoSize));
869
870 name = Form("ModErrorRate_m%02u", uMod);
871 title = Form("Error Counts per second in Module %02u; Time[s] ; Error Counts", uMod);
872 fvhModErrorRate.push_back(new TH1D(name.Data(), title.Data(), fuHistoryHistoSize, 0, fuHistoryHistoSize));
873
874 name = Form("ModErrorRatio_m%02u", uMod);
875 title = Form("Error to data ratio per second in Module %02u; Time[s] ; Error ratio[]", uMod);
876 fvhModErrorRatio.push_back(new TProfile(name.Data(), title.Data(), fuHistoryHistoSize, 0, fuHistoryHistoSize));
877 } // for( UInt_t uMod = 0; uMod < fuNrOfModules; uMod ++ )
878
879 /*******************************************************************/
881 if (kTRUE == fbPulserModeEnable) {
883 for (UInt_t uFeeA = 0; uFeeA < fuNrOfFeePerGdpb * fuNrOfGdpbs; uFeeA++) {
885 for (UInt_t uFeeB = 0; uFeeB < fuNrOfFeePerGdpb * fuNrOfGdpbs; uFeeB++) {
886 if (uFeeA < uFeeB) {
887 UInt_t uGdpbA = uFeeA / (fuNrOfFeePerGdpb);
888 UInt_t uFeeIdA = uFeeA - (fuNrOfFeePerGdpb * uGdpbA);
889 UInt_t uGdpbB = uFeeB / (fuNrOfFeePerGdpb);
890 UInt_t uFeeIdB = uFeeB - (fuNrOfFeePerGdpb * uGdpbB);
891 fvhTimeDiffPulser[uFeeA][uFeeB] =
892 new TH1I(Form("hTimeDiffPulser_g%02u_f%1u_g%02u_f%1u", uGdpbA, uFeeIdA, uGdpbB, uFeeIdB),
893 Form("Time difference for pulser on gDPB %02u FEE %1u and "
894 "gDPB %02u FEE %1u; DeltaT [ps]; Counts",
895 uGdpbA, uFeeIdA, uGdpbB, uFeeIdB),
896 uNbBinsDt, dMinDt, dMaxDt);
897 } // if( uFeeA < uFeeB )
898 else
899 fvhTimeDiffPulser[uFeeA][uFeeB] = NULL;
900 } // for( UInt_t uFeeB = uFeeA; uFeeB < fuNrOfFeePerGdpb * fuNrOfGdpbs - 1; uFeeB++)
901 } // for( UInt_t uFeeA = 0; uFeeA < kuNbChanTest - 1; uFeeA++)
902
903 name = "hTimeMeanPulser";
904 fhTimeMeanPulser = new TH2D(name.Data(), "Time difference Mean for each FEE pairs; FEE A; FEE B ; Mean [ps]",
907
908 name = "hTimeRmsPulser";
909 fhTimeRmsPulser = new TH2D(name.Data(), "Time difference RMS for each FEE pairs; FEE A; FEE B ; RMS [ps]",
912
913 name = "hTimeRmsZoomFitPuls";
914 fhTimeRmsZoomFitPuls = new TH2D(name.Data(),
915 "Time difference RMS after zoom for each "
916 "FEE pairs; FEE A; FEE B ; RMS [ps]",
919
920 name = "hTimeResFitPuls";
921 fhTimeResFitPuls = new TH2D(name.Data(),
922 "Time difference Res from fit for each FEE "
923 "pairs; FEE A; FEE B ; Sigma [ps]",
926
929 for (UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb) {
930 for (UInt_t uGbtx = 0; uGbtx < kuNbGbtxPerGdpb - 1; ++uGbtx) {
931 name = Form("hPulserTimeDiffEvoGdpb%02uGbtx00Gbtx%02u", uGdpb, uGbtx + 1);
932 fvhPulserTimeDiffEvoGbtxGbtx[uGdpb * (kuNbGbtxPerGdpb - 1) + uGbtx] =
933 new TProfile(name.Data(),
934 Form("Time difference of the 1st FEE in the 1st GBTx of gDPB %02u "
935 "vs GBTx %02u; time in run [min]; dt [ps]",
936 uGdpb, uGbtx + 1),
938 } // for( UInt_t uGbtx = 0; uGbtx < kuNbGbtxPerGdpb; ++uGbtx )
939
940 fvvhPulserTimeDiffEvoGdpbGdpb[uGdpb].resize(fuNrOfGdpbs, NULL);
941 for (UInt_t uGdpbB = uGdpb + 1; uGdpbB < fuNrOfGdpbs; ++uGdpbB) {
942 name = Form("hPulserTimeDiffEvoGdpb%02uGdpb%02u", uGdpb, uGdpbB);
943 fvvhPulserTimeDiffEvoGdpbGdpb[uGdpb][uGdpbB] =
944 new TProfile(name.Data(),
945 Form("Time difference of the 1st FEE in the 1st GBTx of "
946 "gDPB %02u vs %02u; time in run [min]; dt [ps]",
947 uGdpb, uGdpbB),
949 } // for( UInt_t uGdpbB = uGdpb + 1; uGdpbB < fuNrOfGdpbs; ++uGdpbB )
950 } // for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
951 } // if( kTRUE == fbPulserModeEnable )
952
954 if (kTRUE == fbCoincMapsEnable) {
955 for (UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb) {
956 name = Form("hCoincMapAllChanGdpb%02u", uGdpb);
957 fvhCoincMapAllChanGdpb.push_back(new TH2D(name.Data(),
958 Form("Coincidence map of all channels of gDPB %02u; Chan A "
959 "[]; Chan B[]; Coinc. []",
960 uGdpb),
963
964 name = Form("hCoincMeanAllChanGdpb%02u", uGdpb);
965 fvhCoincMeanAllChanGdpb.push_back(new TProfile2D(name.Data(),
966 Form("Coincidence mean of all channels of gDPB %02u; "
967 "Chan A []; Chan B[]; Mean dt [ps]",
968 uGdpb),
971 } // for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
972 } // if( kTRUE == fbCoincMapsEnable )
973
974 if (server) {
975 server->Register("/TofRaw", fhMessType);
976 server->Register("/TofRaw", fhSysMessType);
977 server->Register("/TofRaw", fhGet4MessType);
978 server->Register("/TofRaw", fhGet4ChanScm);
979 server->Register("/TofRaw", fhGet4ChanErrors);
980 server->Register("/TofRaw", fhGet4EpochFlags);
981
982 server->Register("/TofRaw", fhGdpbMessType);
983 server->Register("/TofRaw", fhGdpbSysMessType);
984 server->Register("/TofRaw", fhGdpbSysMessPattType);
985 server->Register("/TofRaw", fhGdpbEpochFlags);
986 server->Register("/TofRaw", fhGdpbEpochSyncEvo);
987 server->Register("/TofRaw", fhGdpbEpochMissEvo);
988
989 server->Register("/TofRaw", fhScmScalerCounters);
990 server->Register("/TofRaw", fhScmDeadtimeCounters);
991 server->Register("/TofRaw", fhScmSeuCounters);
992 server->Register("/TofRaw", fhScmSeuCountersEvo);
993
994 server->Register("/TofRaw", fhPatternMissmatch);
995 server->Register("/TofRaw", fhPatternEnable);
996 server->Register("/TofRaw", fhPatternResync);
997
998 for (UInt_t uTotPlot = 0; uTotPlot < fvhRawTot_gDPB.size(); ++uTotPlot)
999 server->Register("/TofRaw", fvhRawTot_gDPB[uTotPlot]);
1000
1001 for (UInt_t uTotPlot = 0; uTotPlot < fvhRemapTot_gDPB.size(); ++uTotPlot)
1002 server->Register("/TofRaw", fvhRemapTot_gDPB[uTotPlot]);
1003
1004 for (UInt_t uMod = 0; uMod < fuNrOfModules; uMod++) {
1005 server->Register("/TofRaw", fvhRemapTotSideA_mod[uMod]);
1006 server->Register("/TofRaw", fvhRemapTotSideB_mod[uMod]);
1007 server->Register("/TofMod", fvhModRate[uMod]);
1008 server->Register("/TofMod", fvhModErrorRate[uMod]);
1009 server->Register("/TofMod", fvhModErrorRatio[uMod]);
1010 } // for( UInt_t uMod = 0; uMod < fuNrOfModules; uMod ++ )
1011
1012 for (UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb) {
1013 server->Register("/TofRaw", fvhGdpbGet4MessType[uGdpb]);
1014 server->Register("/TofRaw", fvhGdpbGet4ChanScm[uGdpb]);
1015 server->Register("/TofRaw", fvhGdpbGet4ChanErrors[uGdpb]);
1016
1017 server->Register("/TofRaw", fvhGdpbPatternMissmatchEvo[uGdpb]);
1018 server->Register("/TofRaw", fvhGdpbPatternEnableEvo[uGdpb]);
1019 server->Register("/TofRaw", fvhGdpbPatternResyncEvo[uGdpb]);
1020
1021 server->Register("/TofRaw", fvhRawFt_gDPB[uGdpb]);
1022 server->Register("/TofRaw", fvhChCount_gDPB[uGdpb]);
1023 server->Register("/TofRates", fvhChannelRate_gDPB[uGdpb]);
1024 server->Register("/TofRaw", fvhRemapChCount_gDPB[uGdpb]);
1025 server->Register("/TofRates", fvhRemapChRate_gDPB[uGdpb]);
1026
1027 for (UInt_t uFee = 0; uFee < fuNrOfFeePerGdpb; ++uFee) {
1028 server->Register("/TofRates", fvhFeeRate_gDPB[uGdpb * fuNrOfFeePerGdpb + uFee]);
1029 server->Register("/TofRates", fvhFeeErrorRate_gDPB[uGdpb * fuNrOfFeePerGdpb + uFee]);
1030 server->Register("/TofRates", fvhFeeErrorRatio_gDPB[uGdpb * fuNrOfFeePerGdpb + uFee]);
1031 server->Register("/TofRates", fvhFeeRateLong_gDPB[uGdpb * fuNrOfFeePerGdpb + uFee]);
1032 server->Register("/TofRates", fvhFeeErrorRateLong_gDPB[uGdpb * fuNrOfFeePerGdpb + uFee]);
1033 server->Register("/TofRates", fvhFeeErrorRatioLong_gDPB[uGdpb * fuNrOfFeePerGdpb + uFee]);
1034 } // for (UInt_t uFee = 0; uFee < fuNrOfFeePerGdpb; ++ uFee)
1035
1036 server->Register("/StarRaw", fvhTokenMsgType[uGdpb]);
1037 server->Register("/StarRaw", fvhTriggerRate[uGdpb]);
1038 server->Register("/StarRaw", fvhCmdDaqVsTrig[uGdpb]);
1039 server->Register("/StarRaw", fvhStarTokenEvo[uGdpb]);
1040 server->Register("/StarRaw", fvhStarTrigGdpbTsEvo[uGdpb]);
1041 server->Register("/StarRaw", fvhStarTrigStarTsEvo[uGdpb]);
1042 } // for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
1043
1044 if (kTRUE == fbPulserModeEnable) {
1045 for (UInt_t uFeeA = 0; uFeeA < fuNrOfFeePerGdpb * fuNrOfGdpbs; uFeeA++)
1046 for (UInt_t uFeeB = 0; uFeeB < fuNrOfFeePerGdpb * fuNrOfGdpbs; uFeeB++)
1047 if (NULL != fvhTimeDiffPulser[uFeeA][uFeeB]) server->Register("/TofDt", fvhTimeDiffPulser[uFeeA][uFeeB]);
1048
1049 server->Register("/TofRaw", fhTimeMeanPulser);
1050 server->Register("/TofRaw", fhTimeRmsPulser);
1051 server->Register("/TofRaw", fhTimeRmsZoomFitPuls);
1052 server->Register("/TofRaw", fhTimeResFitPuls);
1053
1054 for (UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb) {
1055 for (UInt_t uGbtx = 0; uGbtx < kuNbGbtxPerGdpb - 1; ++uGbtx)
1056 if (NULL != fvhPulserTimeDiffEvoGbtxGbtx[uGdpb * (kuNbGbtxPerGdpb - 1) + uGbtx])
1057 server->Register("/TofDtEvo", fvhPulserTimeDiffEvoGbtxGbtx[uGdpb * (kuNbGbtxPerGdpb - 1) + uGbtx]);
1058
1059 for (UInt_t uGdpbB = uGdpb + 1; uGdpbB < fuNrOfGdpbs; ++uGdpbB)
1060 if (NULL != fvvhPulserTimeDiffEvoGdpbGdpb[uGdpb][uGdpbB])
1061 server->Register("/TofDtEvo", fvvhPulserTimeDiffEvoGdpbGdpb[uGdpb][uGdpbB]);
1062 } // for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
1063 } // if( kTRUE == fbPulserModeEnable )
1064 if (kTRUE == fbCoincMapsEnable) {
1065 for (UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb) {
1066 server->Register("/TofCoinc", fvhCoincMapAllChanGdpb[uGdpb]);
1067 server->Register("/TofCoinc", fvhCoincMeanAllChanGdpb[uGdpb]);
1068 } // for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
1069 } // if( kTRUE == fbCoincMapsEnable )
1070
1071 server->RegisterCommand("/Reset_All_eTOF", "bMcbmMoniTofResetHistos=kTRUE");
1072 server->RegisterCommand("/Save_All_eTof", "bMcbmMoniTofSaveHistos=kTRUE");
1073 server->RegisterCommand("/Update_PulsFit", "bMcbmMoniTofUpdateZoomedFit=kTRUE");
1074 server->RegisterCommand("/Print_Raw_Data", "bMcbmMoniTofRawDataPrint=kTRUE");
1075 server->RegisterCommand("/Print_AllHits", "bMcbmMoniTofPrintAllHitsEna=kTRUE");
1076 server->RegisterCommand("/Print_AllEps", "bMcbmMoniTofPrintAllEpochsEna=kTRUE");
1077
1078 server->Restrict("/Reset_All_eTof", "allow=admin");
1079 server->Restrict("/Save_All_eTof", "allow=admin");
1080 server->Restrict("/Update_PulsFit", "allow=admin");
1081 server->Restrict("/Print_Raw_Data", "allow=admin");
1082 server->Restrict("/Print_AllHits", "allow=admin");
1083 server->Restrict("/Print_AllEps", "allow=admin");
1084 } // if( server )
1085
1087 Double_t w = 10;
1088 Double_t h = 10;
1089 TCanvas* cSummary = new TCanvas("cSummary", "gDPB Monitoring Summary", w, h);
1090 cSummary->Divide(2, 3);
1091
1092 // 1st Column: Messages types
1093 cSummary->cd(1);
1094 gPad->SetLogy();
1095 fhMessType->Draw();
1096
1097 cSummary->cd(2);
1098 gPad->SetLogy();
1099 fhSysMessType->Draw();
1100
1101 cSummary->cd(3);
1102 gPad->SetLogz();
1103 fhGet4MessType->Draw("colz");
1104
1105 // 2nd Column: GET4 Errors + Epoch flags + SCm
1106 cSummary->cd(4);
1107 gPad->SetLogz();
1108 fhGet4ChanErrors->Draw("colz");
1109
1110 cSummary->cd(5);
1111 gPad->SetLogz();
1112 fhGet4EpochFlags->Draw("colz");
1113
1114 cSummary->cd(6);
1115 fhGet4ChanScm->Draw("colz");
1116 /*****************************/
1117
1119 TCanvas* cSummaryGdpb = new TCanvas("cSummaryGdpb", "gDPB Monitoring Summary", w, h);
1120 cSummaryGdpb->Divide(2, 3);
1121
1122 cSummaryGdpb->cd(1);
1123 gPad->SetLogz();
1124 fhGdpbMessType->Draw("colz");
1125
1126 cSummaryGdpb->cd(3);
1127 gPad->SetLogz();
1128 fhGdpbSysMessType->Draw("colz");
1129
1130 cSummaryGdpb->cd(5);
1131 gPad->SetLogz();
1132 fhGdpbEpochFlags->Draw("text colz");
1133
1134 cSummaryGdpb->cd(4);
1135 fhGdpbEpochSyncEvo->Draw("colz");
1136
1137 cSummaryGdpb->cd(6);
1138 gPad->SetLogz();
1139 fhGdpbEpochMissEvo->Draw("colz");
1140 /*****************************/
1141
1143 for (UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb) {
1144 TCanvas* cSumGdpbGet4 =
1145 new TCanvas(Form("cSumGdpbGet4_g%02u", uGdpb), Form("Summary per GET4 or channel for gDPB %02u", uGdpb), w, h);
1146 cSumGdpbGet4->Divide(2, 2);
1147
1148 cSumGdpbGet4->cd(1);
1149 gPad->SetLogz();
1150 fvhGdpbGet4MessType[uGdpb]->Draw("colz");
1151
1152 cSumGdpbGet4->cd(2);
1153 gPad->SetLogz();
1154 fvhGdpbGet4ChanScm[uGdpb]->Draw("colz");
1155
1156 cSumGdpbGet4->cd(3);
1157 gPad->SetLogz();
1158 fvhGdpbGet4ChanErrors[uGdpb]->Draw("colz");
1159 } // for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
1160 /*****************************/
1161
1163 TCanvas* cFeeRates = new TCanvas("cFeeRates", "gDPB Monitoring FEE rates", w, h);
1164 // cFeeRates->Divide(fuNrOfFeePerGdpb, fuNrOfGdpbs );
1165 cFeeRates->Divide(kuNbFeePerGbtx, kuNbGbtxPerGdpb * fuNrOfGdpbs);
1166 for (UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb) {
1167 for (UInt_t uFee = 0; uFee < fuNrOfFeePerGdpb; ++uFee) {
1168 cFeeRates->cd(1 + uGdpb * fuNrOfFeePerGdpb + uFee);
1169 gPad->SetLogy();
1170 fvhFeeRate_gDPB[uGdpb * fuNrOfFeePerGdpb + uFee]->Draw("hist");
1171
1172 fvhFeeErrorRate_gDPB[uGdpb * fuNrOfFeePerGdpb + uFee]->SetLineColor(kRed);
1173 fvhFeeErrorRate_gDPB[uGdpb * fuNrOfFeePerGdpb + uFee]->Draw("same hist");
1174 } // for (UInt_t uFee = 0; uFee < fuNrOfFeePerGdpb; ++uFee )
1175 } // for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
1176 /*****************************/
1177
1179 TCanvas* cFeeErrRatio = new TCanvas("cFeeErrRatio", "gDPB Monitoring FEE error ratios", w, h);
1180 // cFeeErrRatio->Divide(fuNrOfFeePerGdpb, fuNrOfGdpbs );
1181 cFeeErrRatio->Divide(kuNbFeePerGbtx, kuNbGbtxPerGdpb * fuNrOfGdpbs);
1182 for (UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb) {
1183 for (UInt_t uFee = 0; uFee < fuNrOfFeePerGdpb; ++uFee) {
1184 cFeeErrRatio->cd(1 + uGdpb * fuNrOfFeePerGdpb + uFee);
1185 gPad->SetLogy();
1186 fvhFeeErrorRatio_gDPB[uGdpb * fuNrOfFeePerGdpb + uFee]->Draw("hist le0");
1187 } // for (UInt_t uFee = 0; uFee < fuNrOfFeePerGdpb; ++uFee )
1188 } // for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
1189 /*****************************/
1190
1191
1193 TCanvas* cFeeRatesLong = new TCanvas("cFeeRatesLong", "gDPB Monitoring FEE rates", w, h);
1194 // cFeeRatesLong->Divide(fuNrOfFeePerGdpb, fuNrOfGdpbs );
1195 cFeeRatesLong->Divide(kuNbFeePerGbtx, kuNbGbtxPerGdpb * fuNrOfGdpbs);
1196 for (UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb) {
1197 for (UInt_t uFee = 0; uFee < fuNrOfFeePerGdpb; ++uFee) {
1198 cFeeRatesLong->cd(1 + uGdpb * fuNrOfFeePerGdpb + uFee);
1199 gPad->SetLogy();
1200 fvhFeeRateLong_gDPB[uGdpb]->Draw("hist");
1201
1202 fvhFeeErrorRateLong_gDPB[uGdpb * fuNrOfFeePerGdpb + uFee]->SetLineColor(kRed);
1203 fvhFeeErrorRateLong_gDPB[uGdpb * fuNrOfFeePerGdpb + uFee]->Draw("same hist");
1204 } // for (UInt_t uFee = 0; uFee < fuNrOfFeePerGdpb; ++uFee )
1205 } // for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
1206 /*****************************/
1207
1209 TCanvas* cFeeErrRatioLong = new TCanvas("cFeeErrRatioLong", "gDPB Monitoring FEE error ratios", w, h);
1210 // cFeeErrRatioLong->Divide(fuNrOfFeePerGdpb, fuNrOfGdpbs );
1211 cFeeErrRatioLong->Divide(kuNbFeePerGbtx, kuNbGbtxPerGdpb * fuNrOfGdpbs);
1212 for (UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb) {
1213 for (UInt_t uFee = 0; uFee < fuNrOfFeePerGdpb; ++uFee) {
1214 cFeeErrRatioLong->cd(1 + uGdpb * fuNrOfFeePerGdpb + uFee);
1215 gPad->SetLogy();
1216 fvhFeeErrorRatioLong_gDPB[uGdpb * fuNrOfFeePerGdpb + uFee]->Draw("hist le0");
1217 } // for (UInt_t uFee = 0; uFee < fuNrOfFeePerGdpb; ++uFee )
1218 } // for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
1219 /*****************************/
1220
1222 TCanvas* cGdpbChannelCount = new TCanvas("cGdpbChannelCount", "Integrated Get4 channel counts per gDPB", w, h);
1223 cGdpbChannelCount->Divide(1, fuNrOfGdpbs);
1224 for (UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb) {
1225 cGdpbChannelCount->cd(1 + uGdpb);
1226 gPad->SetGridx();
1227 gPad->SetGridy();
1228 gPad->SetLogy();
1229 fvhChCount_gDPB[uGdpb]->Draw();
1230 } // for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
1231 /*****************************/
1232
1234 TCanvas* cGdpbRemapChCount = new TCanvas("cGdpbRemapChCount", "Integrated PADI channel counts per gDPB", w, h);
1235 cGdpbRemapChCount->Divide(1, fuNrOfGdpbs);
1236 for (UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb) {
1237 cGdpbRemapChCount->cd(1 + uGdpb);
1238 gPad->SetGridx();
1239 gPad->SetGridy();
1240 gPad->SetLogy();
1241 fvhRemapChCount_gDPB[uGdpb]->Draw();
1242 } // for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
1243 /*****************************/
1244
1246 TCanvas* cGdpbChannelRate = new TCanvas("cGdpbChannelRate", "Get4 channel rate per gDPB", w, h);
1247 cGdpbChannelRate->Divide(1, fuNrOfGdpbs);
1248 for (UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb) {
1249 cGdpbChannelRate->cd(1 + uGdpb);
1250 gPad->SetGridx();
1251 gPad->SetGridy();
1252 gPad->SetLogz();
1253 fvhChannelRate_gDPB[uGdpb]->Draw("colz");
1254 } // for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
1255 /*****************************/
1256
1258 TCanvas* cGdpbRemapChRate = new TCanvas("cGdpbRemapChRate", "PADI channel rate per gDPB", w, h);
1259 cGdpbRemapChRate->Divide(1, fuNrOfGdpbs);
1260 for (UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb) {
1261 cGdpbRemapChRate->cd(1 + uGdpb);
1262 gPad->SetGridx();
1263 gPad->SetGridy();
1264 gPad->SetLogz();
1265 fvhRemapChRate_gDPB[uGdpb]->Draw("colz");
1266 } // for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
1267 /*****************************/
1268
1270 TCanvas* cTotPnt = NULL;
1271 for (UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb) {
1272 cTotPnt = new TCanvas(Form("cTot_g%02u", uGdpb), Form("gDPB %02u TOT distributions", uGdpb), w, h);
1273
1274 cTotPnt->cd();
1275 gPad->SetGridx();
1276 gPad->SetGridy();
1277 gPad->SetLogz();
1278
1279 fvhRawTot_gDPB[uGdpb]->Draw("colz");
1280 } // for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
1281 cTotPnt = new TCanvas("cTot_all", "TOT distributions", w, h);
1282 cTotPnt->Divide(fuNrOfGdpbs);
1283 for (UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb) {
1284 cTotPnt->cd(1 + uGdpb);
1285 gPad->SetGridx();
1286 gPad->SetGridy();
1287 gPad->SetLogz();
1288
1289 fvhRawTot_gDPB[uGdpb]->Draw("colz");
1290 } // for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
1291
1292 /**************************************************/
1293
1295 cTotPnt = NULL;
1296 for (UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb) {
1297 cTotPnt = new TCanvas(Form("cTotRemap_g%02u", uGdpb), Form("PADI ch gDPB %02u TOT distributions", uGdpb), w, h);
1298
1299 cTotPnt->cd();
1300 gPad->SetGridx();
1301 gPad->SetGridy();
1302 gPad->SetLogz();
1303
1304 fvhRemapTot_gDPB[uGdpb]->Draw("colz");
1305 } // for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
1306 cTotPnt = new TCanvas("cTotRemap_all", "TOT distributions", w, h);
1307 cTotPnt->Divide(fuNrOfGdpbs);
1308 for (UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb) {
1309 cTotPnt->cd(1 + uGdpb);
1310 gPad->SetGridx();
1311 gPad->SetGridy();
1312 gPad->SetLogz();
1313
1314 fvhRemapTot_gDPB[uGdpb]->Draw("colz");
1315 } // for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
1316
1317 /**************************************************/
1318
1320 cTotPnt = NULL;
1321 for (UInt_t uMod = 0; uMod < fuNrOfModules; uMod++) {
1322 cTotPnt =
1323 new TCanvas(Form("cTotRemapSides_m%02u", uMod), Form("Sides ch module %02u TOT distributions", uMod), w, h);
1324 cTotPnt->Divide(1, 2);
1325
1326 cTotPnt->cd(1);
1327 gPad->SetGridx();
1328 gPad->SetGridy();
1329 gPad->SetLogz();
1330
1331 fvhRemapTotSideA_mod[uMod]->Draw("colz");
1332
1333 cTotPnt->cd(2);
1334 gPad->SetGridx();
1335 gPad->SetGridy();
1336 gPad->SetLogz();
1337 fvhRemapTotSideB_mod[uMod]->Draw("colz");
1338 } // for( UInt_t uMod = 0; uMod < fuNrOfModules; uMod ++ )
1339
1340 /**************************************************/
1341
1343 for (UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; uGdpb++) {
1344 TCanvas* cStarToken =
1345 new TCanvas(Form("cStarToken_g%02u", uGdpb), Form("STAR token detection info for gDPB %02u", uGdpb), w, h);
1346 cStarToken->Divide(2, 2);
1347
1348 cStarToken->cd(1);
1349 fvhTriggerRate[uGdpb]->Draw();
1350
1351 cStarToken->cd(2);
1352 fvhCmdDaqVsTrig[uGdpb]->Draw("colz");
1353
1354 cStarToken->cd(3);
1355 fvhStarTokenEvo[uGdpb]->Draw();
1356
1357 cStarToken->cd(4);
1358 fvhStarTrigGdpbTsEvo[uGdpb]->Draw("hist le0");
1359 fvhStarTrigStarTsEvo[uGdpb]->SetLineColor(kRed);
1360 fvhStarTrigStarTsEvo[uGdpb]->Draw("same hist le0");
1361 } // for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; uGdpb ++)
1362 /*****************************/
1363
1364 if (kTRUE == fbPulserModeEnable) {
1366 TCanvas* cPulser =
1367 new TCanvas("cPulser", "Time difference RMS for pulser channels when FEE pulser mode is ON", w, h);
1368 cPulser->Divide(2, 2);
1369
1370 cPulser->cd(1);
1371 gPad->SetGridx();
1372 gPad->SetGridy();
1373 fhTimeRmsPulser->Draw("colz");
1374
1375 cPulser->cd(2);
1376 gPad->SetGridx();
1377 gPad->SetGridy();
1378 fhTimeMeanPulser->Draw("colz");
1379
1380 cPulser->cd(3);
1381 gPad->SetGridx();
1382 gPad->SetGridy();
1383 fhTimeRmsZoomFitPuls->Draw("colz");
1384
1385 cPulser->cd(4);
1386 gPad->SetGridx();
1387 gPad->SetGridy();
1388 fhTimeResFitPuls->Draw("colz");
1389 /*****************************/
1390
1392 TCanvas* cPulserEvo =
1393 new TCanvas("cPulserEvo", "Time difference evolution between 1st FEE of 1st GBTx of gDPB pairs", w, h);
1394 cPulserEvo->Divide(1, fuNrOfGdpbs - 1);
1395 for (UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs - 1; uGdpb++) {
1396 cPulserEvo->cd(1 + uGdpb);
1397 gPad->SetGridx();
1398 gPad->SetGridy();
1399 if (NULL != fvvhPulserTimeDiffEvoGdpbGdpb[uGdpb][uGdpb + 1])
1400 fvvhPulserTimeDiffEvoGdpbGdpb[uGdpb][uGdpb + 1]->Draw();
1401
1402 } // for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs - 1; uGdpb ++)
1403 /*****************************/
1405 for (UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; uGdpb++) {
1406 TCanvas* cPulserEvoGbtx = new TCanvas(Form("cPulserEvoGbtx%02u", uGdpb),
1407 Form("Time difference evolution between 1st FEE of GBTx "
1408 "pairs in gDPB %02u",
1409 uGdpb),
1410 w, h);
1411 cPulserEvoGbtx->Divide(1, kuNbGbtxPerGdpb - 1);
1412
1413 for (UInt_t uGbtx = 0; uGbtx < kuNbGbtxPerGdpb - 1; ++uGbtx) {
1414 cPulserEvoGbtx->cd(1 + uGbtx);
1415 gPad->SetGridx();
1416 gPad->SetGridy();
1417
1418 if (NULL != fvhPulserTimeDiffEvoGbtxGbtx[uGdpb * (kuNbGbtxPerGdpb - 1) + uGbtx])
1419 fvhPulserTimeDiffEvoGbtxGbtx[uGdpb * (kuNbGbtxPerGdpb - 1) + uGbtx]->Draw();
1420 } // for( UInt_t uGbtx = 0; uGbtx < kuNbGbtxPerGdpb - 1; ++uGbtx )
1421 } // for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; uGdpb ++)
1422 /*****************************/
1423 } // if( kTRUE == fbPulserModeEnable )
1424
1426 for (UInt_t uMod = 0; uMod < fuNrOfModules; uMod++) {
1427 TCanvas* cModRates =
1428 new TCanvas(Form("cModRate_m%02u", uMod), Form("Hit and error Rates for module %02u", uMod), w, h);
1429 cModRates->cd();
1430 gPad->SetLogy();
1431 fvhModRate[uMod]->Draw("hist");
1432
1433 fvhModErrorRate[uMod]->SetLineColor(kRed);
1434 fvhModErrorRate[uMod]->Draw("same hist");
1435 } // for( UInt_t uMod = 0; uMod < fuNrOfModules; uMod ++ )
1436 /*****************************/
1437
1439 // Try to recover canvas in case it was created already by another monitor
1440 // If not existing, create it
1441 fcMsSizeAll = dynamic_cast<TCanvas*>(gROOT->FindObject("cMsSizeAll"));
1442 if (NULL == fcMsSizeAll) {
1443 fcMsSizeAll = new TCanvas("cMsSizeAll", "Evolution of MS size in last 300 s", w, h);
1444 fcMsSizeAll->Divide(4, 3);
1445 LOG(info) << "Created MS size canvas in TOF monitor";
1446 } // if( NULL == fcMsSizeAll )
1447 else
1448 LOG(info) << "Recovered MS size canvas in TOF monitor";
1449
1450 LOG(info) << "Leaving CreateHistograms";
1451}
1452
1453Bool_t CbmStar2019MonitorTof::DoUnpack(const fles::Timeslice& ts, size_t component)
1454{
1456 LOG(info) << "Reset eTOF STAR histos ";
1458 bMcbmMoniTofResetHistos = kFALSE;
1459 } // if( bMcbmMoniTofResetHistos )
1461 LOG(info) << "Start saving eTOF STAR histos ";
1462 SaveAllHistos("data/histos_Shift_StarTof.root");
1463 bMcbmMoniTofSaveHistos = kFALSE;
1464 } // if( bSaveStsHistos )
1468 } // if (bMcbmMoniTofUpdateZoomedFit)
1471 bMcbmMoniTofRawDataPrint = kFALSE;
1472 } // if( bMcbmMoniTofRawDataPrint )
1476 } // if( bMcbmMoniTofPrintAllHitsEna )
1480 } // if( bMcbmMoniTofPrintAllEpochsEna )
1481
1483 std::chrono::time_point<std::chrono::system_clock> timeCurrent = std::chrono::system_clock::now();
1484 std::chrono::duration<double> elapsed_seconds = timeCurrent - fTimeLastHistoSaving;
1485 if (0 == fTimeLastHistoSaving.time_since_epoch().count()) {
1486 fTimeLastHistoSaving = timeCurrent;
1487 // fulNbBuiltSubEventLastPrintout = fulNbBuiltSubEvent;
1488 // fulNbStarSubEventLastPrintout = fulNbStarSubEvent;
1489 } // if( 0 == fTimeLastHistoSaving.time_since_epoch().count() )
1490 else if (300 < elapsed_seconds.count()) {
1491 std::time_t cTimeCurrent = std::chrono::system_clock::to_time_t(timeCurrent);
1492 char tempBuff[80];
1493 std::strftime(tempBuff, 80, "%F %T", localtime(&cTimeCurrent));
1494 /*
1495 LOG(info) << "CbmTofStarEventBuilder2018::DoUnpack => " << tempBuff
1496 << " Total number of Built events: " << std::setw(9) << fulNbBuiltSubEvent
1497 << ", " << std::setw(9) << (fulNbBuiltSubEvent - fulNbBuiltSubEventLastPrintout)
1498 << " events in last " << std::setw(4) << elapsed_seconds.count() << " s";
1499 fTimeLastPrintoutNbStarEvent = timeCurrent;
1500 fulNbBuiltSubEventLastPrintout = fulNbBuiltSubEvent;
1501*/
1502 fTimeLastHistoSaving = timeCurrent;
1503 SaveAllHistos("data/histos_shift.root");
1504 } // else if( 300 < elapsed_seconds.count() )
1505
1506 LOG(debug1) << "Timeslice contains " << ts.num_microslices(component) << "microslices.";
1507
1509 UInt_t uNbMsLoop = fuNbCoreMsPerTs;
1510 if (kFALSE == fbIgnoreOverlapMs) uNbMsLoop += fuNbOverMsPerTs;
1511
1512 Int_t messageType = -111;
1513 Double_t dTsStartTime = -1;
1514
1516 for (UInt_t uMsIdx = 0; uMsIdx < uNbMsLoop; uMsIdx++) {
1517 if (fuMsAcceptsPercent < uMsIdx) continue;
1518
1519 fuCurrentMs = uMsIdx;
1520
1521 if (0 == fulCurrentTsIndex && 0 == uMsIdx) {
1522 for (UInt_t uMsCompIdx = 0; uMsCompIdx < fvMsComponentsList.size(); ++uMsCompIdx) {
1523 UInt_t uMsComp = fvMsComponentsList[uMsCompIdx];
1524 auto msDescriptor = ts.descriptor(uMsComp, uMsIdx);
1525 LOG(info) << "hi hv eqid flag si sv idx/start crc size offset";
1526 LOG(info) << Form(
1527 "%02x %02x %04x %04x %02x %02x %016lx %08x %08x %016lx", static_cast<unsigned int>(msDescriptor.hdr_id),
1528 static_cast<unsigned int>(msDescriptor.hdr_ver), msDescriptor.eq_id, msDescriptor.flags,
1529 static_cast<unsigned int>(msDescriptor.sys_id), static_cast<unsigned int>(msDescriptor.sys_ver),
1530 msDescriptor.idx, msDescriptor.crc, msDescriptor.size, msDescriptor.offset);
1531 } // for( UInt_t uMsCompIdx = 0; uMsCompIdx < fvMsComponentsList.size(); ++uMsCompIdx )
1532 } // if( 0 == fulCurrentTsIndex && 0 == uMsIdx )
1533
1535 for (UInt_t uMsCompIdx = 0; uMsCompIdx < fvMsComponentsList.size(); ++uMsCompIdx) {
1536 constexpr uint32_t kuBytesPerMessage = 8;
1537
1538 UInt_t uMsComp = fvMsComponentsList[uMsCompIdx];
1539 auto msDescriptor = ts.descriptor(uMsComp, uMsIdx);
1540 fiEquipmentId = msDescriptor.eq_id;
1541 fdMsIndex = static_cast<double>(msDescriptor.idx);
1542 const uint8_t* msContent = reinterpret_cast<const uint8_t*>(ts.content(uMsComp, uMsIdx));
1543
1544 uint32_t size = msDescriptor.size;
1545 // fulLastMsIdx = msDescriptor.idx;
1546 if (size > 0) LOG(debug) << "Microslice: " << msDescriptor.idx << " has size: " << size;
1547 /*
1548 if( numCompMsInTs - fuOverlapMsNb <= m )
1549 {
1550// LOG(info) << "Ignore overlap Microslice: " << msDescriptor.idx;
1551 continue;
1552 } // if( numCompMsInTs - fuOverlapMsNb <= m )
1553*/
1554 if (0 == uMsIdx && 0 == uMsCompIdx) dTsStartTime = (1e-9) * fdMsIndex;
1555
1556 if (fdStartTimeMsSz < 0) fdStartTimeMsSz = (1e-9) * fdMsIndex;
1557 fvhMsSzPerLink[uMsComp]->Fill(size);
1558 if (2 * fuHistoryHistoSize < (1e-9) * fdMsIndex - fdStartTimeMsSz) {
1559 // Reset the evolution Histogram and the start time when we reach the end of the range
1560 fvhMsSzTimePerLink[uMsComp]->Reset();
1561 fdStartTimeMsSz = (1e-9) * fdMsIndex;
1562 } // if( 2 * fuHistoryHistoSize < (1e-9) * fdMsIndex - fdStartTimeMsSz )
1563 fvhMsSzTimePerLink[uMsComp]->Fill((1e-9) * fdMsIndex - fdStartTimeMsSz, size);
1564
1565 // If not integer number of message in input buffer, print warning/error
1566 if (0 != (size % kuBytesPerMessage))
1567 LOG(error) << "The input microslice buffer does NOT "
1568 << "contain only complete nDPB messages!";
1569
1570 // Compute the number of complete messages in the input microslice buffer
1571 uint32_t uNbMessages = (size - (size % kuBytesPerMessage)) / kuBytesPerMessage;
1572
1573 // Get the gDPB ID from the MS header
1575
1577 auto it = fGdpbIdIndexMap.find(fuGdpbId);
1578 if (it == fGdpbIdIndexMap.end()) {
1579 LOG(info) << "---------------------------------------------------------------";
1580 LOG(info) << "hi hv eqid flag si sv idx/start crc size offset";
1581 LOG(info) << Form(
1582 "%02x %02x %04x %04x %02x %02x %016lx %08x %08x %016lx", static_cast<unsigned int>(msDescriptor.hdr_id),
1583 static_cast<unsigned int>(msDescriptor.hdr_ver), msDescriptor.eq_id, msDescriptor.flags,
1584 static_cast<unsigned int>(msDescriptor.sys_id), static_cast<unsigned int>(msDescriptor.sys_ver),
1585 msDescriptor.idx, msDescriptor.crc, msDescriptor.size, msDescriptor.offset);
1586 LOG(fatal) << "Could not find the gDPB index for AFCK id 0x" << std::hex << fuGdpbId << std::dec
1587 << " in timeslice " << fulCurrentTsIndex << " in microslice " << fdMsIndex << " component "
1588 << uMsCompIdx << "\n"
1589 << "If valid this index has to be added in the TOF "
1590 "parameter file in the RocIdArray field";
1591 continue;
1592 } // if( it == fGdpbIdIndexMap.end() )
1593 else
1595
1596 // Prepare variables for the loop on contents
1597 const uint64_t* pInBuff = reinterpret_cast<const uint64_t*>(msContent);
1598 for (uint32_t uIdx = 0; uIdx < uNbMessages; uIdx++) {
1599 // Fill message
1600 uint64_t ulData = static_cast<uint64_t>(pInBuff[uIdx]);
1601
1603 if (0 == uIdx && kFALSE == fbOldFwData) {
1604 ProcessEpochCycle(ulData);
1605 continue;
1606 } // if( 0 == uIdx && kFALSE == fbOldFwData )
1607
1608 gdpbv100::Message mess(ulData);
1609
1610 if (fuRawDataPrintMsgIdx < fuRawDataPrintMsgNb || fair::Logger::Logging(fair::Severity::debug2)) {
1611 mess.printDataCout();
1613 }
1614
1615 // Increment counter for different message types
1616 // and fill the corresponding histogram
1617 messageType = mess.getMessageType();
1618 fviMsgCounter[messageType]++;
1619 fhMessType->Fill(messageType);
1620 fhGdpbMessType->Fill(messageType, fuGdpbNr);
1621
1625 if (fuGdpbNr == fuDiamondDpbIdx) fuGet4Id = mess.getGdpbGenChipId();
1627
1629 LOG(warning) << "Message with Get4 ID too high: " << fuGet4Id << " VS " << fuNrOfGet4PerGdpb
1630 << " set in parameters.";
1631
1632 switch (messageType) {
1633 case gdpbv100::MSG_HIT: {
1634 if (mess.getGdpbHitIs24b()) {
1635 fhGet4MessType->Fill(fuGet4Nr, 4);
1637 PrintGenInfo(mess);
1638 } // if( getGdpbHitIs24b() )
1639 else {
1640 fhGet4MessType->Fill(fuGet4Nr, 0);
1642 fvmEpSupprBuffer[fuGet4Nr].push_back(mess);
1643 } // else of if( getGdpbHitIs24b() )
1644 break;
1645 } // case gdpbv100::MSG_HIT:
1646 case gdpbv100::MSG_EPOCH: {
1648 if (1 == mess.getGdpbEpSync()) {
1649 fhGdpbEpochFlags->Fill(fuGdpbNr, 0);
1651 } // if (1 == mess.getGdpbEpSync())
1652
1653 if (1 == mess.getGdpbEpDataLoss()) fhGdpbEpochFlags->Fill(fuGdpbNr, 1);
1654
1655 if (1 == mess.getGdpbEpEpochLoss()) fhGdpbEpochFlags->Fill(fuGdpbNr, 2);
1656
1657 if (1 == mess.getGdpbEpMissmatch()) {
1658 fhGdpbEpochFlags->Fill(fuGdpbNr, 3);
1660 } // if (1 == mess.getGdpbEpMissmatch())
1661
1662 for (uint32_t uGet4Index = 0; uGet4Index < fuNrOfGet4PerGdpb; uGet4Index++) {
1663 fuGet4Id = uGet4Index;
1665 gdpbv100::Message tmpMess(mess);
1666 tmpMess.setGdpbGenChipId(uGet4Index);
1667
1668 fhGet4MessType->Fill(fuGet4Nr, 1);
1670 FillEpochInfo(tmpMess);
1671 } // for( uint32_t uGet4Index = 0; uGet4Index < fuNrOfGet4PerGdpb; uGetIndex ++ )
1672
1673 if (kTRUE == fbPrintAllEpochsEnable) {
1674 LOG(info) << "Epoch: " << Form("0x%08x ", fuGdpbId) << ", Merg"
1675 << ", Link " << std::setw(1) << mess.getGdpbEpLinkId() << ", epoch " << std::setw(8)
1676 << mess.getGdpbEpEpochNb() << ", Sync " << std::setw(1) << mess.getGdpbEpSync()
1677 << ", Data loss " << std::setw(1) << mess.getGdpbEpDataLoss() << ", Epoch loss "
1678 << std::setw(1) << mess.getGdpbEpEpochLoss() << ", Epoch miss " << std::setw(1)
1679 << mess.getGdpbEpMissmatch();
1680 } // if( kTRUE == fbPrintAllEpochsEnable )
1681 } // if this epoch message is a merged one valid for all chips
1682 else {
1683 fhGet4MessType->Fill(fuGet4Nr, 1);
1686 // FillEpochInfo(mess);
1687
1688 if (kTRUE == fbPrintAllEpochsEnable) {
1689 LOG(info) << "Epoch: " << Form("0x%08x ", fuGdpbId) << ", " << std::setw(4) << fuGet4Nr << ", Link "
1690 << std::setw(1) << mess.getGdpbEpLinkId() << ", epoch " << std::setw(8)
1691 << mess.getGdpbEpEpochNb() << ", Sync " << std::setw(1) << mess.getGdpbEpSync()
1692 << ", Data loss " << std::setw(1) << mess.getGdpbEpDataLoss() << ", Epoch loss "
1693 << std::setw(1) << mess.getGdpbEpEpochLoss() << ", Epoch miss " << std::setw(1)
1694 << mess.getGdpbEpMissmatch();
1695 } // if( kTRUE == fbPrintAllEpochsEnable )
1696 } // if single chip epoch message
1697 break;
1698 } // case gdpbv100::MSG_EPOCH:
1699 case gdpbv100::MSG_SLOWC: {
1700 fhGet4MessType->Fill(fuGet4Nr, 2);
1702 PrintSlcInfo(mess);
1703 break;
1704 } // case gdpbv100::MSG_SLOWC:
1705 case gdpbv100::MSG_SYST: {
1706 fhSysMessType->Fill(mess.getGdpbSysSubType());
1709 fhGet4MessType->Fill(fuGet4Nr, 3);
1711
1712 UInt_t uFeeNr = (fuGet4Id / fuNrOfGet4PerFee);
1713 if (0 <= fdStartTime) {
1717 1e-9 * (mess.getMsgFullTimeD(fvulCurrentEpoch[fuGet4Nr]) - fdStartTime), 1, 1);
1718
1719 UInt_t uGbtxNr = (uFeeNr / kuNbFeePerGbtx);
1720 UInt_t uGbtxNrInSys = fuGdpbNr * kuNbGbtxPerGdpb + uGbtxNr;
1721 fvhModErrorRate[fviModuleId[uGbtxNrInSys]]->Fill(
1723 fvhModErrorRatio[fviModuleId[uGbtxNrInSys]]->Fill(
1724 1e-9 * (mess.getMsgFullTimeD(fvulCurrentEpoch[fuGet4Nr]) - fdStartTime), 1, 1);
1725 } // if (0 <= fdStartTime)
1726 if (0 <= fdStartTimeLong) {
1728 1e-9 / 60.0 * (mess.getMsgFullTimeD(fvulCurrentEpoch[fuGet4Nr]) - fdStartTimeLong), 1 / 60.0);
1730 1e-9 / 60.0 * (mess.getMsgFullTimeD(fvulCurrentEpoch[fuGet4Nr]) - fdStartTimeLong), 1, 1 / 60.0);
1731 } // if (0 <= fdStartTime)
1732
1735 switch (mess.getGdpbSysErrData()) {
1737 fhGet4ChanErrors->Fill(dFullChId, 0);
1738 fvhGdpbGet4ChanErrors[fuGdpbNr]->Fill(dGdpbChId, 0);
1739 break;
1741 fhGet4ChanErrors->Fill(dFullChId, 1);
1742 fvhGdpbGet4ChanErrors[fuGdpbNr]->Fill(dGdpbChId, 1);
1743 break;
1745 fhGet4ChanErrors->Fill(dFullChId, 2);
1746 fvhGdpbGet4ChanErrors[fuGdpbNr]->Fill(dGdpbChId, 2);
1747 break;
1749 fhGet4ChanErrors->Fill(dFullChId, 3);
1750 fvhGdpbGet4ChanErrors[fuGdpbNr]->Fill(dGdpbChId, 3);
1751 break;
1753 fhGet4ChanErrors->Fill(dFullChId, 4);
1754 fvhGdpbGet4ChanErrors[fuGdpbNr]->Fill(dGdpbChId, 4);
1755 break;
1757 fhGet4ChanErrors->Fill(dFullChId, 5);
1758 fvhGdpbGet4ChanErrors[fuGdpbNr]->Fill(dGdpbChId, 5);
1759 break;
1761 fhGet4ChanErrors->Fill(dFullChId, 6);
1762 fvhGdpbGet4ChanErrors[fuGdpbNr]->Fill(dGdpbChId, 6);
1763 break;
1765 fhGet4ChanErrors->Fill(dFullChId, 7);
1766 fvhGdpbGet4ChanErrors[fuGdpbNr]->Fill(dGdpbChId, 7);
1767 break;
1769 fhGet4ChanErrors->Fill(dFullChId, 8);
1770 fvhGdpbGet4ChanErrors[fuGdpbNr]->Fill(dGdpbChId, 8);
1771 break;
1773 fhGet4ChanErrors->Fill(dFullChId, 9);
1774 fvhGdpbGet4ChanErrors[fuGdpbNr]->Fill(dGdpbChId, 9);
1775 break;
1777 fhGet4ChanErrors->Fill(dFullChId, 10);
1778 fvhGdpbGet4ChanErrors[fuGdpbNr]->Fill(dGdpbChId, 10);
1779 break;
1781 fhGet4ChanErrors->Fill(dFullChId, 11);
1782 fvhGdpbGet4ChanErrors[fuGdpbNr]->Fill(dGdpbChId, 11);
1783 break;
1785 fhGet4ChanErrors->Fill(dFullChId, 12);
1786 fvhGdpbGet4ChanErrors[fuGdpbNr]->Fill(dGdpbChId, 12);
1787 break;
1789 fhGet4ChanErrors->Fill(dFullChId, 13);
1790 fvhGdpbGet4ChanErrors[fuGdpbNr]->Fill(dGdpbChId, 13);
1791 break;
1793 fhGet4ChanErrors->Fill(dFullChId, 14);
1794 fvhGdpbGet4ChanErrors[fuGdpbNr]->Fill(dGdpbChId, 14);
1795 break;
1797 fhGet4ChanErrors->Fill(dFullChId, 15);
1798 fvhGdpbGet4ChanErrors[fuGdpbNr]->Fill(dGdpbChId, 15);
1799 break;
1801 fhGet4ChanErrors->Fill(dFullChId, 16);
1802 fvhGdpbGet4ChanErrors[fuGdpbNr]->Fill(dGdpbChId, 16);
1803 break;
1805 fhGet4ChanErrors->Fill(dFullChId, 17);
1806 fvhGdpbGet4ChanErrors[fuGdpbNr]->Fill(dGdpbChId, 17);
1807 break;
1809 fhGet4ChanErrors->Fill(dFullChId, 18);
1810 fvhGdpbGet4ChanErrors[fuGdpbNr]->Fill(dGdpbChId, 18);
1811 break;
1813 fhGet4ChanErrors->Fill(dFullChId, 19);
1814 fvhGdpbGet4ChanErrors[fuGdpbNr]->Fill(dGdpbChId, 19);
1815 break;
1817 fhGet4ChanErrors->Fill(dFullChId, 20);
1818 fvhGdpbGet4ChanErrors[fuGdpbNr]->Fill(dGdpbChId, 20);
1819 break;
1820 default: // Corrupt error or not yet supported error
1821 fhGet4ChanErrors->Fill(dFullChId, 21);
1822 fvhGdpbGet4ChanErrors[fuGdpbNr]->Fill(dGdpbChId, 21);
1823 break;
1824 } // Switch( mess.getGdpbSysErrData() )
1825 } // if( gdpbv100::SYSMSG_GET4_EVENT == mess.getGdpbSysSubType() )
1828 FillPattInfo(mess);
1829 } // if( gdpbv100::SYS_PATTERN == mess.getGdpbSysSubType() )
1830 PrintSysInfo(mess);
1831 break;
1832 } // case gdpbv100::MSG_SYST:
1837 // fhGet4MessType->Fill(fuGet4Nr, 5);
1838 FillStarTrigInfo(mess);
1839 break;
1840 default:
1841 LOG(error) << "Message type " << std::hex << std::setw(2) << static_cast<uint16_t>(messageType)
1842 << " not included in Get4 unpacker.";
1843 } // switch( mess.getMessageType() )
1844 } // for (uint32_t uIdx = 0; uIdx < uNbMessages; uIdx ++)
1845 } // for( UInt_t uMsCompIdx = 0; uMsCompIdx < fvMsComponentsList.size(); ++uMsCompIdx )
1846
1848 if (kTRUE == fbPulserModeEnable) {
1850 for (UInt_t uFeeA = 0; uFeeA < fuNrOfFeePerGdpb * fuNrOfGdpbs; ++uFeeA) {
1851 if (1 != fvuFeeNbHitsLastMs[uFeeA]) {
1853 fvuFeeNbHitsLastMs[uFeeA] = 0;
1854 continue;
1855 } // if( 1 != fvuFeeNbHitsLastMs[ uFeeA ] )
1856
1857 UInt_t uGdpbNr = uFeeA / fuNrOfFeePerGdpb;
1858 UInt_t uGbtxNr = uFeeA / kuNbFeePerGbtx;
1860 for (UInt_t uFeeB = uFeeA + 1; uFeeB < fuNrOfFeePerGdpb * fuNrOfGdpbs; ++uFeeB) {
1861 if (1 != fvuFeeNbHitsLastMs[uFeeB]) continue;
1862
1863 if (NULL != fvhTimeDiffPulser[uFeeA][uFeeB]) {
1864 Double_t dTimeDiff = 1e3 * (fdTsLastPulserHit[uFeeB] - fdTsLastPulserHit[uFeeA]);
1865 if (TMath::Abs(dTimeDiff) < kdMaxDtPulserPs) {
1866 fvhTimeDiffPulser[uFeeA][uFeeB]->Fill(dTimeDiff);
1867
1869 if (0 == uFeeA % kuNbFeePerGbtx && 0 == uFeeB % kuNbFeePerGbtx) {
1871 if (uGdpbNr == uFeeB / fuNrOfFeePerGdpb) {
1872 if (0 == uGbtxNr) {
1873 UInt_t uPlotIdx =
1874 uGdpbNr * (kuNbGbtxPerGdpb - 1) + (uFeeB - uGdpbNr * fuNrOfFeePerGdpb) / kuNbFeePerGbtx - 1;
1875 fvhPulserTimeDiffEvoGbtxGbtx[uPlotIdx]->Fill(1e-9 / 60.0 * (fdTsLastPulserHit[uFeeA] - fdStartTime),
1876 dTimeDiff);
1877 } // if( 0 == uGbtxNr )
1878 } // if( uGdpbNr == uFeeB / fuNrOfFeePerGdpb )
1879 else // if( NULL != fvvhPulserTimeDiffEvoGdpbGdpb[ uFeeB / fuNrOfFeePerGdpb ][ uGdpbNr ] )
1880 {
1882 if (0 == uGbtxNr && 0 == uFeeB / kuNbFeePerGbtx)
1883 fvvhPulserTimeDiffEvoGdpbGdpb[uGdpbNr][uFeeB / fuNrOfFeePerGdpb]->Fill(
1884 1e-9 / 60.0 * (fdTsLastPulserHit[uFeeA] - fdStartTime), dTimeDiff);
1885 } // else of if( uGdpbNr == uFeeB / fuNrOfFeePerGdpb )
1886 } // if( 0 == uFeeA % kuNbFeePerGbtx && 0 == uFeeB % kuNbFeePerGbtx )
1887 } // if( TMath::Abs( dTimeDiff ) < kdMaxDtPulserPs )
1888 // else if( 10 == uFeeA && 20 == uFeeB )
1889 // LOG(info) << "new in 10 " << dTimeDiff;
1890 } // if( NULL != fvhTimeDiffPulser[uFeeA][uFeeB] )
1891 } // for( UInt_t uFeeB = uFeeA + 1; uFeeB < fuNrOfFeePerGdpb * fuNrOfGdpbs; ++uFeeB)
1892
1894 fvuFeeNbHitsLastMs[uFeeA] = 0;
1895 } // for( UInt_t uFeeA = 0; uFeeA < fuNrOfFeePerGdpb * fuNrOfGdpbs; ++uFeeA)
1896 } // if( kTRUE == fbPulserModeEnable )
1897
1898 if (kTRUE == fbCoincMapsEnable) {
1899 for (UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb) {
1900 for (UInt_t uChanA = 0; uChanA < fuNrOfChannelsPerGdpb; ++uChanA) {
1901 if (0 == fvuCoincNbHitsLastMs[uGdpb][uChanA]) {
1903 continue;
1904 } // if( 0 == fvuCoincNbHitsLastMs[ uGdpb ][ uChanA ] )
1905
1906 for (UInt_t uChanB = uChanA + 1; uChanB < fuNrOfChannelsPerGdpb; ++uChanB) {
1907 if (0 == fvuCoincNbHitsLastMs[uGdpb][uChanB]) {
1909 continue;
1910 } // if( 0 == fvuCoincNbHitsLastMs[ uGdpb ][ uChanB ] )
1911
1912 Double_t dTimeDiff = 1e3 * (fvdCoincTsLastHit[uGdpb][uChanB] - fvdCoincTsLastHit[uGdpb][uChanA]);
1913
1914 fvhCoincMeanAllChanGdpb[uGdpb]->Fill(uChanA, uChanB, dTimeDiff);
1915
1916 if (TMath::Abs(dTimeDiff) < kdMaxDtPulserPs) {
1917 fvhCoincMapAllChanGdpb[uGdpb]->Fill(uChanA, uChanB);
1918 } // if( TMath::Abs( dTimeDiff ) < kdMaxDtPulserPs )
1919 } // for( UInt_t uChanA = 0; uChanA < fuNrOfChannelsPerGdpb; ++uChan A )
1920
1922 fvuCoincNbHitsLastMs[uGdpb][uChanA] = 0;
1923 } // for( UInt_t uChanA = 0; uChanA < fuNrOfChannelsPerGdpb; ++uChan A )
1924 } // for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
1925 } // if( kTRUE == fbCoincMapsEnable )
1926 } // for( UInt_t uMsIdx = 0; uMsIdx < uNbMsLoop; uMsIdx ++ )
1927
1928 // Update RMS plots only every 10s in data
1929 if (kTRUE == fbPulserModeEnable) {
1930 if (10.0 < dTsStartTime - fdLastRmsUpdateTime) {
1931 // Reset summary histograms for safety
1932 fhTimeMeanPulser->Reset();
1933 fhTimeRmsPulser->Reset();
1934
1935 for (UInt_t uFeeA = 0; uFeeA < fuNrOfFeePerGdpb * fuNrOfGdpbs; uFeeA++)
1936 for (UInt_t uFeeB = 0; uFeeB < fuNrOfFeePerGdpb * fuNrOfGdpbs; uFeeB++)
1937 if (NULL != fvhTimeDiffPulser[uFeeA][uFeeB]) {
1938 fhTimeMeanPulser->Fill(uFeeA, uFeeB, fvhTimeDiffPulser[uFeeA][uFeeB]->GetMean());
1939 fhTimeRmsPulser->Fill(uFeeA, uFeeB, fvhTimeDiffPulser[uFeeA][uFeeB]->GetRMS());
1940 } // for( UInt_t uChan = 0; uChan < kuNbChanTest - 1; uChan++)
1942 } // if( 10.0 < dTsStartTime - fdLastRmsUpdateTime )
1943 } // if( kTRUE == fbPulserModeEnable )
1944
1946
1947 return kTRUE;
1948}
1949
1951{
1952 uint64_t ulEpochCycleVal = ulCycleData & gdpbv100::kulEpochCycleFieldSz;
1953
1954 if (fuRawDataPrintMsgIdx < fuRawDataPrintMsgNb || fair::Logger::Logging(fair::Severity::debug2)) {
1955 LOG(info) << "CbmMcbm2018MonitorTof::ProcessEpochCyle => "
1956 << Form(" TS %5lu MS %3lu In data 0x%016lX Cycle 0x%016lX", fulCurrentTsIndex, fuCurrentMs, ulCycleData,
1957 ulEpochCycleVal);
1959 }
1960
1961 for (uint32_t uGet4Index = 0; uGet4Index < fuNrOfGet4PerGdpb; uGet4Index++) {
1962 fuGet4Id = uGet4Index;
1964 /*
1965 if( !( ulEpochCycleVal == fvulCurrentEpochCycle[fuGet4Nr] ||
1966 ulEpochCycleVal == fvulCurrentEpochCycle[fuGet4Nr] + 1 ) )
1967 LOG(error) << "CbmMcbm2018MonitorTof::ProcessEpochCyle => "
1968 << " Missmatch in epoch cycles detected, probably fake cycles due to epoch index corruption! "
1969 << Form( " Current cycle 0x%09X New cycle 0x%09X", fvulCurrentEpochCycle[fuGet4Nr], ulEpochCycleVal );
1970*/
1971 fvulCurrentEpochCycle[fuGet4Nr] = ulEpochCycleVal;
1972 } // for( uint32_t uGet4Index = 0; uGet4Index < fuNrOfGet4PerGdpb; uGet4Index ++ )
1973 return;
1974}
1975
1977{
1978 UInt_t uChannel = mess.getGdpbHitChanId();
1979 UInt_t uTot = mess.getGdpbHit32Tot();
1980 UInt_t uFts = mess.getGdpbHitFineTs();
1981
1982 ULong64_t ulCurEpochGdpbGet4 = fvulCurrentEpoch[fuGet4Nr];
1983
1984 // In Ep. Suppr. Mode, receive following epoch instead of previous
1985 if (0 < ulCurEpochGdpbGet4) ulCurEpochGdpbGet4--;
1986 else
1987 ulCurEpochGdpbGet4 = gdpbv100::kuEpochCounterSz; // Catch epoch cycle!
1988
1989 UInt_t uChannelNr = fuGet4Id * fuNrOfChannelsPerGet4 + uChannel;
1990 UInt_t uChannelNrInFee = (fuGet4Id % fuNrOfGet4PerFee) * fuNrOfChannelsPerGet4 + uChannel;
1991 UInt_t uFeeNr = (fuGet4Id / fuNrOfGet4PerFee);
1992 UInt_t uFeeNrInSys = fuGdpbNr * fuNrOfFeePerGdpb + uFeeNr;
1993 UInt_t uRemappedChannelNr = uFeeNr * fuNrOfChannelsPerFee + fvuGet4ToPadi[uChannelNrInFee];
1994 UInt_t uGbtxNr = (uFeeNr / kuNbFeePerGbtx);
1995 UInt_t uFeeInGbtx = (uFeeNr % kuNbFeePerGbtx);
1996 UInt_t uGbtxNrInSys = fuGdpbNr * kuNbGbtxPerGdpb + uGbtxNr;
1997
1998 ULong_t ulHitTime = mess.getMsgFullTime(ulCurEpochGdpbGet4);
1999 Double_t dHitTime = mess.getMsgFullTimeD(ulCurEpochGdpbGet4);
2000
2001 // In 32b mode the coarse counter is already computed back to 112 FTS bins
2002 // => need to hide its contribution from the Finetime
2003 // => FTS = Fullt TS modulo 112
2004 uFts = mess.getGdpbHitFullTs() % 112;
2005
2006 fvhChCount_gDPB[fuGdpbNr]->Fill(uChannelNr);
2007 fvhRawFt_gDPB[fuGdpbNr]->Fill(uChannelNr, uFts);
2008 fvhRawTot_gDPB[fuGdpbNr]->Fill(uChannelNr, uTot);
2009
2011 fvhRemapChCount_gDPB[fuGdpbNr]->Fill(uRemappedChannelNr);
2012 fvhRemapTot_gDPB[fuGdpbNr]->Fill(uRemappedChannelNr, uTot);
2013
2014 if (uGbtxNrInSys < fuNrOfGbtx) {
2015 UInt_t uOffset = uGbtxNrInSys * kuNbFeeSide * fuNrOfChannelsPerFee;
2016 if (fviRpcSide[uGbtxNrInSys])
2017 fvhRemapTotSideB_mod[fviModuleId[uGbtxNrInSys]]->Fill(uRemappedChannelNr - uOffset, uTot);
2018 else
2019 fvhRemapTotSideA_mod[fviModuleId[uGbtxNrInSys]]->Fill(uRemappedChannelNr - uOffset, uTot);
2020 } // if( uGbtxNrInSys < fuNrOfGbtx )
2021 /*
2022 switch( ( uRemappedChannelNr / fuNrOfChannelsPerFee ) / kuNbFeeSide )
2023 {
2024 case 0: // Module 1 Side A
2025 fvhRemapTotSideA_mod[ fuGdpbNr ]->Fill( uRemappedChannelNr , uTot);
2026 case 1: // Module 1 Side B
2027 fvhRemapTotSideB_mod[ fuGdpbNr ]->Fill( uRemappedChannelNr - 160 , uTot);
2028 case 2: // Module 2 Side A
2029 fvhRemapTotSideA_mod[ fuGdpbNr ]->Fill( uRemappedChannelNr - 160 , uTot);
2030 case 3: // Module 2 Side B
2031 fvhRemapTotSideB_mod[ fuGdpbNr ]->Fill( uRemappedChannelNr - 320, uTot);
2032 } // switch( ( uRemappedChannelNr / fuNrOfChannelsPerFee ) % kuNbFeeSide )
2033*/
2034 // In Run rate evolution
2035 if (fdStartTime < 0) fdStartTime = dHitTime;
2036
2037 // Reset the evolution Histogram and the start time when we reach the end of the range
2038 if (fuHistoryHistoSize < 1e-9 * (dHitTime - fdStartTime)) {
2040 fdStartTime = dHitTime;
2041 } // if( fuHistoryHistoSize < 1e-9 * (dHitTime - fdStartTime) )
2042
2043 // In Run rate evolution
2044 if (fdStartTimeLong < 0) fdStartTimeLong = dHitTime;
2045
2046 // Reset the evolution Histogram and the start time when we reach the end of the range
2047 if (fuHistoryHistoSizeLong < 1e-9 * (dHitTime - fdStartTimeLong) / 60.0) {
2049 fdStartTimeLong = dHitTime;
2050 } // if( fuHistoryHistoSize < 1e-9 * (dHitTime - fdStartTime) / 60.0 )
2051
2053 if (kTRUE == fbPulserModeEnable) {
2056 if (gdpbv100::kuFeePulserChannel == uChannelNrInFee) {
2057 fdTsLastPulserHit[uFeeNrInSys] = dHitTime;
2058 fvuFeeNbHitsLastMs[uFeeNrInSys]++;
2059 /*
2061 for( UInt_t uFeeB = 0; uFeeB < uFeeNrInSys; uFeeB++)
2062 if( NULL != fvhTimeDiffPulser[uFeeB][uFeeNrInSys] )
2063 {
2064 Double_t dTimeDiff = 1e3 * ( fdTsLastPulserHit[ uFeeNrInSys ] - fdTsLastPulserHit[ uFeeB ] );
2065 if( TMath::Abs( dTimeDiff ) < kdMaxDtPulserPs )
2066 {
2067 fvhTimeDiffPulser[uFeeB][uFeeNrInSys]->Fill( dTimeDiff );
2068
2070 if( 0 == uFeeInGbtx && 0 == uFeeB % kuNbFeePerGbtx )
2071 {
2073 if( fuGdpbNr == uFeeB / fuNrOfFeePerGdpb )
2074 {
2075 if( 0 == uFeeB / kuNbFeePerGbtx )
2076 {
2077 UInt_t uPlotIdx = fuGdpbNr * ( kuNbGbtxPerGdpb - 1) + uGbtxNr - 1;
2078 fvhPulserTimeDiffEvoGbtxGbtx[ uPlotIdx ]->Fill( 1e-9 / 60.0 * (dHitTime - fdStartTime), dTimeDiff );
2079 }
2080 } // if( fuGdpbNr == uFeeB / fuNrOfFeePerGdpb )
2081 else // if( NULL != fvvhPulserTimeDiffEvoGdpbGdpb[ uFeeB / fuNrOfFeePerGdpb ][ fuGdpbNr ] )
2082 {
2084 if( 0 == uGbtxNr && 0 == uFeeB / kuNbFeePerGbtx )
2085 fvvhPulserTimeDiffEvoGdpbGdpb[ uFeeB / fuNrOfFeePerGdpb ][ fuGdpbNr ]->Fill(
2086 1e-9 / 60.0 * (dHitTime - fdStartTime), dTimeDiff );
2087 } // else of if( fuGdpbNr == uFeeB / fuNrOfFeePerGdpb )
2088 } // if( 0 == uFeeInGbtx && 0 == uFeeB % kuNbFeePerGbtx )
2089 } // if( TMath::Abs( dTimeDiff ) < kdMaxDtPulserPs )
2090 // else if( 10 == uFeeB && 20 == uFeeNrInSys )
2091 // LOG(info) << "new in 20 " << dTimeDiff;
2092 } // if( NULL != fvhTimeDiffPulser[uFeeB][uFeeB] )
2093
2095 for( UInt_t uFeeB = uFeeNrInSys + 1; uFeeB < fuNrOfFeePerGdpb * fuNrOfGdpbs; uFeeB++)
2096 if( NULL != fvhTimeDiffPulser[uFeeNrInSys][uFeeB] )
2097 {
2098 Double_t dTimeDiff = 1e3 * ( fdTsLastPulserHit[ uFeeB ] - fdTsLastPulserHit[ uFeeNrInSys ] );
2099 if( TMath::Abs( dTimeDiff ) < kdMaxDtPulserPs )
2100 {
2101 fvhTimeDiffPulser[uFeeNrInSys][uFeeB]->Fill( dTimeDiff );
2102
2104 if( 0 == uFeeInGbtx && 0 == uFeeB % kuNbFeePerGbtx )
2105 {
2107 if( fuGdpbNr == uFeeB / fuNrOfFeePerGdpb )
2108 {
2109 if( 0 == uGbtxNr )
2110 {
2111 UInt_t uPlotIdx = fuGdpbNr * ( kuNbGbtxPerGdpb - 1)
2112 + ( uFeeB - fuGdpbNr * fuNrOfFeePerGdpb) / kuNbFeePerGbtx - 1;
2113 fvhPulserTimeDiffEvoGbtxGbtx[ uPlotIdx ]->Fill( 1e-9 / 60.0 * (dHitTime - fdStartTime), dTimeDiff );
2114 } // if( 0 == uGbtxNr )
2115 } // if( fuGdpbNr == uFeeB / fuNrOfFeePerGdpb )
2116 else // if( NULL != fvvhPulserTimeDiffEvoGdpbGdpb[ uFeeB / fuNrOfFeePerGdpb ][ fuGdpbNr ] )
2117 {
2119 if( 0 == uGbtxNr && 0 == uFeeB / kuNbFeePerGbtx )
2120 fvvhPulserTimeDiffEvoGdpbGdpb[ fuGdpbNr ][ uFeeB / fuNrOfFeePerGdpb ]->Fill(
2121 1e-9 / 60.0 * (dHitTime - fdStartTime), dTimeDiff );
2122 } // else of if( fuGdpbNr == uFeeB / fuNrOfFeePerGdpb )
2123 } // if( 0 == uFeeInGbtx && 0 == uFeeB % kuNbFeePerGbtx )
2124 } // if( TMath::Abs( dTimeDiff ) < kdMaxDtPulserPs )
2125 // else if( 10 == uFeeNrInSys && 20 == uFeeB )
2126 // LOG(info) << "new in 10 " << dTimeDiff;
2127 } // if( NULL != fvhTimeDiffPulser[uFeeNrInSys][uFeeB] )
2128*/
2129 } // if( gdpbv100::kuFeePulserChannel == uChannelNrInFee )
2131 else if (fuGdpbNr == fuDiamondDpbIdx && 0 == uChannelNrInFee) {
2132 fdTsLastPulserHit[uFeeNrInSys] = dHitTime;
2133 fvuFeeNbHitsLastMs[uFeeNrInSys]++;
2134 } // if( fuGdpbNr == fuDiamondDpbIdx && 0 == uChannelNrInFee )
2135 } // if( kTRUE == fbPulserModeEnable )
2136
2138 if (kTRUE == fbCoincMapsEnable) {
2139 fvdCoincTsLastHit[fuGdpbNr][uRemappedChannelNr] = dHitTime;
2140 fvuCoincNbHitsLastMs[fuGdpbNr][uRemappedChannelNr]++;
2141 } // if( kTRUE == fbCoincMapsEnable )
2142
2143 if (0 <= fdStartTime) {
2144 fvhChannelRate_gDPB[fuGdpbNr]->Fill(1e-9 * (dHitTime - fdStartTime), uChannelNr);
2145 fvhRemapChRate_gDPB[fuGdpbNr]->Fill(1e-9 * (dHitTime - fdStartTime), uRemappedChannelNr);
2146 fvhFeeRate_gDPB[(fuGdpbNr * fuNrOfFeePerGdpb) + uFeeNr]->Fill(1e-9 * (dHitTime - fdStartTime));
2147 fvhFeeErrorRatio_gDPB[(fuGdpbNr * fuNrOfFeePerGdpb) + uFeeNr]->Fill(1e-9 * (dHitTime - fdStartTime), 0, 1);
2148
2149 fvhModRate[fviModuleId[uGbtxNrInSys]]->Fill(1e-9 * (dHitTime - fdStartTime));
2150 fvhModErrorRatio[fviModuleId[uGbtxNrInSys]]->Fill(1e-9 * (dHitTime - fdStartTime), 0, 1);
2151 } // if (0 <= fdStartTime)
2152
2153 if (0 <= fdStartTimeLong) {
2154 fvhFeeRateLong_gDPB[(fuGdpbNr * fuNrOfFeePerGdpb) + uFeeNr]->Fill(1e-9 / 60.0 * (dHitTime - fdStartTimeLong),
2155 1 / 60.0);
2156 fvhFeeErrorRatioLong_gDPB[(fuGdpbNr * fuNrOfFeePerGdpb) + uFeeNr]->Fill(1e-9 / 60.0 * (dHitTime - fdStartTimeLong),
2157 0, 1 / 60.0);
2158 } // if (0 <= fdStartTimeLong)
2159
2160 if (kTRUE == fbPrintAllHitsEnable) {
2161 LOG(info) << "Hit: " << Form("0x%08x ", fuGdpbId) << ", " << std::setw(2) << fuGet4Nr << ", " << std::setw(3)
2162 << uChannel << ", " << std::setw(3) << uTot << ", epoch " << std::setw(3) << ulCurEpochGdpbGet4
2163 << ", FullTime Clk " << Form("%12lu ", ulHitTime) << ", FullTime s " << Form("%12.9f ", dHitTime / 1e9)
2164 << ", FineTime " << uFts;
2165 } // if( kTRUE == fbPrintAllHitsEnable )
2166}
2167
2169{
2170 ULong64_t ulEpochNr = mess.getGdpbEpEpochNb();
2171 /*
2172 if( fvulCurrentEpoch[fuGet4Nr] < ulEpochNr )
2173 fvulCurrentEpochCycle[fuGet4Nr]++;
2174*/
2175 fvulCurrentEpoch[fuGet4Nr] = ulEpochNr;
2177
2178 if (1 == mess.getGdpbEpSync()) fhGet4EpochFlags->Fill(fuGet4Nr, 0);
2179 if (1 == mess.getGdpbEpDataLoss()) fhGet4EpochFlags->Fill(fuGet4Nr, 1);
2180 if (1 == mess.getGdpbEpEpochLoss()) fhGet4EpochFlags->Fill(fuGet4Nr, 2);
2181 if (1 == mess.getGdpbEpMissmatch()) fhGet4EpochFlags->Fill(fuGet4Nr, 3);
2182
2183 fulCurrentEpochTime = mess.getMsgFullTime(ulEpochNr);
2184
2187 if (0 < ulEpochNr) mess.setGdpbEpEpochNb(ulEpochNr - 1);
2188 else
2190
2191 Int_t iBufferSize = fvmEpSupprBuffer[fuGet4Nr].size();
2192 if (0 < iBufferSize) {
2193 LOG(debug) << "Now processing stored messages for for get4 " << fuGet4Nr << " with epoch number "
2194 << (fvulCurrentEpoch[fuGet4Nr] - 1);
2195
2198 std::stable_sort(fvmEpSupprBuffer[fuGet4Nr].begin(), fvmEpSupprBuffer[fuGet4Nr].begin());
2199
2200 for (Int_t iMsgIdx = 0; iMsgIdx < iBufferSize; iMsgIdx++) {
2202 } // for( Int_t iMsgIdx = 0; iMsgIdx < iBufferSize; iMsgIdx++ )
2203
2204 fvmEpSupprBuffer[fuGet4Nr].clear();
2205 } // if( 0 < fvmEpSupprBuffer[fuGet4Nr] )
2206}
2207
2209{
2210 if (fGdpbIdIndexMap.end() != fGdpbIdIndexMap.find(fuGdpbId)) {
2211 UInt_t uChip = mess.getGdpbGenChipId();
2212 UInt_t uChan = mess.getGdpbSlcChan();
2213 UInt_t uEdge = mess.getGdpbSlcEdge();
2214 UInt_t uData = mess.getGdpbSlcData();
2215 UInt_t uType = mess.getGdpbSlcType();
2216 Double_t dGdpbChId = fuGet4Id * fuNrOfChannelsPerGet4 + mess.getGdpbSlcChan() + 0.5 * mess.getGdpbSlcEdge();
2217 Double_t dFullChId = fuGet4Nr * fuNrOfChannelsPerGet4 + mess.getGdpbSlcChan() + 0.5 * mess.getGdpbSlcEdge();
2218 Double_t dMessTime = static_cast<Double_t>(fulCurrentEpochTime) * 1.e-9;
2219
2220 switch (uType) {
2221 case 0: // Scaler counter
2222 {
2223 fhGet4ChanScm->Fill(dFullChId, uType);
2224 fvhGdpbGet4ChanScm[fuGdpbNr]->Fill(dGdpbChId, uType);
2225 fhScmScalerCounters->Fill(uData, dFullChId);
2226 break;
2227 }
2228 case 1: // Deadtime counter
2229 {
2230 fhGet4ChanScm->Fill(dFullChId, uType);
2231 fvhGdpbGet4ChanScm[fuGdpbNr]->Fill(dGdpbChId, uType);
2232 fhScmDeadtimeCounters->Fill(uData, dFullChId);
2233 break;
2234 }
2235 case 2: // SPI message
2236 {
2237 /*
2238 LOG(info) << "GET4 Slow Control message, epoch "
2239 << static_cast<Int_t>(fvulCurrentEpoch[fuGet4Nr]) << ", time " << std::setprecision(9)
2240 << std::fixed << dMessTime << " s "
2241 << " for board ID " << std::hex << std::setw(4) << fuGdpbId
2242 << std::dec << "\n"
2243 << " +++++++ > Chip = "
2244 << std::setw(2) << uChip << ", Chan = "
2245 << std::setw(1) << uChan << ", Edge = "
2246 << std::setw(1) << uEdge << ", Type = "
2247 << std::setw(1) << mess.getGdpbSlcType() << ", Data = 0x"
2248// << std::hex << std::setw(6) << mess.getGdpbSlcData()
2249 << Form( "%06x", uData )
2250 << std::dec << ", CRC = " << uCRC
2251// << " RAW: " << Form( "%08x", mess.getGdpbSlcMess() );
2252*/
2253 fhGet4ChanScm->Fill(dFullChId, uType);
2254 fvhGdpbGet4ChanScm[fuGdpbNr]->Fill(dGdpbChId, uType);
2255 break;
2256 }
2257 case 3: // Start message or SEU counter
2258 {
2259 if (0 == mess.getGdpbSlcChan() && 0 == mess.getGdpbSlcEdge()) // START message
2260 {
2261 /*
2262 LOG(info) << std::setprecision(9)
2263 << std::fixed << dMessTime << " s ";
2264 LOG(info) << "GET4 Slow Control message, epoch "
2265 << static_cast<Int_t>(fvulCurrentEpoch[fuGet4Nr]) << ", time " << std::setprecision(9)
2266 << std::fixed << dMessTime << " s "
2267 << " for board ID " << std::hex << std::setw(4) << fuGdpbId
2268 << std::dec << "\n"
2269 << " +++++++ > Chip = "
2270 << std::setw(2) << mess.getGdpbGenChipId() << ", Chan = "
2271 << std::setw(1) << mess.getGdpbSlcChan() << ", Edge = "
2272 << std::setw(1) << mess.getGdpbSlcEdge() << ", Type = "
2273 << std::setw(1) << mess.getGdpbSlcType() << ", Data = 0x"
2274 // << std::hex << std::setw(6) << mess.getGdpbSlcData()
2275 << Form( "%06x", mess.getGdpbSlcData() )
2276 << std::dec << ", CRC = " << mess.getGdpbSlcCrc();
2277*/
2278 fhGet4ChanScm->Fill(dFullChId, uType + 1);
2279 fvhGdpbGet4ChanScm[fuGdpbNr]->Fill(dGdpbChId, uType + 1);
2280 } // if( 0 == mess.getGdpbSlcChan() && 0 == mess.getGdpbSlcEdge() )
2281 else if (0 == mess.getGdpbSlcChan() && 1 == mess.getGdpbSlcEdge()) // SEU counter message
2282 {
2283 /*
2284 LOG(info) << "GET4 Slow Control message, epoch "
2285 << static_cast<Int_t>(fvulCurrentEpoch[fuGet4Nr]) << ", time " << std::setprecision(9)
2286 << std::fixed << dMessTime << " s "
2287 << " for board ID " << std::hex << std::setw(4) << fuGdpbId
2288 << std::dec << "\n"
2289 << " +++++++ > Chip = "
2290 << std::setw(2) << mess.getGdpbGenChipId() << ", Chan = "
2291 << std::setw(1) << mess.getGdpbSlcChan() << ", Edge = "
2292 << std::setw(1) << mess.getGdpbSlcEdge() << ", Type = "
2293 << std::setw(1) << mess.getGdpbSlcType() << ", Data = 0x"
2294// << std::hex << std::setw(6) << mess.getGdpbSlcData()
2295 << Form( "%06x", mess.getGdpbSlcData() )
2296 << std::dec << ", CRC = " << mess.getGdpbSlcCrc();
2297*/
2298 fhGet4ChanScm->Fill(dFullChId, uType);
2299 fvhGdpbGet4ChanScm[fuGdpbNr]->Fill(dGdpbChId, uType);
2300 fhScmSeuCounters->Fill(uData, dFullChId);
2301 fhScmSeuCountersEvo->Fill(dMessTime - fdStartTime * 1.e-9, uData, dFullChId);
2302 } // else if( 0 == mess.getGdpbSlcChan() && 1 == mess.getGdpbSlcEdge() )
2303 break;
2304 }
2305 default: // Should never happen
2306 break;
2307 } // switch( mess.getGdpbSlcType() )
2308 }
2309}
2310
2312{
2313 Int_t mType = mess.getMessageType();
2314 Int_t channel = mess.getGdpbHitChanId();
2315 uint64_t uData = mess.getData();
2316
2317 LOG(debug) << "Get4 MSG type " << mType << " from gdpbId " << fuGdpbId << ", getId " << fuGet4Id << ", (hit channel) "
2318 << channel << " data " << std::hex << uData;
2319}
2320
2322{
2323 if (fGdpbIdIndexMap.end() != fGdpbIdIndexMap.find(fuGdpbId))
2324 LOG(debug) << "GET4 System message, epoch " << static_cast<Int_t>(fvulCurrentEpoch[fuGet4Nr]) << ", time "
2325 << std::setprecision(9) << std::fixed << Double_t(fulCurrentEpochTime) * 1.e-9 << " s "
2326 << " for board ID " << std::hex << std::setw(4) << fuGdpbId << std::dec;
2327
2328 switch (mess.getGdpbSysSubType()) {
2330 uint32_t uData = mess.getGdpbSysErrData();
2334 LOG(debug) << " +++++++ > gDPB: " << std::hex << std::setw(4) << fuGdpbId << std::dec
2335 << ", Chip = " << std::setw(2) << mess.getGdpbGenChipId() << ", Chan = " << std::setw(1)
2336 << mess.getGdpbSysErrChanId() << ", Edge = " << std::setw(1) << mess.getGdpbSysErrEdge()
2337 << ", Empt = " << std::setw(1) << mess.getGdpbSysErrUnused() << ", Data = " << std::hex
2338 << std::setw(2) << uData << std::dec << " -- GET4 V1 Error Event";
2339 else
2340 LOG(debug) << " +++++++ >gDPB: " << std::hex << std::setw(4) << fuGdpbId << std::dec
2341 << ", Chip = " << std::setw(2) << mess.getGdpbGenChipId() << ", Chan = " << std::setw(1)
2342 << mess.getGdpbSysErrChanId() << ", Edge = " << std::setw(1) << mess.getGdpbSysErrEdge()
2343 << ", Empt = " << std::setw(1) << mess.getGdpbSysErrUnused() << ", Data = " << std::hex
2344 << std::setw(2) << uData << std::dec << " -- GET4 V1 Error Event ";
2345 break;
2346 } // case gdpbv100::SYSMSG_GET4_EVENT
2348 LOG(debug) << "Unknown GET4 message, data: " << std::hex << std::setw(8) << mess.getGdpbSysUnkwData() << std::dec
2349 << " Full message: " << std::hex << std::setw(16) << mess.getData() << std::dec;
2350 break;
2351 } // case gdpbv100::SYS_GDPB_UNKWN:
2353 LOG(debug) << "GET4 synchronization pulse missing";
2354 break;
2355 } // case gdpbv100::SYS_GET4_SYNC_MISS:
2356 case gdpbv100::SYS_PATTERN: {
2357 LOG(debug) << "ASIC pattern for missmatch, disable or resync";
2358 break;
2359 } // case gdpbv100::SYS_PATTERN:
2360 default: {
2361 LOG(debug) << "Crazy system message, subtype " << mess.getGdpbSysSubType();
2362 break;
2363 } // default
2364
2365 } // switch( getGdpbSysSubType() )
2366}
2367
2369{
2370 uint16_t usType = mess.getGdpbSysPattType();
2371 uint16_t usIndex = mess.getGdpbSysPattIndex();
2372 uint32_t uPattern = mess.getGdpbSysPattPattern();
2373
2374 switch (usType) {
2376 LOG(debug) << Form("Missmatch pattern message => Type %d, Index %2d, Pattern 0x%08X", usType, usIndex, uPattern);
2377 for (UInt_t uBit = 0; uBit < 32; ++uBit)
2378 if ((uPattern >> uBit) & 0x1) {
2379 UInt_t uBadAsic = ConvertElinkToGet4(32 * usIndex + uBit);
2381 if (fuGdpbNr == fuDiamondDpbIdx) uBadAsic = 32 * usIndex + uBit;
2382 fhPatternMissmatch->Fill(uBadAsic, fuGdpbNr);
2384 } // if( ( uPattern >> uBit ) & 0x1 )
2385
2386 break;
2387 } // case gdpbv100::PATT_MISSMATCH:
2388 case gdpbv100::PATT_ENABLE: {
2389 for (UInt_t uBit = 0; uBit < 32; ++uBit)
2390 if ((uPattern >> uBit) & 0x1) {
2391 UInt_t uBadAsic = ConvertElinkToGet4(32 * usIndex + uBit);
2393 if (fuGdpbNr == fuDiamondDpbIdx) uBadAsic = 32 * usIndex + uBit;
2394 fhPatternEnable->Fill(uBadAsic, fuGdpbNr);
2396 } // if( ( uPattern >> uBit ) & 0x1 )
2397
2398 break;
2399 } // case gdpbv100::PATT_ENABLE:
2400 case gdpbv100::PATT_RESYNC: {
2401 LOG(debug) << Form("RESYNC pattern message => Type %d, Index %2d, Pattern 0x%08X", usType, usIndex, uPattern);
2402
2403 for (UInt_t uBit = 0; uBit < 32; ++uBit)
2404 if ((uPattern >> uBit) & 0x1) {
2405 UInt_t uBadAsic = ConvertElinkToGet4(32 * usIndex + uBit);
2407 if (fuGdpbNr == fuDiamondDpbIdx) uBadAsic = 32 * usIndex + uBit;
2408 fhPatternResync->Fill(uBadAsic, fuGdpbNr);
2410 } // if( ( uPattern >> uBit ) & 0x1 )
2411
2412 break;
2413 } // case gdpbv100::PATT_RESYNC:
2414 default: {
2415 LOG(debug) << "Crazy pattern message, subtype " << usType;
2416 break;
2417 } // default
2418 } // switch( usType )
2419
2420 return;
2421}
2422
2424{
2425 Int_t iMsgIndex = mess.getStarTrigMsgIndex();
2426
2427 switch (iMsgIndex) {
2428 case 0:
2429 fvhTokenMsgType[fuGdpbNr]->Fill(0);
2431 break;
2432 case 1:
2433 fvhTokenMsgType[fuGdpbNr]->Fill(1);
2436 break;
2437 case 2:
2438 fvhTokenMsgType[fuGdpbNr]->Fill(2);
2440 break;
2441 case 3: {
2442 fvhTokenMsgType[fuGdpbNr]->Fill(3);
2443
2444 ULong64_t ulNewGdpbTsFull = (fvulGdpbTsMsb[fuGdpbNr] << 24) + (fvulGdpbTsLsb[fuGdpbNr]);
2445 ULong64_t ulNewStarTsFull =
2447 UInt_t uNewToken = mess.getStarTokenStarD();
2448 UInt_t uNewDaqCmd = mess.getStarDaqCmdStarD();
2449 UInt_t uNewTrigCmd = mess.getStarTrigCmdStarD();
2450
2451 if ((uNewToken == fvuStarTokenLast[fuGdpbNr]) && (ulNewGdpbTsFull == fvulGdpbTsFullLast[fuGdpbNr])
2452 && (ulNewStarTsFull == fvulStarTsFullLast[fuGdpbNr]) && (uNewDaqCmd == fvuStarDaqCmdLast[fuGdpbNr])
2453 && (uNewTrigCmd == fvuStarTrigCmdLast[fuGdpbNr])) {
2454 UInt_t uTrigWord = ((fvuStarTrigCmdLast[fuGdpbNr] & 0x00F) << 16)
2455 + ((fvuStarDaqCmdLast[fuGdpbNr] & 0x00F) << 12) + ((fvuStarTokenLast[fuGdpbNr] & 0xFFF));
2456 LOG(warning) << "Possible error: identical STAR tokens found twice in "
2457 "a row => ignore 2nd! "
2458 << " TS " << fulCurrentTsIndex << " gDBB #" << fuGdpbNr << " "
2459 << Form("token = %5u ", fvuStarTokenLast[fuGdpbNr])
2460 << Form("gDPB ts = %12llu ", fvulGdpbTsFullLast[fuGdpbNr])
2461 << Form("STAR ts = %12llu ", fvulStarTsFullLast[fuGdpbNr])
2462 << Form("DAQ cmd = %2u ", fvuStarDaqCmdLast[fuGdpbNr])
2463 << Form("TRG cmd = %2u ", fvuStarTrigCmdLast[fuGdpbNr]) << Form("TRG Wrd = %5x ", uTrigWord);
2464 return;
2465 } // if exactly same message repeated
2466 /*
2467 if( (uNewToken != fuStarTokenLast[fuGdpbNr] + 1) &&
2468 0 < fvulGdpbTsFullLast[fuGdpbNr] && 0 < fvulStarTsFullLast[fuGdpbNr] &&
2469 ( 4095 != fvuStarTokenLast[fuGdpbNr] || 1 != uNewToken) )
2470 LOG(warning) << "Possible error: STAR token did not increase by exactly 1! "
2471 << " gDBB #" << fuGdpbNr << " "
2472 << Form("old = %5u vs new = %5u ", fvuStarTokenLast[fuGdpbNr], uNewToken)
2473 << Form("old = %12llu vs new = %12llu ", fvulGdpbTsFullLast[fuGdpbNr], ulNewGdpbTsFull)
2474 << Form("old = %12llu vs new = %12llu ", fvulStarTsFullLast[fuGdpbNr], ulNewStarTsFull)
2475 << Form("old = %2u vs new = %2u ", fvuStarDaqCmdLast[fuGdpbNr], uNewDaqCmd)
2476 << Form("old = %2u vs new = %2u ", fvuStarTrigCmdLast[fuGdpbNr], uNewTrigCmd);
2477*/
2478 // STAR TS counter reset detection
2479 if (ulNewStarTsFull < fvulStarTsFullLast[fuGdpbNr])
2480 LOG(debug) << "Probable reset of the STAR TS: old = " << Form("%16llu", fvulStarTsFullLast[fuGdpbNr])
2481 << " new = " << Form("%16llu", ulNewStarTsFull) << " Diff = -"
2482 << Form("%8llu", fvulStarTsFullLast[fuGdpbNr] - ulNewStarTsFull);
2483
2484 /*
2485 LOG(info) << "Updating trigger token for " << fuGdpbNr
2486 << " " << fuStarTokenLast[fuGdpbNr] << " " << uNewToken;
2487*/
2488 ULong64_t ulGdpbTsDiff = ulNewGdpbTsFull - fvulGdpbTsFullLast[fuGdpbNr];
2489 fvulGdpbTsFullLast[fuGdpbNr] = ulNewGdpbTsFull;
2490 fvulStarTsFullLast[fuGdpbNr] = ulNewStarTsFull;
2491 fvuStarTokenLast[fuGdpbNr] = uNewToken;
2492 fvuStarDaqCmdLast[fuGdpbNr] = uNewDaqCmd;
2493 fvuStarTrigCmdLast[fuGdpbNr] = uNewTrigCmd;
2494
2496 if (fuCurrentMs < fuCoreMs) {
2498 if (0 <= fdStartTime) {
2503 } // if( fuHistoryHistoSize < 1e-9 * (fulGdpbTsFullLast * gdpbv100::kdClockCycleSizeNs - fdStartTime) )
2504
2505 fvhTriggerRate[fuGdpbNr]->Fill(1e-9
2516 } // if( 0 < fdStartTime )
2517 else
2520 } // if( fuCurrentMs < fuCoreMs )
2521
2522 break;
2523 } // case 3
2524 default: LOG(error) << "Unknown Star Trigger messageindex: " << iMsgIndex;
2525 } // switch( iMsgIndex )
2526}
2527
2529
2531{
2532 // Printout some stats on what was unpacked
2533 TString message_type;
2534 for (unsigned int i = 0; i < fviMsgCounter.size(); ++i) {
2535 switch (i) {
2536 case 0: message_type = "NOP"; break;
2537 case 1: message_type = "HIT"; break;
2538 case 2: message_type = "EPOCH"; break;
2539 case 3: message_type = "SYNC"; break;
2540 case 4: message_type = "AUX"; break;
2541 case 5: message_type = "EPOCH2"; break;
2542 case 6: message_type = "GET4"; break;
2543 case 7: message_type = "SYS"; break;
2544 case 8: message_type = "GET4_SLC"; break;
2545 case 9: message_type = "GET4_32B"; break;
2546 case 10: message_type = "GET4_SYS"; break;
2547 default: message_type = "UNKNOWN"; break;
2548 } // switch(i)
2549 LOG(info) << message_type << " messages: " << fviMsgCounter[i];
2550 } // for (unsigned int i=0; i< fviMsgCounter.size(); ++i)
2551
2552 LOG(info) << "-------------------------------------";
2553 for (UInt_t i = 0; i < fuNrOfGdpbs; ++i) {
2554 for (UInt_t j = 0; j < fuNrOfGet4PerGdpb; ++j) {
2555 LOG(info) << "Last epoch for gDPB: " << std::hex << std::setw(4) << i << std::dec << " , GET4 " << std::setw(4)
2556 << j << " => " << fvulCurrentEpoch[GetArrayIndex(i, j)];
2557 } // for (UInt_t j = 0; j < fuNrOfGet4PerGdpb; ++j)
2558 } // for (UInt_t i = 0; i < fuNrOfGdpbs; ++i)
2559 LOG(info) << "-------------------------------------";
2560
2561
2563 if (kTRUE == fbPulserModeEnable) {
2565 fhTimeMeanPulser->Reset();
2566 fhTimeRmsPulser->Reset();
2567
2568 for (UInt_t uFeeA = 0; uFeeA < fuNrOfFeePerGdpb * fuNrOfGdpbs; uFeeA++)
2569 for (UInt_t uFeeB = 0; uFeeB < fuNrOfFeePerGdpb * fuNrOfGdpbs; uFeeB++)
2570 if (NULL != fvhTimeDiffPulser[uFeeA][uFeeB]) {
2571 fhTimeMeanPulser->Fill(uFeeA, uFeeB, fvhTimeDiffPulser[uFeeA][uFeeB]->GetMean());
2572 fhTimeRmsPulser->Fill(uFeeA, uFeeB, fvhTimeDiffPulser[uFeeA][uFeeB]->GetRMS());
2573 } // for( UInt_t uChan = 0; uChan < kuNbChanTest - 1; uChan++)
2574
2577 } // if( kTRUE == fbPulserModeEnable )
2578
2579 SaveAllHistos();
2580}
2581
2583{
2584 TDirectory* oldDir = NULL;
2585 TFile* histoFile = NULL;
2586 if ("" != sFileName) {
2587 // Store current directory position to allow restore later
2588 oldDir = gDirectory;
2589 // open separate histo file in recreate mode
2590 histoFile = new TFile(sFileName, "RECREATE");
2591 histoFile->cd();
2592 } // if( "" != sFileName )
2593
2594 gDirectory->mkdir("Tof_Raw_gDPB");
2595 gDirectory->cd("Tof_Raw_gDPB");
2596
2597 fhMessType->Write();
2598 fhSysMessType->Write();
2599 fhGet4MessType->Write();
2600 fhGet4ChanScm->Write();
2601 fhGet4ChanErrors->Write();
2602 fhGet4EpochFlags->Write();
2603
2604 fhGdpbMessType->Write();
2605 fhGdpbSysMessType->Write();
2606 fhGdpbSysMessPattType->Write();
2607 fhGdpbEpochFlags->Write();
2608 fhGdpbEpochSyncEvo->Write();
2609 fhGdpbEpochMissEvo->Write();
2610
2611 fhScmScalerCounters->Write();
2612 fhScmDeadtimeCounters->Write();
2613 fhScmSeuCounters->Write();
2614 fhScmSeuCountersEvo->Write();
2615
2616 fhPatternMissmatch->Write();
2617 fhPatternEnable->Write();
2618 fhPatternResync->Write();
2619
2620 for (UInt_t uTotPlot = 0; uTotPlot < fvhRawTot_gDPB.size(); ++uTotPlot)
2621 fvhRawTot_gDPB[uTotPlot]->Write();
2622
2623 for (UInt_t uTotPlot = 0; uTotPlot < fvhRemapTot_gDPB.size(); ++uTotPlot)
2624 fvhRemapTot_gDPB[uTotPlot]->Write();
2625
2626 for (UInt_t uMod = 0; uMod < fuNrOfModules; uMod++) {
2627 fvhRemapTotSideA_mod[uMod]->Write();
2628 fvhRemapTotSideB_mod[uMod]->Write();
2629
2630 fvhModRate[uMod]->Write();
2631 fvhModErrorRate[uMod]->Write();
2632 fvhModErrorRatio[uMod]->Write();
2633 } // for( UInt_t uMod = 0; uMod < fuNrOfModules; uMod ++ )
2634
2635 for (UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb) {
2636 fvhGdpbGet4MessType[uGdpb]->Write();
2637 fvhGdpbGet4ChanScm[uGdpb]->Write();
2638 fvhGdpbGet4ChanErrors[uGdpb]->Write();
2639
2640 fvhGdpbPatternMissmatchEvo[uGdpb]->Write();
2641 fvhGdpbPatternEnableEvo[uGdpb]->Write();
2642 fvhGdpbPatternResyncEvo[uGdpb]->Write();
2643
2644 fvhRawFt_gDPB[uGdpb]->Write();
2645 fvhChCount_gDPB[uGdpb]->Write();
2646 fvhChannelRate_gDPB[uGdpb]->Write();
2647 fvhRemapChCount_gDPB[uGdpb]->Write();
2648 fvhRemapChRate_gDPB[uGdpb]->Write();
2649
2650 for (UInt_t uFee = 0; uFee < fuNrOfFeePerGdpb; ++uFee) {
2651 fvhFeeRate_gDPB[uGdpb * fuNrOfFeePerGdpb + uFee]->Write();
2652 fvhFeeErrorRate_gDPB[uGdpb * fuNrOfFeePerGdpb + uFee]->Write();
2653 fvhFeeErrorRatio_gDPB[uGdpb * fuNrOfFeePerGdpb + uFee]->Write();
2654 fvhFeeRateLong_gDPB[uGdpb * fuNrOfFeePerGdpb + uFee]->Write();
2655 fvhFeeErrorRateLong_gDPB[uGdpb * fuNrOfFeePerGdpb + uFee]->Write();
2656 fvhFeeErrorRatioLong_gDPB[uGdpb * fuNrOfFeePerGdpb + uFee]->Write();
2657 } // for (UInt_t uFee = 0; uFee < fuNrOfFeePerGdpb; ++ uFee)
2658
2659 } // for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
2660 if (kTRUE == fbPulserModeEnable) {
2661 fhTimeMeanPulser->Write();
2662 fhTimeRmsPulser->Write();
2663 fhTimeRmsZoomFitPuls->Write();
2664 fhTimeResFitPuls->Write();
2665 } // if( kTRUE == fbPulserModeEnable )
2666 gDirectory->cd("..");
2667
2669 gDirectory->mkdir("Star_Raw");
2670 gDirectory->cd("Star_Raw");
2671 for (UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb) {
2672 fvhTokenMsgType[uGdpb]->Write();
2673 fvhTriggerRate[uGdpb]->Write();
2674 fvhCmdDaqVsTrig[uGdpb]->Write();
2675 fvhStarTokenEvo[uGdpb]->Write();
2676 fvhStarTrigGdpbTsEvo[uGdpb]->Write();
2677 fvhStarTrigStarTsEvo[uGdpb]->Write();
2678 } // for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
2679 gDirectory->cd("..");
2680
2682 if (kTRUE == fbPulserModeEnable) {
2683 gDirectory->mkdir("TofDt");
2684 gDirectory->cd("TofDt");
2685 for (UInt_t uFeeA = 0; uFeeA < fuNrOfFeePerGdpb * fuNrOfGdpbs; uFeeA++)
2686 for (UInt_t uFeeB = 0; uFeeB < fuNrOfFeePerGdpb * fuNrOfGdpbs; uFeeB++)
2687 if (NULL != fvhTimeDiffPulser[uFeeA][uFeeB]) fvhTimeDiffPulser[uFeeA][uFeeB]->Write();
2688 gDirectory->cd("..");
2689
2691 gDirectory->mkdir("TofDtEvo");
2692 gDirectory->cd("TofDtEvo");
2693 for (UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb) {
2694 for (UInt_t uGbtx = 0; uGbtx < kuNbGbtxPerGdpb - 1; ++uGbtx)
2695 fvhPulserTimeDiffEvoGbtxGbtx[uGdpb * (kuNbGbtxPerGdpb - 1) + uGbtx]->Write();
2696
2697 for (UInt_t uGdpbB = uGdpb + 1; uGdpbB < fuNrOfGdpbs; ++uGdpbB)
2698 fvvhPulserTimeDiffEvoGdpbGdpb[uGdpb][uGdpbB]->Write();
2699 } // for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
2700 gDirectory->cd("..");
2701 } // if( kTRUE == fbPulserModeEnable )
2702
2704 if (kTRUE == fbCoincMapsEnable) {
2705 gDirectory->mkdir("TofCoinc");
2706 gDirectory->cd("TofCoinc");
2707 for (UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb) {
2708 fvhCoincMapAllChanGdpb[uGdpb]->Write();
2709 fvhCoincMeanAllChanGdpb[uGdpb]->Write();
2710 } // for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
2711 gDirectory->cd("..");
2712 } // if( kTRUE == fbCoincMapsEnable )
2713
2714 gDirectory->mkdir("Flib_Raw");
2715 gDirectory->cd("Flib_Raw");
2716 for (UInt_t uLinks = 0; uLinks < fvhMsSzPerLink.size(); uLinks++) {
2717 if (NULL == fvhMsSzPerLink[uLinks]) continue;
2718
2719 fvhMsSzPerLink[uLinks]->Write();
2720 fvhMsSzTimePerLink[uLinks]->Write();
2721 } // for( UInt_t uLinks = 0; uLinks < fvhMsSzPerLink.size(); uLinks++ )
2722
2723 TH1* pMissedTsH1 = dynamic_cast<TH1*>(gROOT->FindObjectAny("Missed_TS"));
2724 if (NULL != pMissedTsH1) pMissedTsH1->Write();
2725
2726 TProfile* pMissedTsEvoP = dynamic_cast<TProfile*>(gROOT->FindObjectAny("Missed_TS_Evo"));
2727 if (NULL != pMissedTsEvoP) pMissedTsEvoP->Write();
2728
2729 gDirectory->cd("..");
2730
2731
2732 if ("" != sFileName) {
2733 // Restore original directory position
2734 histoFile->Close();
2735 oldDir->cd();
2736 } // if( "" != sFileName )
2737 if ("" != sFileName) {
2738 // Restore original directory position
2739 histoFile->Close();
2740 oldDir->cd();
2741 } // if( "" != sFileName )
2742}
2743
2745{
2746 LOG(info) << "Reseting all TOF histograms.";
2747
2748 fhMessType->Reset();
2749 fhSysMessType->Reset();
2750 fhGet4MessType->Reset();
2751 fhGet4ChanScm->Reset();
2752 fhGet4ChanErrors->Reset();
2753 fhGet4EpochFlags->Reset();
2754
2755 fhGdpbMessType->Reset();
2756 fhGdpbSysMessType->Reset();
2757 fhGdpbSysMessPattType->Reset();
2758 fhGdpbEpochFlags->Reset();
2759 fhGdpbEpochSyncEvo->Reset();
2760 fhGdpbEpochMissEvo->Reset();
2761
2762 fhScmScalerCounters->Reset();
2763 fhScmDeadtimeCounters->Reset();
2764 fhScmSeuCounters->Reset();
2765 fhScmSeuCountersEvo->Reset();
2766
2767 fhPatternMissmatch->Reset();
2768 fhPatternEnable->Reset();
2769 fhPatternResync->Reset();
2770
2771 for (UInt_t uTotPlot = 0; uTotPlot < fvhRawTot_gDPB.size(); ++uTotPlot)
2772 fvhRawTot_gDPB[uTotPlot]->Reset();
2773
2774 for (UInt_t uTotPlot = 0; uTotPlot < fvhRemapTot_gDPB.size(); ++uTotPlot)
2775 fvhRemapTot_gDPB[uTotPlot]->Reset();
2776
2777 for (UInt_t uMod = 0; uMod < fuNrOfModules; uMod++) {
2778 fvhRemapTotSideA_mod[uMod]->Reset();
2779 fvhRemapTotSideB_mod[uMod]->Reset();
2780
2781 fvhModRate[uMod]->Reset();
2782 fvhModErrorRate[uMod]->Reset();
2783 fvhModErrorRatio[uMod]->Reset();
2784 } // for( UInt_t uMod = 0; uMod < fuNrOfModules; uMod ++ )
2785
2786 for (UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb) {
2787 fvhGdpbGet4MessType[uGdpb]->Reset();
2788 fvhGdpbGet4ChanScm[uGdpb]->Reset();
2789 fvhGdpbGet4ChanErrors[uGdpb]->Reset();
2790
2791 fvhGdpbPatternMissmatchEvo[uGdpb]->Reset();
2792 fvhGdpbPatternEnableEvo[uGdpb]->Reset();
2793 fvhGdpbPatternResyncEvo[uGdpb]->Reset();
2794
2795 fvhRawFt_gDPB[uGdpb]->Reset();
2796 fvhChCount_gDPB[uGdpb]->Reset();
2797 fvhChannelRate_gDPB[uGdpb]->Reset();
2798 fvhRemapChCount_gDPB[uGdpb]->Reset();
2799 fvhRemapChRate_gDPB[uGdpb]->Reset();
2800
2801 for (UInt_t uFee = 0; uFee < fuNrOfFeePerGdpb; ++uFee) {
2802 fvhFeeRate_gDPB[uGdpb * fuNrOfFeePerGdpb + uFee]->Reset();
2803 fvhFeeErrorRate_gDPB[uGdpb * fuNrOfFeePerGdpb + uFee]->Reset();
2804 fvhFeeErrorRatio_gDPB[uGdpb * fuNrOfFeePerGdpb + uFee]->Reset();
2805 fvhFeeRateLong_gDPB[uGdpb * fuNrOfFeePerGdpb + uFee]->Reset();
2806 fvhFeeErrorRateLong_gDPB[uGdpb * fuNrOfFeePerGdpb + uFee]->Reset();
2807 fvhFeeErrorRatioLong_gDPB[uGdpb * fuNrOfFeePerGdpb + uFee]->Reset();
2808 } // for (UInt_t uFee = 0; uFee < fuNrOfFeePerGdpb; ++ uFee)
2809
2810 } // for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
2811
2813 for (UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb) {
2814 fvhTokenMsgType[uGdpb]->Reset();
2815 fvhTriggerRate[uGdpb]->Reset();
2816 fvhCmdDaqVsTrig[uGdpb]->Reset();
2817 fvhStarTokenEvo[uGdpb]->Reset();
2818 fvhStarTrigGdpbTsEvo[uGdpb]->Reset();
2819 fvhStarTrigStarTsEvo[uGdpb]->Reset();
2820 } // for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
2821
2823 if (kTRUE == fbPulserModeEnable) {
2824 for (UInt_t uFeeA = 0; uFeeA < fuNrOfFeePerGdpb * fuNrOfGdpbs; uFeeA++)
2825 for (UInt_t uFeeB = 0; uFeeB < fuNrOfFeePerGdpb * fuNrOfGdpbs; uFeeB++)
2826 if (NULL != fvhTimeDiffPulser[uFeeA][uFeeB]) fvhTimeDiffPulser[uFeeA][uFeeB]->Reset();
2827
2828 fhTimeMeanPulser->Reset();
2829 fhTimeRmsPulser->Reset();
2830 fhTimeRmsZoomFitPuls->Reset();
2831 fhTimeResFitPuls->Reset();
2832
2833 for (UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb) {
2834 for (UInt_t uGbtx = 0; uGbtx < kuNbGbtxPerGdpb - 1; ++uGbtx)
2835 fvhPulserTimeDiffEvoGbtxGbtx[uGdpb * (kuNbGbtxPerGdpb - 1) + uGbtx]->Reset();
2836
2837 for (UInt_t uGdpbB = uGdpb + 1; uGdpbB < fuNrOfGdpbs; ++uGdpbB)
2838 fvvhPulserTimeDiffEvoGdpbGdpb[uGdpb][uGdpbB]->Reset();
2839 } // for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
2840 } // if( kTRUE == fbPulserModeEnable )
2841
2843 if (kTRUE == fbCoincMapsEnable) {
2844 for (UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb) {
2845 fvhCoincMapAllChanGdpb[uGdpb]->Reset();
2846 fvhCoincMeanAllChanGdpb[uGdpb]->Reset();
2847 } // for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
2848 } // if( kTRUE == fbCoincMapsEnable )
2849
2850 for (UInt_t uLinks = 0; uLinks < fvhMsSzPerLink.size(); uLinks++) {
2851 if (NULL == fvhMsSzPerLink[uLinks]) continue;
2852
2853 fvhMsSzPerLink[uLinks]->Reset();
2854 fvhMsSzTimePerLink[uLinks]->Reset();
2855 } // for( UInt_t uLinks = 0; uLinks < fvhMsSzPerLink.size(); uLinks++ )
2856
2857 fdStartTime = -1;
2858 fdStartTimeLong = -1;
2859 fdStartTimeMsSz = -1;
2860}
2862{
2863 for (UInt_t uGdpbLoop = 0; uGdpbLoop < fuNrOfGdpbs; ++uGdpbLoop) {
2864 fvhChannelRate_gDPB[uGdpbLoop]->Reset();
2865 fvhRemapChRate_gDPB[uGdpbLoop]->Reset();
2866 for (UInt_t uFeeLoop = 0; uFeeLoop < fuNrOfFeePerGdpb; ++uFeeLoop) {
2867 fvhFeeRate_gDPB[(uGdpbLoop * fuNrOfFeePerGdpb) + uFeeLoop]->Reset();
2868 fvhFeeErrorRate_gDPB[(uGdpbLoop * fuNrOfFeePerGdpb) + uFeeLoop]->Reset();
2869 fvhFeeErrorRatio_gDPB[(uGdpbLoop * fuNrOfFeePerGdpb) + uFeeLoop]->Reset();
2870 } // for( UInt_t uFeeLoop = 0; uFeeLoop < fuNrOfFeePerGdpb; ++uFeeLoop )
2871 fvhTriggerRate[uGdpbLoop]->Reset();
2872 fvhStarTokenEvo[uGdpbLoop]->Reset();
2873 fvhStarTrigGdpbTsEvo[uGdpbLoop]->Reset();
2874 fvhStarTrigStarTsEvo[uGdpbLoop]->Reset();
2875 } // for( UInt_t uGdpbLoop = 0; uGdpbLoop < fuNrOfGdpbs; ++uGdpbLoop )
2876
2877 for (UInt_t uMod = 0; uMod < fuNrOfModules; uMod++) {
2878 fvhModRate[uMod]->Reset();
2879 fvhModErrorRate[uMod]->Reset();
2880 fvhModErrorRatio[uMod]->Reset();
2881 } // for( UInt_t uMod = 0; uMod < fuNrOfModules; uMod ++ )
2882
2883 fdStartTime = -1;
2884}
2886{
2887 for (UInt_t uGdpbLoop = 0; uGdpbLoop < fuNrOfGdpbs; uGdpbLoop++) {
2888 for (UInt_t uFeeLoop = 0; uFeeLoop < fuNrOfFeePerGdpb; uFeeLoop++) {
2889 fvhFeeRateLong_gDPB[(uGdpbLoop * fuNrOfFeePerGdpb) + uFeeLoop]->Reset();
2890 fvhFeeErrorRateLong_gDPB[(uGdpbLoop * fuNrOfFeePerGdpb) + uFeeLoop]->Reset();
2891 fvhFeeErrorRatioLong_gDPB[(uGdpbLoop * fuNrOfFeePerGdpb) + uFeeLoop]->Reset();
2892 } // for (UInt_t uFeeLoop = 0; uFeeLoop < fuNrOfFeePerGdpb; uFeeLoop++)
2893 } // for (UInt_t uFeeLoop = 0; uFeeLoop < fuNrOfFeePerGdpb; uFeeLoop++)
2894
2895 fdStartTimeLong = -1;
2896}
2897
2899{
2900 if (kFALSE == fbPulserModeEnable) {
2901 LOG(error) << "CbmStar2019MonitorTof::UpdateZoomedFit => "
2902 << "Pulser mode not enabled in root macro, doinb nothing !!! ";
2903 return;
2904 } // if( kFALSE == fbPulserModeEnable )
2905
2906 // Only do something is the user defined the width he want for the zoom
2907 if (0.0 < fdFitZoomWidthPs) {
2908 // Reset summary histograms for safety
2909 fhTimeRmsZoomFitPuls->Reset();
2910 fhTimeResFitPuls->Reset();
2911
2912 Double_t dRes = 0;
2914
2915 for (UInt_t uFeeA = 0; uFeeA < fuNrOfFeePerGdpb * fuNrOfGdpbs; uFeeA++)
2916 for (UInt_t uFeeB = 0; uFeeB < fuNrOfFeePerGdpb * fuNrOfGdpbs; uFeeB++)
2917 if (NULL != fvhTimeDiffPulser[uFeeA][uFeeB]) {
2918 // Check that we have at least 1 entry
2919 if (0 == fvhTimeDiffPulser[uFeeA][uFeeB]->GetEntries()) {
2920 fhTimeRmsZoomFitPuls->Fill(uFeeA, uFeeB, 0.0);
2921 fhTimeResFitPuls->Fill(uFeeA, uFeeB, 0.0);
2922 LOG(info) << "CbmStar2019MonitorTof::UpdateZoomedFit => Empty input "
2923 << "for FEE pair " << uFeeA << " and " << uFeeB << " !!! ";
2924 continue;
2925 } // if( 0 == fvhTimeDiffPulser[uFeeA][uFeeB]->GetEntries() )
2926
2927 // Read the peak position (bin with max counts) + total nb of entries
2928 Int_t iBinWithMax = fvhTimeDiffPulser[uFeeA][uFeeB]->GetMaximumBin();
2929 Double_t dNbCounts = fvhTimeDiffPulser[uFeeA][uFeeB]->Integral();
2930
2931 // Zoom the X axis to +/- ZoomWidth around the peak position
2932 Double_t dPeakPos = fvhTimeDiffPulser[uFeeA][uFeeB]->GetXaxis()->GetBinCenter(iBinWithMax);
2933 fvhTimeDiffPulser[uFeeA][uFeeB]->GetXaxis()->SetRangeUser(dPeakPos - fdFitZoomWidthPs,
2934 dPeakPos + fdFitZoomWidthPs);
2935
2936 // Read integral and check how much we lost due to the zoom (% loss allowed)
2937 Double_t dZoomCounts = fvhTimeDiffPulser[uFeeA][uFeeB]->Integral();
2938 if ((dZoomCounts / dNbCounts) < 0.99) {
2939 fhTimeRmsZoomFitPuls->Fill(uFeeA, uFeeB, 0.0);
2940 fhTimeResFitPuls->Fill(uFeeA, uFeeB, 0.0);
2941 LOG(warning) << "CbmStar2019MonitorTof::UpdateZoomedFit => Zoom too strong, "
2942 << "more than 1% loss for FEE pair " << uFeeA << " and " << uFeeB << " !!! ";
2943 continue;
2944 } // if( ( dZoomCounts / dNbCounts ) < 0.99 )
2945
2946 // Fill new RMS after zoom into summary histo
2947 fhTimeRmsZoomFitPuls->Fill(uFeeA, uFeeB, fvhTimeDiffPulser[uFeeA][uFeeB]->GetRMS());
2948
2949
2950 // Fit using zoomed boundaries + starting gaussian width, store into summary histo
2951 dRes = 0;
2952 fitFuncPairs[uFeeA][uFeeB] = new TF1(Form("fPair_%02d_%02d", uFeeA, uFeeB), "gaus",
2953 dPeakPos - fdFitZoomWidthPs, dPeakPos + fdFitZoomWidthPs);
2954 // Fix the Mean fit value around the Histogram Mean
2955 fitFuncPairs[uFeeA][uFeeB]->SetParameter(0, dZoomCounts);
2956 fitFuncPairs[uFeeA][uFeeB]->SetParameter(1, dPeakPos);
2957 fitFuncPairs[uFeeA][uFeeB]->SetParameter(2, 200.0); // Hardcode start with ~4*BinWidth, do better later
2958 // Using integral instead of bin center seems to lead to unrealistic values => no "I"
2959 fvhTimeDiffPulser[uFeeA][uFeeB]->Fit(Form("fPair_%02d_%02d", uFeeA, uFeeB), "QRM0");
2960 // Get Sigma
2961 dRes = fitFuncPairs[uFeeA][uFeeB]->GetParameter(2);
2962 // Cleanup memory
2963 delete fitFuncPairs[uFeeA][uFeeB];
2964 // Fill summary
2965 fhTimeResFitPuls->Fill(uFeeA, uFeeB, dRes / TMath::Sqrt2());
2966
2967
2968 LOG(info) << "CbmStar2019MonitorTof::UpdateZoomedFit => "
2969 << "For FEE pair " << uFeeA << " and " << uFeeB
2970 << " we have zoomed RMS = " << fvhTimeDiffPulser[uFeeA][uFeeB]->GetRMS() << " and a resolution of "
2971 << dRes / TMath::Sqrt2();
2972
2973 // Restore original axis state?
2974 fvhTimeDiffPulser[uFeeA][uFeeB]->GetXaxis()->UnZoom();
2975 } // loop on uFeeA and uFeeB + check if corresponding fvhTimeDiffPulser exists
2976 } // if( 0.0 < fdFitZoomWidthPs )
2977 else {
2978 LOG(error) << "CbmStar2019MonitorTof::UpdateZoomedFit => Zoom width not defined, "
2979 << "please use SetFitZoomWidthPs, e.g. in macro, before trying this "
2980 "update !!!";
2981 } // else of if( 0.0 < fdFitZoomWidthPs )
2982}
2983
ClassImp(CbmConverterManager)
Histogram manager.
Bool_t bMcbmMoniTofPrintAllEpochsEna
Bool_t bMcbmMoniTofSaveHistos
Bool_t bMcbmMoniTofResetHistos
Bool_t bMcbmMoniTofUpdateZoomedFit
Bool_t bMcbmMoniTofRawDataPrint
Bool_t bMcbmMoniTofPrintAllHitsEna
static double dTsStartTime
static constexpr size_t size()
Definition KfSimdPseudo.h:2
int Int_t
bool Bool_t
void FillStarTrigInfo(gdpbv100::Message)
std::vector< ULong64_t > fvulStarTsMsb
std::vector< TH2 * > fvhCmdDaqVsTrig
std::vector< UInt_t > fvuStarDaqCmdLast
std::vector< TH2 * > fvhStarTokenEvo
std::vector< ULong64_t > fvulStarTsFullLast
std::vector< ULong64_t > fvulGdpbTsMsb
std::vector< TProfile * > fvhStarTrigStarTsEvo
std::vector< UInt_t > fvuStarTrigCmdLast
std::vector< ULong64_t > fvulGdpbTsLsb
std::vector< TH1 * > fvhTriggerRate
std::vector< UInt_t > fvuStarTokenLast
std::vector< TH1 * > fvhTokenMsgType
std::vector< ULong64_t > fvulStarTsMid
std::vector< TProfile * > fvhStarTrigGdpbTsEvo
std::vector< ULong64_t > fvulGdpbTsFullLast
std::vector< ULong64_t > fvulStarTsMsb
UInt_t ConvertElinkToGet4(UInt_t uElinkIdx)
ULong64_t fulCurrentEpochTime
Epoch + Epoch Cycle.
std::vector< TProfile * > fvhStarTrigGdpbTsEvo
std::vector< TProfile * > fvhStarTrigStarTsEvo
std::vector< UInt_t > fvuFeeNbHitsLastMs
Buffer for pulser channels.
std::vector< TH1 * > fvhFeeErrorRateLong_gDPB
std::vector< TH2 * > fvhGdpbGet4ChanScm
std::vector< Int_t > fviNrOfRpc
std::vector< ULong64_t > fvulStarTsMid
std::vector< TProfile2D * > fvhCoincMeanAllChanGdpb
void FillHitInfo(gdpbv100::Message)
virtual void AddMsComponentToList(size_t component, UShort_t usDetectorId)
virtual void SetNbMsInTs(size_t uCoreMsNb, size_t uOverlapMsNb)
void PrintSysInfo(gdpbv100::Message)
std::vector< TH1 * > fvhModRate
module plots
std::vector< TH2 * > fvhGdpbPatternEnableEvo
TH2 * fhGet4MessType
Per GET4 in system.
std::vector< TH2 * > fvhRawFt_gDPB
TODO: Channel rate plots!
std::vector< TH2 * > fvhCmdDaqVsTrig
void SaveAllHistos(TString sFileName="")
std::vector< TProfile * > fvhPulserTimeDiffEvoGbtxGbtx
std::vector< std::vector< TH1 * > > fvhTimeDiffPulser
std::vector< TProfile * > fvhFeeErrorRatioLong_gDPB
std::vector< Int_t > fviRpcSide
std::vector< Bool_t > fvbFirstEpochSeen
const UInt_t kuNbBinsDt
[ fuNrOfGdpbs ][ fuNrOfChannelsPerGdpb ]
std::vector< TH2 * > fvhRawTot_gDPB
std::vector< std::vector< UInt_t > > fvuCoincNbHitsLastMs
[ fuFeeNr ]
std::vector< TH2 * > fvhGdpbGet4MessType
Per GET4 in gDPB.
void ProcessEpochCycle(uint64_t ulCycleData)
std::vector< TProfile * > fvhModErrorRatio
std::vector< TH2 * > fvhRemapChRate_gDPB
std::vector< ULong64_t > fvulCurrentEpochCycle
std::vector< TH2 * > fvhRemapTotSideB_mod
std::vector< ULong64_t > fvulCurrentEpochFull
Epoch cycle from the Ms Start message and Epoch counter flip.
std::vector< UInt_t > fvuStarTrigCmdLast
std::map< UInt_t, UInt_t > fGdpbIdIndexMap
Map of ID to index for the gDPBs.
void FillEpochInfo(gdpbv100::Message)
void FillPattInfo(gdpbv100::Message)
std::vector< TH1 * > fvhMsSzPerLink
static const UInt_t kuNbGet4PerGbtx
TH2 * fhScmScalerCounters
Slow control messages.
std::vector< std::vector< Double_t > > fvdCoincTsLastHit
[ fuNrOfGdpbs ][ fuNrOfChannelsPerGdpb ]
Int_t GetArrayIndex(Int_t gdpbId, Int_t get4Id)
void FillStarTrigInfo(gdpbv100::Message)
TH2 * fhPatternMissmatch
Pattern messages per gDPB.
std::vector< TH1 * > fvhTriggerRate
std::vector< TH1 * > fvhFeeRate_gDPB
std::vector< Double_t > fdTsLastPulserHit
[ fuFeeNr ]
std::vector< TH2 * > fvhStarTokenEvo
std::vector< TH2 * > fvhGdpbGet4ChanErrors
std::vector< Int_t > fviModuleId
std::vector< size_t > fvMsComponentsList
FLES containers.
std::vector< UInt_t > fvuStarDaqCmdLast
std::vector< ULong64_t > fvulStarTsFullLast
std::vector< std::vector< gdpbv100::Message > > fvmEpSupprBuffer
Buffer for suppressed epoch processing.
std::vector< TH2 * > fvhGdpbPatternMissmatchEvo
Per MS in gDPB.
std::vector< TH2 * > fvhChannelRate_gDPB
std::vector< UInt_t > fvuPadiToGet4
std::vector< TH2 * > fvhCoincMapAllChanGdpb
std::vector< TH2 * > fvhRemapTotSideA_mod
std::vector< UInt_t > fvuElinkToGet4
5 FEE with 8 GET4 each
std::vector< TH1 * > fvhChCount_gDPB
std::vector< ULong64_t > fvulCurrentEpoch
std::vector< ULong64_t > fvulGdpbTsMsb
std::vector< ULong64_t > fvulGdpbTsLsb
std::vector< TH1 * > fvhModErrorRate
std::vector< ULong64_t > fvulGdpbTsFullLast
std::vector< UInt_t > fvuGet4ToPadi
std::vector< int > fviMsgCounter
void PrintGenInfo(gdpbv100::Message)
CbmStar2019TofPar * fUnpackPar
std::vector< UInt_t > fvuStarTokenLast
std::chrono::time_point< std::chrono::system_clock > fTimeLastHistoSaving
std::vector< TH1 * > fvhFeeErrorRate_gDPB
std::vector< TH1 * > fvhFeeRateLong_gDPB
void PrintSlcInfo(gdpbv100::Message)
std::vector< UInt_t > fvuGet4ToElink
std::vector< TProfile * > fvhMsSzTimePerLink
std::vector< Int_t > fviRpcType
std::vector< TH2 * > fvhGdpbPatternResyncEvo
std::vector< TH2 * > fvhRemapTot_gDPB
std::vector< TH1 * > fvhTokenMsgType
size_t fuMsAcceptsPercent
/‍** Ignore Overlap Ms: all fuOverlapMsNb MS at the end of timeslice **‍/
std::vector< std::vector< TProfile * > > fvvhPulserTimeDiffEvoGdpbGdpb
std::vector< TH1 * > fvhRemapChCount_gDPB
std::vector< TProfile * > fvhFeeErrorRatio_gDPB
virtual Bool_t DoUnpack(const fles::Timeslice &ts, size_t component)
Data class with information on a STS local track.
uint16_t getGdpbHitIs24b() const
uint16_t getStarTrigMsgIndex() const
uint32_t getGdpbEpEpochNb() const
uint16_t getGdpbHit32Tot() const
uint16_t getGdpbSysPattType() const
uint64_t getStarTsMidStarC() const
uint16_t getGdpbSysSubType() const
uint64_t getGdpbTsLsbStarB() const
uint32_t getGdpbSlcChan() const
bool getGdpbEpMissmatch() const
double getMsgFullTimeD(uint64_t epoch) const
Returns expanded and adjusted time of message in double (in ns)
uint16_t getGdpbSysPattIndex() const
void setGdpbGenChipId(uint32_t v)
uint16_t getGdpbHitFineTs() const
uint32_t getGdpbSlcData() const
void setGdpbEpEpochNb(uint32_t v)
bool getGdpbEpSync() const
uint64_t getStarTsMsbStarB() const
uint32_t getGdpbHitFullTs() const
uint32_t getStarTrigCmdStarD() const
uint32_t getGdpbSysPattPattern() const
uint32_t getGdpbSlcType() const
bool isStarTrigger() const
Returns true is message type is MSG_STAR_TRI_A, _B, _C, _D (STAR Trigger message)
bool getGdpbEpEpochLoss() const
uint32_t getStarDaqCmdStarD() const
uint64_t getGdpbTsMsbStarA() const
uint16_t getGdpbGenChipId() const
uint16_t getGdpbSysErrUnused() const
void printDataCout(unsigned kind=msg_print_Prefix|msg_print_Data, uint32_t epoch=0) const
Print message in human readable format to cout.
bool getGdpbEpDataLoss() const
uint32_t getStarTokenStarD() const
bool getGdpbEpLinkId() const
uint32_t getGdpbSlcEdge() const
uint64_t getStarTsLsbStarD() const
uint8_t getMessageType() const
Returns the message type. Valid for all message types. 4 bit.
uint32_t getGdpbSysUnkwData() const
uint64_t getData() const
uint16_t getGdpbSysErrData() const
bool getGdpbSysErrEdge() const
uint16_t getGdpbSysErrChanId() const
uint64_t getMsgFullTime(uint64_t epoch) const
Returns expanded and adjusted time of message (in ns)
uint16_t getGdpbHitChanId() const
const double kdClockCycleSizeNs
@ GET4_V2X_ERR_LOST_EVT
@ GET4_V2X_ERR_TOT_RANGE
@ GET4_V2X_ERR_FIFO_WRITE
@ GET4_V2X_ERR_DLL_LOCK
@ GET4_V2X_ERR_EPOCH_OVERF
@ GET4_V2X_ERR_EP_CNT_SYNC
@ GET4_V2X_ERR_UNKNOWN
@ GET4_V2X_ERR_READ_INIT
@ GET4_V2X_ERR_DLL_RESET
@ GET4_V2X_ERR_ADD_RIS_EDG
@ GET4_V2X_ERR_TOK_RING_ST
@ GET4_V2X_ERR_TOKEN
@ GET4_V2X_ERR_TOT_OVERWRT
@ GET4_V2X_ERR_READOUT_ERR
@ GET4_V2X_ERR_EVT_DISCARD
@ GET4_V2X_ERR_UNPAIR_FALL
@ GET4_V2X_ERR_CHAN_STATE
@ GET4_V2X_ERR_SYNC
@ GET4_V2X_ERR_SEQUENCE_ER
const uint32_t kuChipIdMergedEpoch
const uint32_t kuFeePulserChannel
const uint32_t kuEpochCounterSz
const uint64_t kulEpochCycleFieldSz
@ SYS_GET4_SYNC_MISS