CbmRoot
Loading...
Searching...
No Matches
CbmMcbm2018MonitorSts.cxx
Go to the documentation of this file.
1/* Copyright (C) 2018-2020 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// ----- CbmMcbm2018MonitorSts -----
8// ----- Created 11/05/18 by P.-A. Loizeau -----
9// ----- -----
10// -----------------------------------------------------------------------------
11
13
15
16// Data
17
18// CbmRoot
19#include "CbmHistManager.h"
20#include "CbmMcbm2018StsPar.h"
21
22// FairRoot
23#include "FairRootManager.h"
24#include "FairRun.h"
25#include "FairRunOnline.h"
26#include "FairRuntimeDb.h"
27#include <Logger.h>
28
29// Root
30#include "TClonesArray.h"
31#include "THttpServer.h"
32#include "TMath.h"
33#include "TROOT.h"
34#include "TRandom.h"
35#include "TString.h"
36#include "TStyle.h"
37#include <TFile.h>
38
39// C++11
40#include <bitset>
41
42// C/C++
43#include <cstdint>
44#include <iomanip>
45#include <iostream>
46
50
57 , fbIgnoreOverlapMs(kFALSE)
58 , fUnpackParSts(NULL)
59 , fuNbModules(0)
62 , fuNrOfDpbs(0)
65 , fuNbFebs(0)
66 , fuNbStsXyters(0)
69 , fviFebType()
72 ,
73 /*
74 fuNrOfDpbs(0),
75 fDpbIdIndexMap(),
76 fuNbStsXyters(0),
77 fUnpackParSts->GetNbChanPerAsic()(0),
78 fuNbFebs(0),
79*/
80 fsHistoFileFullname("data/SetupHistos.root")
81 , fbPrintMessages(kFALSE)
82 , fPrintMessCtrl(stsxyter::MessagePrintMask::msg_print_Human)
83 ,
84 // fbEnableCoincidenceMaps( kFALSE ),
88 , fmMsgCounter()
90 , fuCurrDpbId(0)
91 , fuCurrDpbIdx(0)
101 , fvdPrevMsTime()
102 , fvdMsTime()
106 ,
107 // fvmChanHitsInTs(),
108 fdStartTime(-1.0)
109 , fdStartTimeMsSz(-1.0)
110 , ftStartTimeUnix(std::chrono::steady_clock::now())
111 , fvmHitsInMs()
114 , fuMaxNbMicroslices(100)
119 , fbLongHistoEnable(kFALSE)
123 , fdCoincCenter(0.0)
124 , fdCoincBorder(50.0)
127 , fHM(new CbmHistManager())
128 , fhStsMessType(NULL)
129 , fhStsSysMessType(NULL)
130 , fhStsMessTypePerDpb(NULL)
132 , fhStsStatusMessType(NULL)
136 , fhStsAllFebsHitRateEvo(nullptr)
137 , fhStsAllAsicsHitRateEvo(nullptr)
138 , fhStsFebAsicHitCounts(nullptr)
144 ,
145 // fhStsFebChanAdcCal(),
146 // fhStsFebChanAdcCalProf(),
163 ,
164 /*
165 fhStsFebChanDtCoinc(),
166 fhStsFebChanCoinc(),
167 fhStsModulePNCoincDt(),
168 fhStsModulePNCoincDtAsicP(),
169 fhStsModulePNCoincDtAsicN(),
170 fhStsModulePNCoincChan(),
171 fhStsModulePNCoincAdc(),
172 fhStsModuleCoincAdcChanP(),
173 fhStsModuleCoincAdcChanN(),
174 fhStsModuleCoincMap(),
175*/
195 , fcMsSizeAll(NULL)
196{
197}
198
200
202{
203 LOG(info) << "Initializing flib StsXyter unpacker for STS";
204
205 FairRootManager* ioman = FairRootManager::Instance();
206 if (ioman == NULL) { LOG(fatal) << "No FairRootManager instance"; }
207
208 return kTRUE;
209}
210
212{
213 LOG(info) << "Setting parameter containers for " << GetName();
214 fUnpackParSts = (CbmMcbm2018StsPar*) (FairRun::Instance()->GetRuntimeDb()->getContainer("CbmMcbm2018StsPar"));
215}
216
217
219{
220 LOG(info) << "Init parameter containers for " << GetName();
221
222 Bool_t bInit = InitStsParameters();
223 if (kTRUE == bInit) CreateHistograms();
224
225 return bInit;
226}
227
229{
230 LOG(info) << "ReInit parameter containers for " << GetName();
231
232 return InitStsParameters();
233}
234
236{
237
238
239 fuNbModules = fUnpackParSts->GetNbOfModules();
240 LOG(info) << "Nr. of STS Modules: " << fuNbModules;
241
244 for (UInt_t uModIdx = 0; uModIdx < fuNbModules; ++uModIdx) {
245 fviModuleType[uModIdx] = fUnpackParSts->GetModuleType(uModIdx);
246 fviModAddress[uModIdx] = fUnpackParSts->GetModuleAddress(uModIdx);
247 LOG(info) << "Module #" << std::setw(2) << uModIdx << " Type " << std::setw(4) << fviModuleType[uModIdx]
248 << " Address 0x" << std::setw(8) << std::hex << fviModAddress[uModIdx] << std::dec;
249 } // for( UInt_t uModIdx = 0; uModIdx < fuNbModules; ++uModIdx)
250
251 fuNrOfDpbs = fUnpackParSts->GetNrOfDpbs();
252 LOG(info) << "Nr. of STS DPBs: " << fuNrOfDpbs;
253
254 fDpbIdIndexMap.clear();
255 for (UInt_t uDpb = 0; uDpb < fuNrOfDpbs; ++uDpb) {
256 fDpbIdIndexMap[fUnpackParSts->GetDpbId(uDpb)] = uDpb;
257 LOG(info) << "Eq. ID for DPB #" << std::setw(2) << uDpb << " = 0x" << std::setw(4) << std::hex
258 << fUnpackParSts->GetDpbId(uDpb) << std::dec << " => " << fDpbIdIndexMap[fUnpackParSts->GetDpbId(uDpb)];
259 } // for( UInt_t uDpb = 0; uDpb < fuNrOfDpbs; ++uDpb )
260
261 fuNbFebs = fUnpackParSts->GetNrOfFebs();
262 LOG(info) << "Nr. of FEBs: " << fuNbFebs;
263
264 fuNbStsXyters = fUnpackParSts->GetNrOfAsics();
265 LOG(info) << "Nr. of StsXyter ASICs: " << fuNbStsXyters;
266
270 fviFebType.resize(fuNrOfDpbs);
273 for (UInt_t uDpb = 0; uDpb < fuNrOfDpbs; ++uDpb) {
274 fvbCrobActiveFlag[uDpb].resize(fUnpackParSts->GetNbCrobsPerDpb());
275 fviFebModuleIdx[uDpb].resize(fUnpackParSts->GetNbCrobsPerDpb());
276 fviFebModuleSide[uDpb].resize(fUnpackParSts->GetNbCrobsPerDpb());
277 fviFebType[uDpb].resize(fUnpackParSts->GetNbCrobsPerDpb());
278 fvdFebAdcGain[uDpb].resize(fUnpackParSts->GetNbCrobsPerDpb());
279 fvdFebAdcOffs[uDpb].resize(fUnpackParSts->GetNbCrobsPerDpb());
280 for (UInt_t uCrobIdx = 0; uCrobIdx < fUnpackParSts->GetNbCrobsPerDpb(); ++uCrobIdx) {
281 fvbCrobActiveFlag[uDpb][uCrobIdx] = fUnpackParSts->IsCrobActive(uDpb, uCrobIdx);
282
283 fviFebModuleIdx[uDpb][uCrobIdx].resize(fUnpackParSts->GetNbFebsPerCrob());
284 fviFebModuleSide[uDpb][uCrobIdx].resize(fUnpackParSts->GetNbFebsPerCrob());
285 fviFebType[uDpb][uCrobIdx].resize(fUnpackParSts->GetNbFebsPerCrob(), -1);
286 fvdFebAdcGain[uDpb][uCrobIdx].resize(fUnpackParSts->GetNbFebsPerCrob(), 0.0);
287 fvdFebAdcOffs[uDpb][uCrobIdx].resize(fUnpackParSts->GetNbFebsPerCrob(), 0.0);
288 for (UInt_t uFebIdx = 0; uFebIdx < fUnpackParSts->GetNbFebsPerCrob(); ++uFebIdx) {
289 fviFebModuleIdx[uDpb][uCrobIdx][uFebIdx] = fUnpackParSts->GetFebModuleIdx(uDpb, uCrobIdx, uFebIdx);
290 fviFebModuleSide[uDpb][uCrobIdx][uFebIdx] = fUnpackParSts->GetFebModuleSide(uDpb, uCrobIdx, uFebIdx);
291 fvdFebAdcGain[uDpb][uCrobIdx][uFebIdx] = fUnpackParSts->GetFebAdcGain(uDpb, uCrobIdx, uFebIdx);
292 fvdFebAdcOffs[uDpb][uCrobIdx][uFebIdx] = fUnpackParSts->GetFebAdcOffset(uDpb, uCrobIdx, uFebIdx);
293
294 if (0 <= fviFebModuleSide[uDpb][uCrobIdx][uFebIdx]
295 && static_cast<UInt_t>(fviFebModuleIdx[uDpb][uCrobIdx][uFebIdx]) < fuNbModules) {
296 switch (fviModuleType[fviFebModuleIdx[uDpb][uCrobIdx][uFebIdx]]) {
297 case 0: // FEB-8-1 with ZIF connector on the right
298 {
299 // P side (0) has type A (0)
300 // N side (1) has type B (1)
301 fviFebType[uDpb][uCrobIdx][uFebIdx] = fviFebModuleSide[uDpb][uCrobIdx][uFebIdx];
302 break;
303 } // case 0: // FEB-8-1 with ZIF connector on the right
304 case 1: // FEB-8-1 with ZIF connector on the left
305 {
306 // P side (0) has type B (1)
307 // N side (1) has type A (0)
308 fviFebType[uDpb][uCrobIdx][uFebIdx] = !(fviFebModuleSide[uDpb][uCrobIdx][uFebIdx]);
309 break;
310 } // case 1: // FEB-8-1 with ZIF connector on the left
311 default: break;
312 } // switch( fviModuleType[ fviFebModuleIdx[ uDpb ][ uCrobIdx ][ uFebIdx ] ] )
313 } // FEB active and module index OK
314 } // for( UInt_t uFebIdx = 0; uFebIdx < fUnpackParSts->GetNbFebsPerCrob(); ++ uFebIdx )
315 } // for( UInt_t uCrobIdx = 0; uCrobIdx < fUnpackParSts->GetNbCrobsPerDpb(); ++uCrobIdx )
316 } // for( UInt_t uDpb = 0; uDpb < fuNrOfDpbs; ++uDpb )
317
318 for (UInt_t uDpb = 0; uDpb < fuNrOfDpbs; ++uDpb) {
319 for (UInt_t uCrobIdx = 0; uCrobIdx < fUnpackParSts->GetNbCrobsPerDpb(); ++uCrobIdx) {
320 LOG(info) << Form("DPB #%02u CROB #%02u Active: ", uDpb, uCrobIdx) << fvbCrobActiveFlag[uDpb][uCrobIdx];
321 } // for( UInt_t uCrobIdx = 0; uCrobIdx < fUnpackParSts->GetNbCrobsPerDpb(); ++uCrobIdx )
322 } // for( UInt_t uDpb = 0; uDpb < fuNrOfDpbs; ++uDpb )
323
324 for (UInt_t uDpb = 0; uDpb < fuNrOfDpbs; ++uDpb) {
325 for (UInt_t uCrobIdx = 0; uCrobIdx < fUnpackParSts->GetNbCrobsPerDpb(); ++uCrobIdx) {
326 LOG(info) << Form("DPB #%02u CROB #%u: ", uDpb, uCrobIdx);
327 for (UInt_t uFebIdx = 0; uFebIdx < fUnpackParSts->GetNbFebsPerCrob(); ++uFebIdx)
328 if (0 <= fviFebModuleIdx[uDpb][uCrobIdx][uFebIdx]) {
329 LOG(info) << Form(" FEB #%02u: Mod. Idx = %03d Side %c (%2d) Type %c (%2d) ADC "
330 "gain %4.0f e- ADC Offs %5.0f e-",
331 uFebIdx, fviFebModuleIdx[uDpb][uCrobIdx][uFebIdx],
332 1 == fviFebModuleSide[uDpb][uCrobIdx][uFebIdx] ? 'N' : 'P',
333 fviFebModuleSide[uDpb][uCrobIdx][uFebIdx],
334 1 == fviFebType[uDpb][uCrobIdx][uFebIdx] ? 'B' : 'A', fviFebType[uDpb][uCrobIdx][uFebIdx],
335 fvdFebAdcGain[uDpb][uCrobIdx][uFebIdx], fvdFebAdcOffs[uDpb][uCrobIdx][uFebIdx]);
336 } // for( UInt_t uFebIdx = 0; uFebIdx < fUnpackParSts->GetNbFebsPerCrob(); ++ uFebIdx )
337 } // for( UInt_t uCrobIdx = 0; uCrobIdx < fUnpackParSts->GetNbCrobsPerDpb(); ++uCrobIdx )
338 } // for( UInt_t uDpb = 0; uDpb < fuNrOfDpbs; ++uDpb )
339
340 if (fbBinningFw) LOG(info) << "Unpacking data in bin sorter FW mode";
341 else
342 LOG(info) << "Unpacking data in full time sorter FW mode (legacy)";
343
344 // Internal status initialization
350 for (UInt_t uDpb = 0; uDpb < fuNrOfDpbs; ++uDpb) {
351 fvulCurrentTsMsb[uDpb] = 0;
352 fvuCurrentTsMsbCycle[uDpb] = 0;
353 fvuInitialHeaderDone[uDpb] = kFALSE;
355 } // for( UInt_t uDpb = 0; uDpb < fuNrOfDpbs; ++uDpb )
356
358
366
367 for (UInt_t uXyterIdx = 0; uXyterIdx < fuNbStsXyters; ++uXyterIdx) {
368 fvulChanLastHitTime[uXyterIdx].resize(fUnpackParSts->GetNbChanPerAsic());
369 fvdChanLastHitTime[uXyterIdx].resize(fUnpackParSts->GetNbChanPerAsic());
370 fvuChanNbHitsInMs[uXyterIdx].resize(fUnpackParSts->GetNbChanPerAsic());
371 fvdChanLastHitTimeInMs[uXyterIdx].resize(fUnpackParSts->GetNbChanPerAsic());
372 fvusChanLastHitAdcInMs[uXyterIdx].resize(fUnpackParSts->GetNbChanPerAsic());
373 fvmAsicHitsInMs[uXyterIdx].clear();
374
375 for (UInt_t uChan = 0; uChan < fUnpackParSts->GetNbChanPerAsic(); ++uChan) {
376 fvulChanLastHitTime[uXyterIdx][uChan] = 0;
377 fvdChanLastHitTime[uXyterIdx][uChan] = -1.0;
378
379 fvuChanNbHitsInMs[uXyterIdx][uChan].resize(fuMaxNbMicroslices);
380 fvdChanLastHitTimeInMs[uXyterIdx][uChan].resize(fuMaxNbMicroslices);
381 fvusChanLastHitAdcInMs[uXyterIdx][uChan].resize(fuMaxNbMicroslices);
382 for (UInt_t uMsIdx = 0; uMsIdx < fuMaxNbMicroslices; ++uMsIdx) {
383 fvuChanNbHitsInMs[uXyterIdx][uChan][uMsIdx] = 0;
384 fvdChanLastHitTimeInMs[uXyterIdx][uChan][uMsIdx] = -1.0;
385 fvusChanLastHitAdcInMs[uXyterIdx][uChan][uMsIdx] = 0;
386 } // for( UInt_t uMsIdx = 0; uMsIdx < fuMaxNbMicroslices; ++uMsIdx )
387 } // for( UInt_t uChan = 0; uChan < fUnpackParSts->GetNbChanPerAsic(); ++uChan )
388 } // for( UInt_t uXyterIdx = 0; uXyterIdx < fuNbStsXyters; ++uXyterIdx )
389
390 LOG(info) << "CbmMcbm2018MonitorSts::ReInitContainers => Changed "
391 "fvuChanNbHitsInMs size "
392 << fvuChanNbHitsInMs.size() << " VS " << fuNbStsXyters;
393 LOG(info) << "CbmMcbm2018MonitorSts::ReInitContainers => Changed "
394 "fvuChanNbHitsInMs size "
395 << fvuChanNbHitsInMs[0].size() << " VS " << fUnpackParSts->GetNbChanPerAsic();
396 LOG(info) << "CbmMcbm2018MonitorSts::ReInitContainers => Changed "
397 "fvuChanNbHitsInMs size "
398 << fvuChanNbHitsInMs[0][0].size() << " VS " << fuMaxNbMicroslices;
399
400 fvmFebHitsInMs.resize(fuNbFebs);
405 for (UInt_t uFebIdx = 0; uFebIdx < fuNbFebs; ++uFebIdx) {
406 fvmFebHitsInMs[uFebIdx].clear();
407 fvdFebChanCountsSinceLastRateUpdate[uFebIdx].resize(fUnpackParSts->GetNbChanPerFeb(), 0.0);
408 fdStsFebChanLastTimeForDist[uFebIdx].resize(fUnpackParSts->GetNbChanPerFeb(), -1.0);
409 } // for( UInt_t uFebIdx = 0; uFebIdx < fuNbFebs; ++uFebIdx )
410
414
415 return kTRUE;
416}
417
418void CbmMcbm2018MonitorSts::AddMsComponentToList(size_t component, UShort_t /*usDetectorId*/)
419{
421 for (UInt_t uCompIdx = 0; uCompIdx < fvMsComponentsList.size(); ++uCompIdx)
422 if (component == fvMsComponentsList[uCompIdx]) return;
423
425 if (kiMaxNbFlibLinks <= component) {
426 LOG(error) << "CbmMcbm2018MonitorSts::AddMsComponentToList => "
427 << "Ignored the addition of component " << component << " as it is above the hadcoded limit of "
428 << static_cast<Int_t>(kiMaxNbFlibLinks) << " !!!!!!!!! "
429 << "\n"
430 << " To change this behavior check kiMaxNbFlibLinks in "
431 "CbmMcbm2018MonitorSts.cxx";
432 return;
433 } // if( kiMaxNbFlibLinks <= component )
434
435
437 fvMsComponentsList.push_back(component);
438 LOG(info) << "CbmMcbm2018MonitorSts::AddMsComponentToList => Added component: " << component;
439
441 if (NULL == fhMsSz[component]) {
442 TString sMsSzName = Form("MsSz_link_%02lu", component);
443 TString sMsSzTitle = Form("Size of MS for nDPB of link %02lu; Ms Size [bytes]", component);
444 fhMsSz[component] = new TH1F(sMsSzName.Data(), sMsSzTitle.Data(), 30000, 0., 30000.);
445 fHM->Add(sMsSzName.Data(), fhMsSz[component]);
446
447 sMsSzName = Form("MsSzTime_link_%02lu", component);
448 sMsSzTitle = Form("Size of MS vs time for gDPB of link %02lu; Time[s] ; Ms Size [bytes]", component);
449 fhMsSzTime[component] = new TProfile(sMsSzName.Data(), sMsSzTitle.Data(), 15000, 0., 300.);
450 fHM->Add(sMsSzName.Data(), fhMsSzTime[component]);
451
452 if (NULL != fcMsSizeAll) {
453 fcMsSizeAll->cd(1 + component);
454 gPad->SetLogy();
455 fhMsSzTime[component]->Draw("hist le0");
456 } // if( NULL != fcMsSizeAll )
457 LOG(info) << "Added MS size histo for component: " << component << " (DPB)";
458
459 THttpServer* server = FairRunOnline::Instance()->GetHttpServer();
460 if (server) {
461 server->Register("/FlibRaw", fhMsSz[component]);
462 server->Register("/FlibRaw", fhMsSzTime[component]);
463 } // if( server )
464 } // if( NULL == fhMsSz[ component ] )
465}
466void CbmMcbm2018MonitorSts::SetNbMsInTs(size_t uCoreMsNb, size_t uOverlapMsNb)
467{
468 fuNbCoreMsPerTs = uCoreMsNb;
469 fuNbOverMsPerTs = uOverlapMsNb;
470
471 UInt_t uNbMsTotal = fuNbCoreMsPerTs + fuNbOverMsPerTs;
472
473 if (fuMaxNbMicroslices < uNbMsTotal) {
474 fuMaxNbMicroslices = uNbMsTotal;
475
480 for (UInt_t uXyterIdx = 0; uXyterIdx < fuNbStsXyters; ++uXyterIdx) {
481 fvuChanNbHitsInMs[uXyterIdx].resize(fUnpackParSts->GetNbChanPerAsic());
482 fvdChanLastHitTimeInMs[uXyterIdx].resize(fUnpackParSts->GetNbChanPerAsic());
483 fvusChanLastHitAdcInMs[uXyterIdx].resize(fUnpackParSts->GetNbChanPerAsic());
484 for (UInt_t uChan = 0; uChan < fUnpackParSts->GetNbChanPerAsic(); ++uChan) {
485 fvuChanNbHitsInMs[uXyterIdx][uChan].resize(fuMaxNbMicroslices);
486 fvdChanLastHitTimeInMs[uXyterIdx][uChan].resize(fuMaxNbMicroslices);
487 fvusChanLastHitAdcInMs[uXyterIdx][uChan].resize(fuMaxNbMicroslices);
488 for (UInt_t uMsIdx = 0; uMsIdx < fuMaxNbMicroslices; ++uMsIdx) {
489 fvuChanNbHitsInMs[uXyterIdx][uChan][uMsIdx] = 0;
490 fvdChanLastHitTimeInMs[uXyterIdx][uChan][uMsIdx] = -1.0;
491 fvusChanLastHitAdcInMs[uXyterIdx][uChan][uMsIdx] = 0;
492 } // for( UInt_t uMsIdx = 0; uMsIdx < fuMaxNbMicroslices; ++uMsIdx )
493 } // for( UInt_t uChan = 0; uChan < fUnpackParSts->GetNbChanPerAsic(); ++uChan )
494 } // for( UInt_t uXyterIdx = 0; uXyterIdx < fuNbStsXyters; ++uXyterIdx )
495 LOG(info) << "CbmMcbm2018MonitorSts::DoUnpack => Changed fvuChanNbHitsInMs size " << fvuChanNbHitsInMs.size()
496 << " VS " << fuNbStsXyters;
497 LOG(info) << "CbmMcbm2018MonitorSts::DoUnpack => Changed fvuChanNbHitsInMs size " << fvuChanNbHitsInMs[0].size()
498 << " VS " << fUnpackParSts->GetNbChanPerAsic();
499 LOG(info) << "CbmMcbm2018MonitorSts::DoUnpack => Changed fvuChanNbHitsInMs size " << fvuChanNbHitsInMs[0][0].size()
500 << " VS " << fuMaxNbMicroslices;
501 } // if( fuMaxNbMicroslices < uNbMsTotal )
502}
503/*
504void CbmMcbm2018MonitorSts::SetCoincidenceBorder( Double_t dCenterPos, Double_t dBorderVal )
505{
506 fdCoincCenter = dCenterPos;
507 fdCoincBorder = dBorderVal;
508 fdCoincMin = dCenterPos - dBorderVal;
509 fdCoincMax = dCenterPos + dBorderVal;
510}
511*/
513{
514 TString sHistName {""};
515 TString title {""};
516
517 sHistName = "hPulserMessageType";
518 title = "Nb of message for each type; Type";
519 fhStsMessType = new TH1I(sHistName, title, 6, 0., 6.);
520 fhStsMessType->GetXaxis()->SetBinLabel(1, "Dummy");
521 fhStsMessType->GetXaxis()->SetBinLabel(2, "Hit");
522 fhStsMessType->GetXaxis()->SetBinLabel(3, "TsMsb");
523 fhStsMessType->GetXaxis()->SetBinLabel(4, "Epoch");
524 fhStsMessType->GetXaxis()->SetBinLabel(5, "Status");
525 fhStsMessType->GetXaxis()->SetBinLabel(6, "Empty");
526
527 /* *** Missing int + MessType OP!!!! ****
528 fhStsMessType->GetXaxis()->SetBinLabel(1 + stsxyter::MessType::Dummy, "Dummy");
529 fhStsMessType->GetXaxis()->SetBinLabel(1 + stsxyter::MessType::Hit, "Hit");
530 fhStsMessType->GetXaxis()->SetBinLabel(1 + stsxyter::MessType::TsMsb, "TsMsb");
531 fhStsMessType->GetXaxis()->SetBinLabel(1 + stsxyter::MessType::ReadDataAck, "ReadDataAck");
532 fhStsMessType->GetXaxis()->SetBinLabel(1 + stsxyter::MessType::Ack, "Ack");
533*/
534
535 sHistName = "hPulserSysMessType";
536 title = "Nb of system message for each type; System Type";
537 fhStsSysMessType = new TH1I(sHistName, title, 17, 0., 17.);
538 /*
539 hSysMessType->GetXaxis()->SetBinLabel(1 + ngdpb::SYSMSG_DAQ_START, "DAQ START");
540 hSysMessType->GetXaxis()->SetBinLabel(1 + ngdpb::SYSMSG_DAQ_FINISH, "DAQ FINISH");
541 hSysMessType->GetXaxis()->SetBinLabel(1 + 16, "GET4 Hack 32B");
542*/
543
544 sHistName = "hPulserMessageTypePerDpb";
545 title = "Nb of message of each type for each DPB; DPB; Type";
546 fhStsMessTypePerDpb = new TH2I(sHistName, title, fuNrOfDpbs, 0, fuNrOfDpbs, 6, 0., 6.);
547 fhStsMessTypePerDpb->GetYaxis()->SetBinLabel(1, "Dummy");
548 fhStsMessTypePerDpb->GetYaxis()->SetBinLabel(2, "Hit");
549 fhStsMessTypePerDpb->GetYaxis()->SetBinLabel(3, "TsMsb");
550 fhStsMessTypePerDpb->GetYaxis()->SetBinLabel(4, "Epoch");
551 fhStsMessTypePerDpb->GetYaxis()->SetBinLabel(5, "Status");
552 fhStsMessTypePerDpb->GetYaxis()->SetBinLabel(6, "Empty");
553 /* *** Missing int + MessType OP!!!! ****
554 fhStsMessType->GetYaxis()->SetBinLabel(1 + stsxyter::MessType::Dummy, "Dummy");
555 fhStsMessType->GetYaxis()->SetBinLabel(1 + stsxyter::MessType::Hit, "Hit");
556 fhStsMessType->GetYaxis()->SetBinLabel(1 + stsxyter::MessType::TsMsb, "TsMsb");
557 fhStsMessType->GetYaxis()->SetBinLabel(1 + stsxyter::MessType::ReadDataAck, "ReadDataAck");
558 fhStsMessType->GetYaxis()->SetBinLabel(1 + stsxyter::MessType::Ack, "Ack");
559*/
560
561 sHistName = "hPulserSysMessTypePerDpb";
562 title = "Nb of system message of each type for each DPB; DPB; System Type";
563 fhStsSysMessTypePerDpb = new TH2I(sHistName, title, fuNrOfDpbs, 0, fuNrOfDpbs, 17, 0., 17.);
564 /*
565 hSysMessType->GetYaxis()->SetBinLabel(1 + ngdpb::SYSMSG_DAQ_START, "DAQ START");
566 hSysMessType->GetYaxis()->SetBinLabel(1 + ngdpb::SYSMSG_DAQ_FINISH, "DAQ FINISH");
567 hSysMessType->GetYaxis()->SetBinLabel(1 + 16, "GET4 Hack 32B");
568*/
569
570 sHistName = "hStsStatusMessType";
571 title = "Nb of status message of each type for each DPB; ASIC; Status Type";
572 fhStsStatusMessType = new TH2I(sHistName, title, fuNbStsXyters, 0, fuNbStsXyters, 16, 0., 16.);
573 /*
574 fhStsStatusMessType->GetYaxis()->SetBinLabel( 1, "Dummy");
575 fhStsStatusMessType->GetYaxis()->SetBinLabel( 2, "Hit");
576 fhStsStatusMessType->GetYaxis()->SetBinLabel( 3, "TsMsb");
577 fhStsStatusMessType->GetYaxis()->SetBinLabel( 4, "Epoch");
578*/
579
580 sHistName = "hStsMsStatusFieldType";
581 title = "For each flag in the MS header, ON/OFF counts; Flag bit []; ON/OFF; MS []";
582 fhStsMsStatusFieldType = new TH2I(sHistName, title, 16, -0.5, 15.5, 2, -0.5, 1.5);
583 /*
584 fhStsMsStatusFieldType->GetYaxis()->SetBinLabel( 1, "Dummy");
585 fhStsMsStatusFieldType->GetYaxis()->SetBinLabel( 2, "Hit");
586 fhStsMsStatusFieldType->GetYaxis()->SetBinLabel( 3, "TsMsb");
587 fhStsMsStatusFieldType->GetYaxis()->SetBinLabel( 4, "Epoch");
588*/
589
590 sHistName = "hStsMessTypePerElink";
591 title = "Nb of message of each type for each DPB; DPB; Type";
592 fhStsMessTypePerElink = new TH2I(sHistName, title, fuNrOfDpbs * fUnpackParSts->GetNbElinkPerDpb(), 0,
593 fuNrOfDpbs * fUnpackParSts->GetNbElinkPerDpb(), 6, 0., 6.);
594 fhStsMessTypePerElink->GetYaxis()->SetBinLabel(1, "Dummy");
595 fhStsMessTypePerElink->GetYaxis()->SetBinLabel(2, "Hit");
596 fhStsMessTypePerElink->GetYaxis()->SetBinLabel(3, "TsMsb");
597 fhStsMessTypePerElink->GetYaxis()->SetBinLabel(4, "Epoch");
598 fhStsMessTypePerElink->GetYaxis()->SetBinLabel(5, "Status");
599 fhStsMessTypePerElink->GetYaxis()->SetBinLabel(6, "Empty");
600
601 sHistName = "hStsHitsElinkPerDpb";
602 title = "Nb of hit messages per eLink for each DPB; DPB; eLink; Hits nb []";
603 fhStsHitsElinkPerDpb = new TH2I(sHistName, title, fuNrOfDpbs, 0, fuNrOfDpbs, 42, 0., 42.);
604
606 sHistName = "hStsAllFebsHitRateEvo";
607 title = "Hits per second & FEB; Time [s]; FEB []; Hits []";
608 fhStsAllFebsHitRateEvo = new TH2I(sHistName, title, 1800, 0, 1800, fuNbFebs, -0.5, fuNbFebs - 0.5);
609
611 sHistName = "hStsAllAsicsHitRateEvo";
612 title = "Hits per second & ASIC; Time [s]; ASIC []; Hits []";
613 fhStsAllAsicsHitRateEvo = new TH2I(sHistName, title, 1800, 0, 1800, fuNbStsXyters, -0.5, fuNbStsXyters - 0.5);
614
616 sHistName = "hStsFebAsicHitCounts";
617 title = "Hits per FEB & ASIC; FEB []; ASIC in FEB[]; Hits []";
618 fhStsFebAsicHitCounts = new TH2I(sHistName, title, fuNbFebs, -0.5, fuNbFebs - 0.5, fUnpackParSts->GetNbAsicsPerFeb(),
619 -0.5, fUnpackParSts->GetNbAsicsPerFeb() - 0.5);
620 /*
621 // Number of rate bins =
622 // 9 for the sub-unit decade
623 // + 9 for each unit of each decade * 10 for the subdecade range
624 // + 1 for the closing bin top edge
625 const Int_t iNbDecadesRate = 9;
626 const Int_t iNbStepsDecade = 9;
627 const Int_t iNbSubStepsInStep = 10;
628 const Int_t iNbBinsRate = iNbStepsDecade
629 + iNbStepsDecade * iNbSubStepsInStep * iNbDecadesRate
630 + 1;
631 Double_t dBinsRate[iNbBinsRate];
632 // First fill sub-unit decade
633 for( Int_t iSubU = 0; iSubU < iNbStepsDecade; iSubU ++ )
634 dBinsRate[ iSubU ] = 0.1 * ( 1 + iSubU );
635 std::cout << std::endl;
636 // Then fill the main decades
637 Double_t dSubstepSize = 1.0 / iNbSubStepsInStep;
638 for( Int_t iDecade = 0; iDecade < iNbDecadesRate; iDecade ++)
639 {
640 Double_t dBase = std::pow( 10, iDecade );
641 Int_t iDecadeIdx = iNbStepsDecade
642 + iDecade * iNbStepsDecade * iNbSubStepsInStep;
643 for( Int_t iStep = 0; iStep < iNbStepsDecade; iStep++ )
644 {
645 Int_t iStepIdx = iDecadeIdx + iStep * iNbSubStepsInStep;
646 for( Int_t iSubStep = 0; iSubStep < iNbSubStepsInStep; iSubStep++ )
647 {
648 dBinsRate[ iStepIdx + iSubStep ] = dBase * (1 + iStep)
649 + dBase * dSubstepSize * iSubStep;
650 } // for( Int_t iSubStep = 0; iSubStep < iNbSubStepsInStep; iSubStep++ )
651 } // for( Int_t iStep = 0; iStep < iNbStepsDecade; iStep++ )
652 } // for( Int_t iDecade = 0; iDecade < iNbDecadesRate; iDecade ++)
653 dBinsRate[ iNbBinsRate - 1 ] = std::pow( 10, iNbDecadesRate );
654*/
657 fuLongHistoBinNb = uAlignedLimit / fuLongHistoBinSizeSec;
658
659 // UInt_t uNbBinEvo = (32768 + 1) * 2;
660 // Double_t dMaxEdgeEvo = stsxyter::kdClockCycleNs
661 // * static_cast< Double_t >( uNbBinEvo ) / 2.0;
662 // Double_t dMinEdgeEvo = dMaxEdgeEvo * -1.0;
663
664 // UInt_t uNbBinDt = static_cast<UInt_t>( (fdCoincMax - fdCoincMin )/stsxyter::kdClockCycleNs );
665
669 // fhStsFebChanDtCoinc.resize( fuNbFebs );
670 // fhStsFebChanCoinc.resize( fuNbFebs );
671 for (UInt_t uFebIdx = 0; uFebIdx < fuNbFebs; ++uFebIdx) {
673 sHistName = Form("hStsFebChanCntRaw_%03u", uFebIdx);
674 title = Form("Hits Count per channel, FEB #%03u; Channel; Hits []", uFebIdx);
675 fhStsFebChanCntRaw.push_back(
676 new TH1I(sHistName, title, fUnpackParSts->GetNbChanPerFeb(), -0.5, fUnpackParSts->GetNbChanPerFeb() - 0.5));
677
678 sHistName = Form("hStsFebChanCntRawGood_%03u", uFebIdx);
679 title = Form("Hits Count per channel in good MS (SX2 bug flag off), FEB "
680 "#%03u; Channel; Hits []",
681 uFebIdx);
682 fhStsFebChanCntRawGood.push_back(
683 new TH1I(sHistName, title, fUnpackParSts->GetNbChanPerFeb(), -0.5, fUnpackParSts->GetNbChanPerFeb() - 0.5));
684
686 sHistName = Form("hStsFebChanAdcRaw_%03u", uFebIdx);
687 title = Form("Raw Adc distribution per channel, FEB #%03u; Channel []; Adc "
688 "[]; Hits []",
689 uFebIdx);
690 fhStsFebChanAdcRaw.push_back(new TH2I(sHistName, title, fUnpackParSts->GetNbChanPerFeb(), -0.5,
691 fUnpackParSts->GetNbChanPerFeb() - 0.5, stsxyter::kuHitNbAdcBins, -0.5,
693
695 sHistName = Form("hStsFebChanAdcRawProfc_%03u", uFebIdx);
696 title = Form("Raw Adc prodile per channel, FEB #%03u; Channel []; Adc []", uFebIdx);
697 fhStsFebChanAdcRawProf.push_back(
698 new TProfile(sHistName, title, fUnpackParSts->GetNbChanPerFeb(), -0.5, fUnpackParSts->GetNbChanPerFeb() - 0.5));
699 /*
701 sHistName = Form( "hStsFebChanAdcCal_%03u", uFebIdx );
702 title = Form( "Cal. Adc distribution per channel, FEB #%03u; Channel []; Adc [e-]; Hits []", uFebIdx );
703 fhStsFebChanAdcCal.push_back( new TH2I(sHistName, title,
704 fUnpackParSts->GetNbChanPerFeb(), -0.5, fUnpackParSts->GetNbChanPerFeb() - 0.5,
705 50, 0., 100000. ) );
706
708 sHistName = Form( "hStsFebChanAdcCalProfc_%03u", uFebIdx );
709 title = Form( "Cal. Adc prodile per channel, FEB #%03u; Channel []; Adc [e-]", uFebIdx );
710 fhStsFebChanAdcCalProf.push_back( new TProfile(sHistName, title,
711 fUnpackParSts->GetNbChanPerFeb(), -0.5, fUnpackParSts->GetNbChanPerFeb() - 0.5 ) );
712*/
714 sHistName = Form("hStsFebChanRawTs_%03u", uFebIdx);
715 title = Form("Raw Timestamp distribution per channel, FEB #%03u; Channel "
716 "[]; Ts []; Hits []",
717 uFebIdx);
718 fhStsFebChanRawTs.push_back(new TH2I(sHistName, title, fUnpackParSts->GetNbChanPerFeb(), -0.5,
719 fUnpackParSts->GetNbChanPerFeb() - 0.5, stsxyter::kuHitNbTsBins, -0.5,
721
723 sHistName = Form("hStsFebChanMissEvt_%03u", uFebIdx);
724 title = Form("Missed Event flags per channel, FEB #%03u; Channel []; Miss "
725 "Evt []; Hits []",
726 uFebIdx);
727 fhStsFebChanMissEvt.push_back(new TH2I(sHistName, title, fUnpackParSts->GetNbChanPerFeb(), -0.5,
728 fUnpackParSts->GetNbChanPerFeb() - 0.5, 2, -0.5, 1.5));
729
731 sHistName = Form("hStsFebChanMissEvtEvo_%03u", uFebIdx);
732 title = Form("Missed Evt flags per second & channel in FEB #%03u; Time "
733 "[s]; Channel []; Missed Evt flags []",
734 uFebIdx);
735 fhStsFebChanMissEvtEvo.push_back(new TH2I(sHistName, title, 1800, 0, 1800, fUnpackParSts->GetNbChanPerFeb(), -0.5,
736 fUnpackParSts->GetNbChanPerFeb() - 0.5));
737
739 sHistName = Form("hStsFebAsicMissEvtEvo_%03u", uFebIdx);
740 title = Form("Missed Evt flags per second & StsXyter in FEB #%03u; Time "
741 "[s]; Asic []; Missed Evt flags []",
742 uFebIdx);
743 fhStsFebAsicMissEvtEvo.push_back(new TH2I(sHistName, title, 1800, 0, 1800, fUnpackParSts->GetNbAsicsPerFeb(), -0.5,
744 fUnpackParSts->GetNbAsicsPerFeb() - 0.5));
745
747 sHistName = Form("hStsFebMissEvtEvo_%03u", uFebIdx);
748 title = Form("Missed Evt flags per second & channel in FEB #%03u; Time "
749 "[s]; Missed Evt flags []",
750 uFebIdx);
751 fhStsFebMissEvtEvo.push_back(new TH1I(sHistName, title, 1800, 0, 1800));
752
754 sHistName = Form("hStsFebChanRateEvo_%03u", uFebIdx);
755 title = Form("Hits per second & channel in FEB #%03u; Time [s]; Channel []; Hits []", uFebIdx);
756 fhStsFebChanHitRateEvo.push_back(new TH2I(sHistName, title, 1800, 0, 1800, fUnpackParSts->GetNbChanPerFeb(), -0.5,
757 fUnpackParSts->GetNbChanPerFeb() - 0.5));
758
760 sHistName = Form("hStsFebChanRateProf_%03u", uFebIdx);
761 title = Form("Hits per second for each channel in FEB #%03u; Channel []; Hits/s []", uFebIdx);
762 fhStsFebChanHitRateProf.push_back(
763 new TProfile(sHistName, title, fUnpackParSts->GetNbChanPerFeb(), -0.5, fUnpackParSts->GetNbChanPerFeb() - 0.5));
764
766 sHistName = Form("hStsFebAsicRateEvo_%03u", uFebIdx);
767 title = Form("Hits per second & StsXyter in FEB #%03u; Time [s]; Asic []; Hits []", uFebIdx);
768 fhStsFebAsicHitRateEvo.push_back(new TH2I(sHistName, title, 1800, 0, 1800, fUnpackParSts->GetNbAsicsPerFeb(), -0.5,
769 fUnpackParSts->GetNbAsicsPerFeb() - 0.5));
770
772 sHistName = Form("hStsFebRateEvo_%03u", uFebIdx);
773 title = Form("Hits per second in FEB #%03u; Time [s]; Hits []", uFebIdx);
774 fhStsFebHitRateEvo.push_back(new TH1I(sHistName, title, 1800, 0, 1800));
775
777 sHistName = Form("hStsFebChanRateEvoLong_%03u", uFebIdx);
778 title = Form("Hits per second & channel in FEB #%03u; Time [min]; Channel []; Hits []", uFebIdx);
779 fhStsFebChanHitRateEvoLong.push_back(new TH2D(sHistName, title, fuLongHistoBinNb, -0.5, uAlignedLimit - 0.5,
780 fUnpackParSts->GetNbChanPerFeb(), -0.5,
781 fUnpackParSts->GetNbChanPerFeb() - 0.5));
782
784 sHistName = Form("hStsFebAsicRateEvoLong_%03u", uFebIdx);
785 title = Form("Hits per second & StsXyter in FEB #%03u; Time [min]; Asic []; Hits []", uFebIdx);
786 fhStsFebAsicHitRateEvoLong.push_back(new TH2D(sHistName, title, fuLongHistoBinNb, -0.5, uAlignedLimit - 0.5,
787 fUnpackParSts->GetNbAsicsPerFeb(), -0.5,
788 fUnpackParSts->GetNbAsicsPerFeb() - 0.5));
789
791 sHistName = Form("hStsFebRateEvoLong_%03u", uFebIdx);
792 title = Form("Hits per second in FEB #%03u; Time [min]; Hits []", uFebIdx);
793 fhStsFebHitRateEvoLong.push_back(new TH1D(sHistName, title, fuLongHistoBinNb, -0.5, uAlignedLimit - 0.5));
794
796 sHistName = Form("hStsFebChanDistT_%03u", uFebIdx);
797 title = Form("Time distance between hits on same channel in FEB #%03u; "
798 "Time difference [ns]; Channel []; ",
799 uFebIdx);
800 fhStsFebChanDistT.push_back(new TH2I(sHistName, title, 1000, -0.5, 6250.0 - 0.5, fUnpackParSts->GetNbChanPerFeb(),
801 -0.5, fUnpackParSts->GetNbChanPerFeb() - 0.5));
802
804 sHistName = Form("hStsFebChanCloseHitsCounts_%03u", uFebIdx);
805 title = Form("Hits with too small dt on same channel in FEB #%03u; Channel "
806 "counts []; ",
807 uFebIdx);
808 fhStsFebChanCloseHitsCounts.push_back(new TH2I(sHistName, title, fUnpackParSts->GetNbChanPerFeb(), -0.5,
809 fUnpackParSts->GetNbChanPerFeb() - 0.5, 2, -0.5, 1.5));
810 sHistName = Form("hStsFebChanCloseHitsRatio_%03u", uFebIdx);
811 title = Form("Ratio of Hits with too small dt on same channel in FEB "
812 "#%03u; Ratio []; ",
813 uFebIdx);
815 new TProfile(sHistName, title, fUnpackParSts->GetNbChanPerFeb(), -0.5, fUnpackParSts->GetNbChanPerFeb() - 0.5));
816
817 /*
819 fhStsFebChanDtCoinc[ uFebIdx ].resize( fuNbFebs, nullptr );
820 fhStsFebChanCoinc[ uFebIdx ].resize( fuNbFebs, nullptr );
821 for( UInt_t uFebIdxB = uFebIdx; uFebIdxB < fuNbFebs; ++uFebIdxB )
822 {
823 sHistName = Form( "hStsFebChanDtCoinc_%03u_%03u", uFebIdx, uFebIdxB );
824 title = Form( "Channel coincidences Time diff between FEB #%03u and FEB #%03u; Time difference [ns]",
825 uFebIdx, uFebIdxB );
826 fhStsFebChanDtCoinc[ uFebIdx ][ uFebIdxB ] = new TH1I( sHistName, title, 400, -1250., 1250.);
827
828 sHistName = Form( "hStsFebChanCoinc_%03u_%03u", uFebIdx, uFebIdxB );
829 title = Form( "Channel coincidences between FEB #%03u and FEB #%03u; Channel FEB #%03u []; Channel FEB #%03u []; Coinc. []",
830 uFebIdx, uFebIdxB, uFebIdx, uFebIdxB );
831 fhStsFebChanCoinc[ uFebIdx ][ uFebIdxB ] = new TH2I( sHistName, title,
832 fUnpackParSts->GetNbChanPerFeb(), -0.5, fUnpackParSts->GetNbChanPerFeb() - 0.5,
833 fUnpackParSts->GetNbChanPerFeb(), -0.5, fUnpackParSts->GetNbChanPerFeb() - 0.5 );
834 } // for( UInt_t uFebIdxB = uFebIdx; uFebIdxB < fuNbFebs; ++uFebIdxB )
835*/
836 } // for( UInt_t uFebIdx = 0; uFebIdx < fuNbFebs; ++uFebIdx )
837
839 /*
841 Double_t dSensorMinX = - fUnpackParSts->GetSensorSzX() / 2.0;
842 Double_t dSensorMaxX = fUnpackParSts->GetSensorSzX() / 2.0;
843 Double_t dSensorMinY = - fUnpackParSts->GetSensorSzY() / 2.0;
844 Double_t dSensorMaxY = fUnpackParSts->GetSensorSzY() / 2.0;
845 for( UInt_t uModIdx = 0; uModIdx < fuNbModules; ++ uModIdx )
846 {
848 sHistName = Form( "hStsModulePNCoincDt_%03u", uModIdx );
849 title = Form( "Channel coincidences Time diff between P and N sides on module #%03u; Time difference [ns]",
850 uModIdx );
851 fhStsModulePNCoincDt.push_back( new TH1I( sHistName, title, 400, -1250., 1250.) );
852
853 sHistName = Form( "hStsModulePNCoincDtAsicP_%03u", uModIdx );
854 title = Form( "Channel coincidences Time diff between P and N sides on module #%03u; Time difference [ns]; Asic on P FEB []",
855 uModIdx );
856 fhStsModulePNCoincDtAsicP.push_back( new TH2I( sHistName, title,
857 400, -1250., 1250.,
858 fUnpackParSts->GetNbAsicsPerFeb(), -0.5, fUnpackParSts->GetNbAsicsPerFeb() - 0.5) );
859
860 sHistName = Form( "hStsModulePNCoincDtAsicN_%03u", uModIdx );
861 title = Form( "Channel coincidences Time diff between P and N sides on module #%03u; Time difference [ns]; Asic on P FEB []",
862 uModIdx );
863 fhStsModulePNCoincDtAsicN.push_back( new TH2I( sHistName, title,
864 400, -1250., 1250.,
865 fUnpackParSts->GetNbAsicsPerFeb(), -0.5, fUnpackParSts->GetNbAsicsPerFeb() - 0.5) );
866
867
869 sHistName = Form( "hStsModulePNCoincChan_%03u", uModIdx );
870 title = Form( "P-N channel coincidences in module #%03u; Channel P []; Channel N []; Cnt []", uModIdx );
871 fhStsModulePNCoincChan.push_back( new TH2D( sHistName, title,
872 fUnpackParSts->GetNbChanPerFeb(), -0.5, fUnpackParSts->GetNbChanPerFeb() - 0.5,
873 fUnpackParSts->GetNbChanPerFeb(), -0.5, fUnpackParSts->GetNbChanPerFeb() - 0.5 ) );
874
876 sHistName = Form( "hStsModulePNCoincAdc_%03u", uModIdx );
877 title = Form( "Adc values of P-N coincidences in module #%03u; ADC Channel P [bin]; ADC Channel N [bin]; Cnt []", uModIdx );
878 fhStsModulePNCoincAdc.push_back( new TH2D( sHistName, title,
879 stsxyter::kuHitNbAdcBins, -0.5, stsxyter::kuHitNbAdcBins - 0.5,
880 stsxyter::kuHitNbAdcBins, -0.5, stsxyter::kuHitNbAdcBins - 0.5 ) );
882 sHistName = Form( "hStsModuleCoincAdcChanP_%03u", uModIdx );
883 title = Form( "Adc values of P chan in P-N coincidences in module #%03u; Channel P [bin]; ADC val. [bin]; Cnt []", uModIdx );
884 fhStsModuleCoincAdcChanP.push_back( new TH2D( sHistName, title,
885 fUnpackParSts->GetNbChanPerFeb(), -0.5, fUnpackParSts->GetNbChanPerFeb() - 0.5,
886 stsxyter::kuHitNbAdcBins, -0.5, stsxyter::kuHitNbAdcBins - 0.5 ) );
888 sHistName = Form( "hStsModuleCoincAdcChanN_%03u", uModIdx );
889 title = Form( "Adc values of N chan in P-N coincidences in module #%03u; Channel N [bin]; ADC val. [bin]; Cnt []", uModIdx );
890 fhStsModuleCoincAdcChanN.push_back( new TH2D( sHistName, title,
891 fUnpackParSts->GetNbChanPerFeb(), -0.5, fUnpackParSts->GetNbChanPerFeb() - 0.5,
892 stsxyter::kuHitNbAdcBins, -0.5, stsxyter::kuHitNbAdcBins - 0.5 ) );
893
895 sHistName = Form( "hStsModuleCoincMap_%03u", uModIdx );
896 title = Form( "X-Y map of P-N coincidences in module #%03u; Pos. X [mm]; Pos. Y [mm]; Cnt []", uModIdx );
897 fhStsModuleCoincMap.push_back( new TH2D( sHistName, title,
898 2*fUnpackParSts->GetSensorSzX(), 2*dSensorMinX, 2*dSensorMaxX,
899 2*fUnpackParSts->GetSensorSzY(), 2*dSensorMinY, 2*dSensorMaxY ) );
900 } // for( UInt_t uModIdx = 0; uModIdx < fuNbModules; ++ uModIdx )
901*/
903
906 if (kTRUE == fbEnableCheckBugSmx20) {
907 for (UInt_t uFebIdx = 0; uFebIdx < fuNbFebs; ++uFebIdx) {
909 sHistName = Form("hStsFebSmxErrRatioEvo_%03u", uFebIdx);
910 title = Form("Proportion of uS with SMX logic error in FEB #%03u; Time "
911 "[s]; Error MS fract. []",
912 uFebIdx);
913 fhStsFebSmxErrRatioEvo.push_back(new TProfile(sHistName, title, 1800, 0, 1800));
914
916 sHistName = Form("hStsFebSmxErrRatioEvoAsic_%03u", uFebIdx);
917 title = Form("Proportion of uS with SMX logic error per ASIC in FEB "
918 "#%03u; Time [s]; ASIC []; Error MS fract. []",
919 uFebIdx);
920 fhStsFebSmxErrRatioEvoAsic.push_back(new TProfile2D(sHistName, title, 1800, 0, 1800,
921 fUnpackParSts->GetNbAsicsPerFeb(), -0.5,
922 fUnpackParSts->GetNbAsicsPerFeb() - 0.5));
923
925 sHistName = Form("hStsFebSmxErrRatioCopyEvo_%03u", uFebIdx);
926 title = Form("Proportion of uS with hit copies in FEB #%03u; Time [s]; "
927 "Copies MS fract. []",
928 uFebIdx);
929 fhStsFebSmxErrRatioCopyEvo.push_back(new TProfile(sHistName, title, 1800, 0, 1800));
930
932 sHistName = Form("hStsFebSmxErrRatioCopyEvoAsic_%03u", uFebIdx);
933 title = Form("Proportion of uS with hit copies per ASIC in FEB #%03u; "
934 "Time [s]; ASIC []; Copies MS fract. []",
935 uFebIdx);
936 fhStsFebSmxErrRatioCopyEvoAsic.push_back(new TProfile2D(sHistName, title, 1800, 0, 1800,
937 fUnpackParSts->GetNbAsicsPerFeb(), -0.5,
938 fUnpackParSts->GetNbAsicsPerFeb() - 0.5));
939
941 sHistName = Form("hStsFebSmxErrRatioCopySameAdcEvo_%03u", uFebIdx);
942 title = Form("Proportion of uS with hit full copies in FEB #%03u; Time "
943 "[s]; Copies MS fract. []",
944 uFebIdx);
945 fhStsFebSmxErrRatioCopySameAdcEvo.push_back(new TProfile(sHistName, title, 1800, 0, 1800));
946
948 sHistName = Form("hStsFebSmxErrRatioCopySameAdcEvoAsic_%03u", uFebIdx);
949 title = Form("Proportion of uS with hit full copies per ASIC in FEB "
950 "#%03u; Time [s]; ASIC []; Copies MS fract. []",
951 uFebIdx);
952 fhStsFebSmxErrRatioCopySameAdcEvoAsic.push_back(new TProfile2D(sHistName, title, 1800, 0, 1800,
953 fUnpackParSts->GetNbAsicsPerFeb(), -0.5,
954 fUnpackParSts->GetNbAsicsPerFeb() - 0.5));
955 } // for( UInt_t uFebIdx = 0; uFebIdx < fuNbFebs; ++uFebIdx )
956 } // if( kTRUE == fbEnableCheckBugSmx20 )
958
959 fhMsErrorsEvo = new TH2I("fhMsErrorsEvo", "; MS index [s]; Error type []; Counts []", 600, 0.0, 600.0, 4, -0.5, 3.5);
961
962 // Miscroslice properties histos
963 for (Int_t component = 0; component < kiMaxNbFlibLinks; component++) {
964 fhMsSz[component] = NULL;
965 fhMsSzTime[component] = NULL;
966 } // for( Int_t component = 0; component < kiMaxNbFlibLinks; component ++ )
967
968 // Online histo browser commands
969 THttpServer* server = FairRunOnline::Instance()->GetHttpServer();
970 if (server) {
971 server->Register("/StsRaw", fhStsMessType);
972 server->Register("/StsRaw", fhStsSysMessType);
973 server->Register("/StsRaw", fhStsMessTypePerDpb);
974 server->Register("/StsRaw", fhStsSysMessTypePerDpb);
975 server->Register("/StsRaw", fhStsStatusMessType);
976 server->Register("/StsRaw", fhStsMsStatusFieldType);
977 server->Register("/StsRaw", fhStsMessTypePerElink);
978 server->Register("/StsRaw", fhStsHitsElinkPerDpb);
979 server->Register("/StsRaw", fhStsAllFebsHitRateEvo);
980 server->Register("/StsRaw", fhStsAllAsicsHitRateEvo);
981 server->Register("/StsRaw", fhStsFebAsicHitCounts);
982
983 for (UInt_t uFebIdx = 0; uFebIdx < fuNbFebs; ++uFebIdx) {
984 if (kTRUE == fUnpackParSts->IsFebActive(uFebIdx)) {
985 server->Register("/StsFeb", fhStsFebChanCntRaw[uFebIdx]);
986 server->Register("/StsFeb", fhStsFebChanCntRawGood[uFebIdx]);
987 server->Register("/StsFeb", fhStsFebChanAdcRaw[uFebIdx]);
988 server->Register("/StsFeb", fhStsFebChanAdcRawProf[uFebIdx]);
989 // server->Register("/StsFeb", fhStsFebChanAdcCal[ uFebIdx ] );
990 // server->Register("/StsFeb", fhStsFebChanAdcCalProf[ uFebIdx ] );
991 server->Register("/StsFeb", fhStsFebChanRawTs[uFebIdx]);
992 server->Register("/StsFeb", fhStsFebChanMissEvt[uFebIdx]);
993 server->Register("/StsFeb", fhStsFebChanMissEvtEvo[uFebIdx]);
994 server->Register("/StsFeb", fhStsFebAsicMissEvtEvo[uFebIdx]);
995 server->Register("/StsFeb", fhStsFebMissEvtEvo[uFebIdx]);
996 server->Register("/StsFeb", fhStsFebChanHitRateEvo[uFebIdx]);
997 server->Register("/StsFeb", fhStsFebChanHitRateProf[uFebIdx]);
998 server->Register("/StsFeb", fhStsFebAsicHitRateEvo[uFebIdx]);
999 server->Register("/StsFeb", fhStsFebHitRateEvo[uFebIdx]);
1000 server->Register("/StsFeb", fhStsFebChanHitRateEvoLong[uFebIdx]);
1001 server->Register("/StsFeb", fhStsFebAsicHitRateEvoLong[uFebIdx]);
1002 server->Register("/StsFeb", fhStsFebHitRateEvoLong[uFebIdx]);
1003 server->Register("/StsFeb", fhStsFebChanDistT[uFebIdx]);
1004 server->Register("/StsFeb", fhStsFebChanCloseHitsCounts[uFebIdx]);
1005 server->Register("/StsFeb", fhStsFebChanCloseHitsRatio[uFebIdx]);
1006 /*
1007 for( UInt_t uFebIdxB = uFebIdx; uFebIdxB < fuNbFebs; ++uFebIdxB )
1008 {
1009 server->Register("/StsFeb", fhStsFebChanDtCoinc[ uFebIdx ][ uFebIdxB ] );
1010 server->Register("/StsFeb", fhStsFebChanCoinc[ uFebIdx ][ uFebIdxB ] );
1011 } // for( UInt_t uFebIdxB = uFebIdx; uFebIdxB < fuNbFebs; ++uFebIdxB )
1012*/
1013 if (kTRUE == fbEnableCheckBugSmx20) {
1014 server->Register("/StsSmxErr", fhStsFebSmxErrRatioEvo[uFebIdx]);
1015 server->Register("/StsSmxErr", fhStsFebSmxErrRatioEvoAsic[uFebIdx]);
1016 server->Register("/StsSmxErr", fhStsFebSmxErrRatioCopyEvo[uFebIdx]);
1017 server->Register("/StsSmxErr", fhStsFebSmxErrRatioCopyEvoAsic[uFebIdx]);
1018 server->Register("/StsSmxErr", fhStsFebSmxErrRatioCopySameAdcEvo[uFebIdx]);
1019 server->Register("/StsSmxErr", fhStsFebSmxErrRatioCopySameAdcEvoAsic[uFebIdx]);
1020 } // if( kTRUE == fbEnableCheckBugSmx20 )
1021 } // if( kTRUE == fUnpackParSts->IsFebActive( uFebIdx ) )
1022 } // for( UInt_t uFebIdx = 0; uFebIdx < fuNbFebs; ++uFebIdx )
1023
1024 server->Register("/StsRaw", fhMsErrorsEvo);
1025 /*
1026 for( UInt_t uModIdx = 0; uModIdx < fuNbModules; ++ uModIdx )
1027 {
1028 server->Register("/StsMod", fhStsModulePNCoincDt[ uModIdx ] );
1029 server->Register("/StsMod", fhStsModulePNCoincDtAsicP[ uModIdx ] );
1030 server->Register("/StsMod", fhStsModulePNCoincDtAsicN[ uModIdx ] );
1031 server->Register("/StsMod", fhStsModulePNCoincChan[ uModIdx ] );
1032 server->Register("/StsMod", fhStsModulePNCoincAdc[ uModIdx ] );
1033 server->Register("/StsMod", fhStsModuleCoincAdcChanP[ uModIdx ] );
1034 server->Register("/StsMod", fhStsModuleCoincAdcChanN[ uModIdx ] );
1035 server->Register("/StsMod", fhStsModuleCoincMap[ uModIdx ] );
1036 } // for( UInt_t uModIdx = 0; uModIdx < fuNbModules; ++ uModIdx )
1037*/
1038
1039 server->RegisterCommand("/Reset_All", "bMcbm2018ResetSts=kTRUE");
1040 server->RegisterCommand("/Write_All", "bMcbm2018WriteSts=kTRUE");
1041 server->RegisterCommand("/ScanNoisyCh", "bMcbm2018ScanNoisySts=kTRUE");
1042
1043
1044 server->Restrict("/Reset_All", "allow=admin");
1045 server->Restrict("/Write_All", "allow=admin");
1046 server->Restrict("/ScanNoisyCh", "allow=admin");
1047 } // if( server )
1048
1050 Double_t w = 10;
1051 Double_t h = 10;
1052
1053 // Summary per FEB
1054 fvcStsSumm.resize(fuNbFebs);
1055 fvcStsSmxErr.resize(fuNbFebs);
1056 for (UInt_t uFebIdx = 0; uFebIdx < fuNbFebs; ++uFebIdx) {
1057 if (kTRUE == fUnpackParSts->IsFebActive(uFebIdx)) {
1058 fvcStsSumm[uFebIdx] =
1059 new TCanvas(Form("cStsSum_%03u", uFebIdx), Form("Summary plots for FEB %03u", uFebIdx), w, h);
1060 fvcStsSumm[uFebIdx]->Divide(2, 3);
1061
1062 fvcStsSumm[uFebIdx]->cd(1);
1063 gPad->SetGridx();
1064 gPad->SetGridy();
1065 gPad->SetLogy();
1066 fhStsFebChanCntRaw[uFebIdx]->Draw();
1067
1068 fvcStsSumm[uFebIdx]->cd(2);
1069 gPad->SetGridx();
1070 gPad->SetGridy();
1071 gPad->SetLogy();
1072 fhStsFebChanHitRateProf[uFebIdx]->Draw("e0");
1073
1074 fvcStsSumm[uFebIdx]->cd(3);
1075 gPad->SetGridx();
1076 gPad->SetGridy();
1077 gPad->SetLogz();
1078 fhStsFebChanAdcRaw[uFebIdx]->Draw("colz");
1079
1080 fvcStsSumm[uFebIdx]->cd(4);
1081 gPad->SetGridx();
1082 gPad->SetGridy();
1083 // gPad->SetLogy();
1084 fhStsFebChanAdcRawProf[uFebIdx]->Draw();
1085
1086 fvcStsSumm[uFebIdx]->cd(5);
1087 gPad->SetGridx();
1088 gPad->SetGridy();
1089 gPad->SetLogz();
1090 fhStsFebChanHitRateEvo[uFebIdx]->Draw("colz");
1091
1092 fvcStsSumm[uFebIdx]->cd(6);
1093 gPad->SetGridx();
1094 gPad->SetGridy();
1095 // gPad->SetLogy();
1096 fhStsFebChanMissEvt[uFebIdx]->Draw("colz");
1097 /*
1098 fvcStsSumm[ uFebIdx ]->cd(5);
1099 gPad->SetGridx();
1100 gPad->SetGridy();
1101 gPad->SetLogz();
1102 fhStsFebChanAdcCal[ uFebIdx ]->Draw( "colz" );
1103
1104 fvcStsSumm[ uFebIdx ]->cd(6);
1105 gPad->SetGridx();
1106 gPad->SetGridy();
1107 // gPad->SetLogy();
1108 fhStsFebChanAdcCalProf[ uFebIdx ]->Draw();
1109*/
1110 server->Register("/canvases", fvcStsSumm[uFebIdx]);
1111
1112 if (kTRUE == fbEnableCheckBugSmx20) {
1113 fvcStsSmxErr[uFebIdx] =
1114 new TCanvas(Form("cStsSmxErr_%03u", uFebIdx), Form("SMX logic error plots for FEB %03u", uFebIdx), w, h);
1115 fvcStsSmxErr[uFebIdx]->Divide(2, 3);
1116
1117 fvcStsSmxErr[uFebIdx]->cd(1);
1118 gPad->SetGridx();
1119 gPad->SetGridy();
1120 gPad->SetLogy();
1121 fhStsFebSmxErrRatioEvo[uFebIdx]->Draw();
1122
1123 fvcStsSmxErr[uFebIdx]->cd(2);
1124 gPad->SetGridx();
1125 gPad->SetGridy();
1126 gPad->SetLogz();
1127 fhStsFebSmxErrRatioEvoAsic[uFebIdx]->Draw("colz");
1128
1129 fvcStsSmxErr[uFebIdx]->cd(3);
1130 gPad->SetGridx();
1131 gPad->SetGridy();
1132 gPad->SetLogy();
1133 fhStsFebSmxErrRatioCopyEvo[uFebIdx]->Draw();
1134
1135 fvcStsSmxErr[uFebIdx]->cd(4);
1136 gPad->SetGridx();
1137 gPad->SetGridy();
1138 gPad->SetLogz();
1139 fhStsFebSmxErrRatioCopyEvoAsic[uFebIdx]->Draw("colz");
1140
1141 fvcStsSmxErr[uFebIdx]->cd(5);
1142 gPad->SetGridx();
1143 gPad->SetGridy();
1144 gPad->SetLogy();
1145 fhStsFebSmxErrRatioCopySameAdcEvo[uFebIdx]->Draw();
1146
1147 fvcStsSmxErr[uFebIdx]->cd(6);
1148 gPad->SetGridx();
1149 gPad->SetGridy();
1150 gPad->SetLogz();
1151 fhStsFebSmxErrRatioCopySameAdcEvoAsic[uFebIdx]->Draw("colz");
1152
1153 server->Register("/canvases", fvcStsSmxErr[uFebIdx]);
1154 } // if( kTRUE == fbEnableCheckBugSmx20 )
1155 } // if( kTRUE == fUnpackParSts->IsFebActive( uFebIdx ) )
1156 } // for( UInt_t uFebIdx = 0; uFebIdx < fuNbFebs; ++uFebIdx )
1157
1158 //====================================================================//
1159
1160 //====================================================================//
1162 // Try to recover canvas in case it was created already by another monitor
1163 // If not existing, create it
1164 fcMsSizeAll = dynamic_cast<TCanvas*>(gROOT->FindObject("cMsSizeAll"));
1165 if (NULL == fcMsSizeAll) {
1166 fcMsSizeAll = new TCanvas("cMsSizeAll", "Evolution of MS size in last 300 s", w, h);
1167 fcMsSizeAll->Divide(4, 4);
1168 LOG(info) << "Created MS size canvas in STS monitor";
1169 } // if( NULL == fcMsSizeAll )
1170 else
1171 LOG(info) << "Recovered MS size canvas in STS monitor";
1172 server->Register("/canvases", fcMsSizeAll);
1173 //====================================================================//
1174
1175 /*****************************/
1176}
1177
1178Bool_t CbmMcbm2018MonitorSts::DoUnpack(const fles::Timeslice& ts, size_t component)
1179{
1180 if (bMcbm2018ResetSts) {
1182 bMcbm2018ResetSts = kFALSE;
1183 } // if( bMcbm2018ResetSts )
1184 if (bMcbm2018WriteSts) {
1186 bMcbm2018WriteSts = kFALSE;
1187 } // if( bMcbm2018WriteSts )
1190 bMcbm2018ScanNoisySts = kFALSE;
1191 } // if( bMcbm2018WriteSts )
1192
1193 LOG(debug) << "Timeslice contains " << ts.num_microslices(component) << "microslices.";
1194 fulCurrentTsIdx = ts.index();
1195
1197 if (0 == fulCurrentTsIdx) return kTRUE;
1198
1199 // Ignore overlap ms if flag set by user
1200 UInt_t uNbMsLoop = fuNbCoreMsPerTs;
1201 if (kFALSE == fbIgnoreOverlapMs) uNbMsLoop += fuNbOverMsPerTs;
1202
1203 // Loop over core microslices (and overlap ones if chosen)
1204 for (UInt_t uMsIdx = 0; uMsIdx < uNbMsLoop; uMsIdx++) {
1205 Double_t dMsTime = (1e-9) * static_cast<double>(ts.descriptor(fvMsComponentsList[0], uMsIdx).idx);
1206
1207 if (0 == fulCurrentTsIdx && 0 == uMsIdx) {
1208 for (UInt_t uMsCompIdx = 0; uMsCompIdx < fvMsComponentsList.size(); ++uMsCompIdx) {
1209 UInt_t uMsComp = fvMsComponentsList[uMsCompIdx];
1210 auto msDescriptor = ts.descriptor(uMsComp, uMsIdx);
1211 /*
1212 LOG(info) << "hi hv eqid flag si sv idx/start crc size offset";
1213 LOG(info) << Form( "%02x %02x %04x %04x %02x %02x %016llx %08x %08x %016llx",
1214 static_cast<unsigned int>(msDescriptor.hdr_id),
1215 static_cast<unsigned int>(msDescriptor.hdr_ver), msDescriptor.eq_id, msDescriptor.flags,
1216 static_cast<unsigned int>(msDescriptor.sys_id),
1217 static_cast<unsigned int>(msDescriptor.sys_ver), msDescriptor.idx, msDescriptor.crc,
1218 msDescriptor.size, msDescriptor.offset );
1219*/
1220 LOG(info) << FormatMsHeaderPrintout(msDescriptor);
1221 uint32_t uEqId = static_cast<uint32_t>(msDescriptor.eq_id & 0xFFFF);
1222 auto it = fDpbIdIndexMap.find(uEqId);
1223 if (fDpbIdIndexMap.end() == it) {
1224 LOG(warning) << "Could not find the sDPB index for AFCK id 0x" << std::hex << uEqId << std::dec
1225 << " component " << uMsCompIdx << "\n"
1226 << "If valid this index has to be added in the TOF parameter file "
1227 "in the RocIdArray field"
1228 << "\n"
1229 << "For now we remove it from the list of components analyzed";
1230 } // if( fDpbIdIndexMap.end() == it )
1231 } // for( UInt_t uMsCompIdx = 0; uMsCompIdx < fvMsComponentsList.size(); ++uMsCompIdx )
1232 } // if( 0 == fulCurrentTsIndex && 0 == uMsIdx )
1233
1234 // Loop over registered components
1235 for (UInt_t uMsCompIdx = 0; uMsCompIdx < fvMsComponentsList.size(); ++uMsCompIdx) {
1236 UInt_t uMsComp = fvMsComponentsList[uMsCompIdx];
1237
1238 if (kFALSE == ProcessStsMs(ts, uMsComp, uMsIdx)) return kFALSE;
1239 /*
1240 auto msDescriptor = ts.descriptor( uMsComp, uMsIdx );
1241 fuCurrentEquipmentId = msDescriptor.eq_id;
1242 const uint8_t* msContent = reinterpret_cast<const uint8_t*>( ts.content( uMsComp, uMsIdx ) );
1243
1244 uint32_t uSize = msDescriptor.size;
1245 fulCurrentMsIdx = msDescriptor.idx;
1246 dMsTime = (1e-9) * static_cast<double>(fulCurrentMsIdx);
1247 LOG(debug) << "Microslice: " << fulCurrentMsIdx
1248 << " from EqId " << std::hex << fuCurrentEquipmentId << std::dec
1249 << " has size: " << uSize;
1250
1251 fuCurrDpbId = static_cast< uint32_t >( fuCurrentEquipmentId & 0xFFFF );
1252 fuCurrDpbIdx = fDpbIdIndexMap[ fuCurrDpbId ];
1253
1254 if( uMsComp < kiMaxNbFlibLinks )
1255 {
1256 if( fdStartTimeMsSz < 0 )
1257 fdStartTimeMsSz = dMsTime;
1258 fhMsSz[ uMsComp ]->Fill( uSize );
1259 fhMsSzTime[ uMsComp ]->Fill( dMsTime - fdStartTimeMsSz, uSize);
1260 } // if( uMsComp < kiMaxNbFlibLinks )
1261
1263 if( static_cast<Int_t>( fvdMsTime[ uMsCompIdx ] ) < static_cast<Int_t>( dMsTime ) )
1264 {
1266 UInt_t uFebIdxOffset = fUnpackParSts->GetNbFebsPerDpb() * fuCurrDpbIdx;
1267 for( UInt_t uFebIdx = 0; uFebIdx < fUnpackParSts->GetNbFebsPerDpb(); ++uFebIdx )
1268 {
1269 UInt_t uFebIdxInSyst = uFebIdxOffset + uFebIdx;
1270
1272 if( 0 == fviFebTimeSecLastRateUpdate[uFebIdxInSyst] )
1273 {
1274 fviFebTimeSecLastRateUpdate[uFebIdxInSyst] = static_cast<Int_t>( dMsTime );
1275 fviFebCountsSinceLastRateUpdate[uFebIdxInSyst] = 0;
1276 for( UInt_t uChan = 0; uChan < fUnpackParSts->GetNbChanPerFeb(); ++uChan )
1277 fvdFebChanCountsSinceLastRateUpdate[uFebIdxInSyst][uChan] = 0.0;
1278 continue;
1279 } // if( 0 == fviFebTimeSecLastRateUpdate[uFebIdxInSyst] )
1280
1281 Int_t iTimeInt = static_cast<Int_t>( dMsTime ) - fviFebTimeSecLastRateUpdate[uFebIdxInSyst];
1282 if( fiTimeIntervalRateUpdate <= iTimeInt )
1283 {
1285 if( 0 == fviFebCountsSinceLastRateUpdate[uFebIdxInSyst] )
1286 {
1287 fviFebTimeSecLastRateUpdate[uFebIdxInSyst] = dMsTime;
1288 continue;
1289 } // if( 0 == fviFebCountsSinceLastRateUpdate[uFebIdxInSyst] )
1290
1291 for( UInt_t uChan = 0; uChan < fUnpackParSts->GetNbChanPerFeb(); ++uChan )
1292 {
1293 fhStsFebChanHitRateProf[uFebIdxInSyst]->Fill( uChan,
1294 fvdFebChanCountsSinceLastRateUpdate[uFebIdxInSyst][uChan] / iTimeInt );
1295 fvdFebChanCountsSinceLastRateUpdate[uFebIdxInSyst][uChan] = 0.0;
1296 } // for( UInt_t uChan = 0; uChan < fUnpackParSts->GetNbChanPerFeb(); ++uChan )
1297
1298 fviFebTimeSecLastRateUpdate[uFebIdxInSyst] = dMsTime;
1299 fviFebCountsSinceLastRateUpdate[uFebIdxInSyst] = 0;
1300 } // if( fiTimeIntervalRateUpdate <= iTimeInt )
1301 } // for( UInt_t uFebIdx = 0; uFebIdx < fUnpackParSts->GetNbFebsPerDpb(); ++uFebIdx )
1302 } // if( static_cast<Int_t>( fvdMsTime[ uMsCompIdx ] ) < static_cast<Int_t>( dMsTime ) )
1303
1304 // Store MS time for coincidence plots
1305 fvdMsTime[ uMsCompIdx ] = dMsTime;
1306
1308 UInt_t uTsMsbCycleHeader = std::floor( fulCurrentMsIdx /
1309 ( stsxyter::kulTsCycleNbBins * stsxyter::kdClockCycleNs ) )
1310 - fvuInitialTsMsbCycleHeader[ fuCurrDpbIdx ];
1311 if( kFALSE == fvuInitialHeaderDone[ fuCurrDpbIdx ] )
1312 {
1313 fvuInitialTsMsbCycleHeader[ fuCurrDpbIdx ] = uTsMsbCycleHeader;
1314 fvuInitialHeaderDone[ fuCurrDpbIdx ] = kTRUE;
1315 } // if( kFALSE == fvuInitialHeaderDone[ fuCurrDpbIdx ] )
1316 else if( uTsMsbCycleHeader != fvuCurrentTsMsbCycle[ fuCurrDpbIdx ] &&
1317 4194303 != fvulCurrentTsMsb[fuCurrDpbIdx] )
1318 {
1319 LOG(warning) << "TS MSB cycle from MS header does not match current cycle from data "
1320 << "for TS " << std::setw( 12 ) << fulCurrentTsIdx
1321 << " MS " << std::setw( 12 ) << fulCurrentMsIdx
1322 << " MsInTs " << std::setw( 3 ) << uMsIdx
1323 << " ====> " << fvuCurrentTsMsbCycle[ fuCurrDpbIdx ]
1324 << " VS " << uTsMsbCycleHeader;
1325 fvuCurrentTsMsbCycle[ fuCurrDpbIdx ] = uTsMsbCycleHeader;
1326 }
1327
1328 // If not integer number of message in input buffer, print warning/error
1329 if( 0 != ( uSize % kuBytesPerMessage ) )
1330 LOG(error) << "The input microslice buffer does NOT "
1331 << "contain only complete nDPB messages!";
1332
1333 // Compute the number of complete messages in the input microslice buffer
1334 uint32_t uNbMessages = ( uSize - ( uSize % kuBytesPerMessage ) )
1335 / kuBytesPerMessage;
1336
1337 // Prepare variables for the loop on contents
1338 const uint32_t* pInBuff = reinterpret_cast<const uint32_t*>( msContent );
1339
1340 for( uint32_t uIdx = 0; uIdx < uNbMessages; ++uIdx )
1341 {
1342 // Fill message
1343 uint32_t ulData = static_cast<uint32_t>( pInBuff[uIdx] );
1344
1345 stsxyter::Message mess( static_cast< uint32_t >( ulData & 0xFFFFFFFF ) );
1346
1347 // Print message if requested
1348 if( fbPrintMessages )
1349 mess.PrintMess( std::cout, fPrintMessCtrl );
1350
1351 if( 1000 == fulCurrentTsIdx )
1352 {
1353 mess.PrintMess( std::cout, fPrintMessCtrl );
1354 } // if( 0 == fulCurrentTsIdx )
1355
1356 stsxyter::MessType typeMess = mess.GetMessType();
1357 fmMsgCounter[ typeMess ] ++;
1358 fhStsMessType->Fill( static_cast< uint16_t > (typeMess) );
1359 fhStsMessTypePerDpb->Fill( fuCurrDpbIdx, static_cast< uint16_t > (typeMess) );
1360
1361 switch( typeMess )
1362 {
1363 case stsxyter::MessType::Hit :
1364 {
1365 // Extract the eLink and Asic indices => Should GO IN the fill method now that obly hits are link/asic specific!
1366 UShort_t usElinkIdx = mess.GetLinkIndex();
1367 UInt_t uCrobIdx = usElinkIdx / fUnpackParSts->GetNbElinkPerCrob();
1368 Int_t uFebIdx = fUnpackParSts->ElinkIdxToFebIdx( usElinkIdx );
1369 if( -1 == uFebIdx )
1370 {
1371 LOG(warning) << "CbmMcbm2018MonitorSts::DoUnpack => "
1372 << "Wrong elink Idx!";
1373 continue;
1374 } // if( -1 == uFebIdx )
1375
1376 UInt_t uAsicIdx = ( fuCurrDpbIdx * fUnpackParSts->GetNbCrobsPerDpb() + uCrobIdx
1377 ) * fUnpackParSts->GetNbAsicsPerCrob()
1378 + fUnpackParSts->ElinkIdxToAsicIdx( 1 == fviFebType[ fuCurrDpbIdx ][ uCrobIdx ][ uFebIdx ],
1379 usElinkIdx );
1380
1381 FillHitInfo( mess, usElinkIdx, uAsicIdx, uMsIdx );
1382 break;
1383 } // case stsxyter::MessType::Hit :
1384 case stsxyter::MessType::TsMsb :
1385 {
1386 FillTsMsbInfo( mess, uIdx, uMsIdx );
1387 break;
1388 } // case stsxyter::MessType::TsMsb :
1389 case stsxyter::MessType::Epoch :
1390 {
1391 // The first message in the TS is a special ones: EPOCH
1392 FillEpochInfo( mess );
1393
1394 if( 0 < uIdx )
1395 LOG(info) << "CbmMcbm2018MonitorSts::DoUnpack => "
1396 << "EPOCH message at unexpected position in MS: message "
1397 << uIdx << " VS message 0 expected!";
1398 break;
1399 } // case stsxyter::MessType::TsMsb :
1400 case stsxyter::MessType::Empty :
1401 {
1402// FillTsMsbInfo( mess );
1403 break;
1404 } // case stsxyter::MessType::Empty :
1405 case stsxyter::MessType::Dummy :
1406 {
1407 break;
1408 } // case stsxyter::MessType::Dummy / ReadDataAck / Ack :
1409 default:
1410 {
1411 LOG(fatal) << "CbmMcbm2018MonitorSts::DoUnpack => "
1412 << "Unknown message type, should never happen, stopping here!";
1413 }
1414 } // switch( mess.GetMessType() )
1415 } // for( uint32_t uIdx = 0; uIdx < uNbMessages; ++uIdx )
1416*/
1417 } // for( UInt_t uMsComp = 0; uMsComp < fvMsComponentsList.size(); ++uMsComp )
1418
1420 // Sort the buffer of hits
1421 std::sort(fvmHitsInMs.begin(), fvmHitsInMs.end());
1422
1423 // Time differences plotting using the fully time sorted hits
1424 if (0 < fvmHitsInMs.size()) {
1425 // ULong64_t ulLastHitTime = ( *( fvmHitsInMs.rbegin() ) ).GetTs();
1426 std::vector<stsxyter::FinalHit>::iterator itA;
1427 // comment unused variable, FU, 18.01.21 std::vector<stsxyter::FinalHit>::iterator itB;
1428
1429 // std::chrono::steady_clock::time_point tNow = std::chrono::steady_clock::now();
1430 // Double_t dUnixTimeInRun = std::chrono::duration_cast< std::chrono::seconds >(tNow - ftStartTimeUnix).count();
1431
1432 for (itA = fvmHitsInMs.begin(); itA != fvmHitsInMs.end();
1433 // itA != fvmHitsInMs.end() && (*itA).GetTs() < ulLastHitTime - 320; // 320 * 3.125 ns = 1000 ns
1434 ++itA) {
1435 UShort_t usAsicIdx = (*itA).GetAsic();
1436 // UShort_t usChanIdx = (*itA).GetChan();
1437 // ULong64_t ulHitTs = (*itA).GetTs();
1438 // UShort_t usHitAdc = (*itA).GetAdc();
1439 UShort_t usFebIdx = usAsicIdx / fUnpackParSts->GetNbAsicsPerFeb();
1440 // UShort_t usAsicInFeb = usAsicIdx % fUnpackParSts->GetNbAsicsPerFeb();
1441
1442 // Double_t dTimeSinceStartSec = (ulHitTs * stsxyter::kdClockCycleNs - fdStartTime)* 1e-9;
1443
1444 fvmAsicHitsInMs[usAsicIdx].push_back((*itA));
1445 fvmFebHitsInMs[usFebIdx].push_back((*itA));
1446 /*
1447 if( 1000 == fulCurrentTsIdx )
1448 {
1449 LOG(info) << Form( "FEB %02u ASIC %u Chan %03u TS %12u ADC %2u Time %8.3f",
1450 usFebIdx, usAsicInFeb, usChanIdx, ulHitTs, usHitAdc,
1451 ulHitTs* stsxyter::kdClockCycleNs );
1452 } // if( 0 == fulCurrentTsIdx )
1453*/
1454 } // loop on time sorted hits and split per asic/feb
1455
1456 // Remove all hits which were already used
1457 fvmHitsInMs.erase(fvmHitsInMs.begin(), itA);
1458
1460 Bool_t bHitCopyInThisMs[fuNbStsXyters];
1461 Bool_t bHitCopySameAdcInThisMs[fuNbStsXyters];
1462 Bool_t bFlagOnInThisMs[fuNbStsXyters];
1463 if (kTRUE == fbEnableCheckBugSmx20)
1464 for (UInt_t uAsic = 0; uAsic < fuNbStsXyters; uAsic++) {
1465 bHitCopyInThisMs[uAsic] = kFALSE;
1466 bHitCopySameAdcInThisMs[uAsic] = kFALSE;
1467 bFlagOnInThisMs[uAsic] = kFALSE;
1468 } // for( UInt_t uAsic = 0; uAsic < fuNbStsXyters; uAsic++)
1470
1471 for (UInt_t uAsic = 0; uAsic < fuNbStsXyters; uAsic++) {
1472 UInt_t uFebIdx = uAsic / fUnpackParSts->GetNbAsicsPerFeb();
1473 UInt_t uAsicInFeb = uAsic % fUnpackParSts->GetNbAsicsPerFeb();
1474
1475 std::vector<ULong64_t> vulLastHitTs(fUnpackParSts->GetNbChanPerAsic(), 0);
1476 std::vector<UShort_t> vusLastHitAdc(fUnpackParSts->GetNbChanPerAsic(), 0);
1477
1478 for (itA = fvmAsicHitsInMs[uAsic].begin(); itA != fvmAsicHitsInMs[uAsic].end(); ++itA) {
1479 // UShort_t usAsicIdx = (*itA).GetAsic();
1480 UShort_t usChanIdx = (*itA).GetChan();
1481 ULong64_t ulHitTs = (*itA).GetTs();
1482 UShort_t usHitAdc = (*itA).GetAdc();
1483
1484 UInt_t uChanInFeb = usChanIdx + fUnpackParSts->GetNbChanPerAsic() * uAsicInFeb;
1485
1487 if (kTRUE == fbEnableCheckBugSmx20) {
1489 Bool_t bIsNotCopy = kTRUE;
1490 if (vulLastHitTs[usChanIdx] == ulHitTs) {
1491 bIsNotCopy = kFALSE;
1492 bHitCopyInThisMs[uAsic] = kTRUE;
1493 if (vusLastHitAdc[usChanIdx] == usHitAdc) bHitCopySameAdcInThisMs[uAsic] = kTRUE;
1494 } // if( vulLastHitTs[ usChanIdx ] == ulHitTs)
1495
1496 vulLastHitTs[usChanIdx] = ulHitTs;
1497 vusLastHitAdc[usChanIdx] = usHitAdc;
1498
1499 if (bIsNotCopy) {
1500 fhStsFebChanCntRawGood[uFebIdx]->Fill(uChanInFeb);
1501 bFlagOnInThisMs[uAsic] |= SmxErrCheckCoinc(uFebIdx, uAsicInFeb, ulHitTs * stsxyter::kdClockCycleNs);
1502 } // if( bIsNotCopy )
1503 } // if( kTRUE == fbEnableCheckBugSmx20 )
1505 } // for( it = fvmAsicHitsInMs[ uAsic ].begin(); it != fvmAsicHitsInMs[ uAsic ].end(); ++it )
1506
1508 fvmAsicHitsInMs[uAsic].clear();
1509 } // for( UInt_t uAsic = 0; uAsic < fuNbStsXyters; uAsic++)
1510
1512 if (kTRUE == fbEnableCheckBugSmx20) {
1513 std::vector<Bool_t> vbCopyOnAnyAsicMs(fuNbFebs, kFALSE);
1514 std::vector<Bool_t> vbCopySameAdcOnAnyAsicMs(fuNbFebs, kFALSE);
1515 std::vector<Bool_t> vbFlagOnAnyAsicMs(fuNbFebs, kFALSE);
1516 // Bool_t bCopyOnAnyMs = kFALSE;
1517 // Bool_t bCopySameAdcOnAnyMs = kFALSE;
1518 // Bool_t bFlagOnAnyMs = kFALSE;
1519 for (UInt_t uAsic = 0; uAsic < fuNbStsXyters; uAsic++) {
1520 UInt_t uFebIdx = uAsic / fUnpackParSts->GetNbAsicsPerFeb();
1521 UInt_t uAsicInFeb = uAsic % fUnpackParSts->GetNbAsicsPerFeb();
1522
1523 fhStsFebSmxErrRatioCopyEvoAsic[uFebIdx]->Fill(dMsTime - fdStartTimeMsSz, uAsicInFeb,
1524 bHitCopyInThisMs[uAsic] ? 1.0 : 0.0);
1525
1526 fhStsFebSmxErrRatioCopySameAdcEvoAsic[uFebIdx]->Fill(dMsTime - fdStartTimeMsSz, uAsicInFeb,
1527 bHitCopySameAdcInThisMs[uAsic] ? 1.0 : 0.0);
1528
1529 fhStsFebSmxErrRatioEvoAsic[uFebIdx]->Fill(dMsTime - fdStartTimeMsSz, uAsicInFeb,
1530 bFlagOnInThisMs[uAsic] ? 1.0 : 0.0);
1531
1532 vbCopyOnAnyAsicMs[uFebIdx] = vbCopyOnAnyAsicMs[uFebIdx] || bHitCopyInThisMs[uAsic];
1533 vbCopySameAdcOnAnyAsicMs[uFebIdx] = vbCopySameAdcOnAnyAsicMs[uFebIdx] || bHitCopySameAdcInThisMs[uAsic];
1534 vbFlagOnAnyAsicMs[uFebIdx] = vbFlagOnAnyAsicMs[uFebIdx] || bFlagOnInThisMs[uAsic];
1535
1536 // bCopyOnAnyMs |= bHitCopyInThisMs[uAsic];
1537 // bCopySameAdcOnAnyMs |= bHitCopySameAdcInThisMs[uAsic];
1538 // bFlagOnAnyMs |= bFlagOnInThisMs[uAsic];
1539 } // for( UInt_t uAsic = 0; uAsic < fuNbStsXyters; uAsic++)
1540
1541 for (UInt_t uFebIdx = 0; uFebIdx < fuNbFebs; ++uFebIdx) {
1542 fhStsFebSmxErrRatioCopyEvo[uFebIdx]->Fill(dMsTime - fdStartTimeMsSz, vbCopyOnAnyAsicMs[uFebIdx] ? 1.0 : 0.0);
1543 fhStsFebSmxErrRatioCopySameAdcEvo[uFebIdx]->Fill(dMsTime - fdStartTimeMsSz,
1544 vbCopySameAdcOnAnyAsicMs[uFebIdx] ? 1.0 : 0.0);
1545 fhStsFebSmxErrRatioEvo[uFebIdx]->Fill(dMsTime - fdStartTimeMsSz, vbFlagOnAnyAsicMs[uFebIdx] ? 1.0 : 0.0);
1546 } // for( UInt_t uFebIdx = 0; uFebIdx < fuNbFebs; ++uFebIdx )
1547 } // if( kTRUE == fbEnableCheckBugSmx20 )
1549
1550 for (UInt_t uFebIdx = 0; uFebIdx < fuNbFebs; ++uFebIdx) {
1551 /*
1552 if( kTRUE == fbEnableCoincidenceMaps )
1553 {
1554 UInt_t uFebA = uFebIdx % fUnpackParSts->GetNbFebsPerCrob();
1555 UInt_t uCrobIdxA = ( uFebIdx / fUnpackParSts->GetNbFebsPerCrob() ) % fUnpackParSts->GetNbCrobsPerDpb();
1556 UInt_t uDpbIdxA = ( uFebIdx / fUnpackParSts->GetNbFebsPerCrob() ) / fUnpackParSts->GetNbCrobsPerDpb();
1557
1558 for( itA = fvmFebHitsInMs[ uFebIdx ].begin(); itA != fvmFebHitsInMs[ uFebIdx ].end(); ++itA )
1559 {
1560 UShort_t usAsicIdxA = (*itA).GetAsic();
1561 UShort_t usAsicInFebA = usAsicIdxA % fUnpackParSts->GetNbAsicsPerFeb();
1562 UShort_t usChanIdxA = (*itA).GetChan();
1563 UInt_t uChanInFebA = usChanIdxA + fUnpackParSts->GetNbChanPerAsic() * usAsicInFebA;
1564 ULong64_t ulHitTsA = (*itA).GetTs();
1565 Double_t dHitTsA = ulHitTsA * stsxyter::kdClockCycleNs;
1566
1568 if( -1 < fdStsFebChanLastTimeForDist[ uFebIdx ][ uChanInFebA ] )
1569 {
1570 fhStsFebChanDistT[ uFebIdx ]->Fill( dHitTsA - fdStsFebChanLastTimeForDist[ uFebIdx ][ uChanInFebA ],
1571 uChanInFebA );
1572 } // if( -1 < fdStsFebChanLastTimeForDist[ uFebIdx ][ uChanInFebA ] )
1573 fdStsFebChanLastTimeForDist[ uFebIdx ][ uChanInFebA ] = dHitTsA;
1574
1575 for( UInt_t uFebIdxB = uFebIdx; uFebIdxB < fuNbFebs; ++uFebIdxB )
1576 {
1578 Bool_t bSameModulePNsides = kFALSE;
1579 Double_t dBestDtMatch = 1e9; // ns
1580 std::vector< stsxyter::FinalHit >::iterator itBestMatch;
1581 UInt_t uFebB = uFebIdxB % fUnpackParSts->GetNbFebsPerCrob();
1582 UInt_t uCrobIdxB = ( uFebIdxB / fUnpackParSts->GetNbFebsPerCrob() ) % fUnpackParSts->GetNbCrobsPerDpb();
1583 UInt_t uDpbIdxB = ( uFebIdxB / fUnpackParSts->GetNbFebsPerCrob() ) / fUnpackParSts->GetNbCrobsPerDpb();
1584 if( fviFebModuleIdx[ uDpbIdxA ][ uCrobIdxA ][ uFebA ] == fviFebModuleIdx[ uDpbIdxB ][ uCrobIdxB ][ uFebB ] &&
1585 fviFebModuleSide[ uDpbIdxA ][ uCrobIdxA ][ uFebA ] != fviFebModuleSide[ uDpbIdxB ][ uCrobIdxB ][ uFebB ] )
1586 bSameModulePNsides = kTRUE;
1587
1588 for( itB = fvmFebHitsInMs[ uFebIdxB ].begin(); itB != fvmFebHitsInMs[ uFebIdxB ].end(); ++itB )
1589 {
1590 UShort_t usAsicIdxB = (*itB).GetAsic();
1591 UShort_t usChanIdxB = (*itB).GetChan();
1592 UInt_t uChanInFebB = usChanIdxB + fUnpackParSts->GetNbChanPerAsic() * (usAsicIdxB % fUnpackParSts->GetNbAsicsPerFeb());
1593
1594 if( uFebIdx == uFebIdxB && uChanInFebA == uChanInFebB )
1595 continue;
1596
1597 ULong64_t ulHitTsB = (*itB).GetTs();
1598 Double_t dHitTsB = ulHitTsB * stsxyter::kdClockCycleNs;
1599 Double_t dDtClk = static_cast< Double_t >( ulHitTsB ) - static_cast< Double_t >( ulHitTsA );
1600 Double_t dDt = dDtClk * stsxyter::kdClockCycleNs;
1601
1602 fhStsFebChanDtCoinc[ uFebIdx ][ uFebIdxB ]->Fill( dDt );
1603
1605 if( -1.0 * fdFebChanCoincidenceLimit < dDt )
1606 {
1608 if( fdFebChanCoincidenceLimit < dDt )
1609 break;
1610
1611 fhStsFebChanCoinc[ uFebIdx ][ uFebIdxB ]->Fill( uChanInFebA, uChanInFebB );
1612
1614 if( kTRUE == bSameModulePNsides )
1615 {
1617 if( TMath::Abs( dDt ) < TMath::Abs( dBestDtMatch ) )
1618 {
1619 itBestMatch = itB;
1620 dBestDtMatch = dDt;
1621 } // if( dDt < dBestDtMatch )
1622 } // if same module and opposite sides
1623 } // if( -1.0 * fdFebChanCoincidenceLimit < dDt )
1624 } // for( itB = fvmFebHitsInMs[ uFebIdxB ].begin(); itB != fvmFebHitsInMs[ uFebIdxB ].end(); ++itB )
1625
1627 if( kTRUE == bSameModulePNsides && dBestDtMatch < fdFebChanCoincidenceLimit )
1628 {
1629 UInt_t uModIdx = fviFebModuleIdx[ uDpbIdxA ][ uCrobIdxA ][ uFebA ];
1630
1631 UShort_t usAsicIdxB = (*itB).GetAsic();
1632 UShort_t usAsicInFebB = usAsicIdxB % fUnpackParSts->GetNbAsicsPerFeb();
1633 UShort_t usChanIdxB = (*itB).GetChan();
1634 UInt_t uChanInFebB = usChanIdxB + fUnpackParSts->GetNbChanPerAsic() * usAsicInFebB;
1635
1636 UShort_t usAdcA = (*itA).GetAdc();
1637 UShort_t usAdcB = (*itBestMatch).GetAdc();
1638 Double_t dPosX = 0.0;
1639 Double_t dPosY = 0.0;
1640
1641 fhStsModulePNCoincDt[ uModIdx ]->Fill( dBestDtMatch );
1642 if( 0 == fviFebModuleSide[ uDpbIdxA ][ uCrobIdxA ][ uFebA ] )
1643 {
1645 fUnpackParSts->ComputeModuleCoordinates( uModIdx, uChanInFebB, uChanInFebA, dPosX, dPosY );
1646
1648 fhStsModulePNCoincChan[ uModIdx ]->Fill( uChanInFebA, uChanInFebB );
1649 fhStsModulePNCoincAdc[ uModIdx ]->Fill( usAdcA, usAdcB );
1650 fhStsModuleCoincAdcChanP[ uModIdx ]->Fill( uChanInFebA, usAdcA );
1651 fhStsModuleCoincAdcChanN[ uModIdx ]->Fill( uChanInFebB, usAdcB );
1652 fhStsModulePNCoincDtAsicP[ uModIdx ]->Fill( dBestDtMatch, usAsicInFebA );
1653 fhStsModulePNCoincDtAsicN[ uModIdx ]->Fill( dBestDtMatch, usAsicInFebB );
1654 fhStsModuleCoincMap[ uModIdx ]->Fill( dPosX, dPosY );
1655 } // if( 0 == fviFebModuleSide[ uDpbIdxA ][ uCrobIdxA ][ uFebA ] )
1656 else
1657 {
1659 fUnpackParSts->ComputeModuleCoordinates( uModIdx, uChanInFebA, uChanInFebB, dPosX, dPosY );
1660
1662 fhStsModulePNCoincChan[ uModIdx ]->Fill( uChanInFebB, uChanInFebA );
1663 fhStsModulePNCoincAdc[ uModIdx ]->Fill( usAdcB, usAdcA );
1664 fhStsModuleCoincAdcChanP[ uModIdx ]->Fill( uChanInFebB, usAdcB );
1665 fhStsModuleCoincAdcChanN[ uModIdx ]->Fill( uChanInFebA, usAdcA );
1666 fhStsModulePNCoincDtAsicP[ uModIdx ]->Fill( dBestDtMatch, usAsicInFebB );
1667 fhStsModulePNCoincDtAsicN[ uModIdx ]->Fill( dBestDtMatch, usAsicInFebA );
1668 fhStsModuleCoincMap[ uModIdx ]->Fill( dPosX, dPosY );
1669 } // else of if( 0 == fviFebModuleSide[ uDpbIdxA ][ uCrobIdxA ][ uFebA ] )
1670
1672 fhStsModuleCoincMap[ uModIdx ]->Fill( dPosX, dPosY );
1673 } // if same module and opposite sides
1674 } // for( UInt_t uFebIdxB = uFebIdx; uFebIdxB < fuNbFebs; ++uFebIdxB )
1675 } // for( itA = fvmFebHitsInMs[ uFebIdx ].begin(); itA != fvmFebHitsInMs[ uFebIdx ].end(); ++itA )
1676 } // if( kTRUE == fbEnableCoincidenceMaps )
1677*/
1679 fvmFebHitsInMs[uFebIdx].clear();
1680 } // for( UInt_t uFebIdx = 0; uFebIdx < fuNbFebs; ++uFebIdx )
1681 } // if( 0 < fvmHitsInMs.size() )
1682 } // for( UInt_t uMsIdx = 0; uMsIdx < uNbMsLoop; uMsIdx ++ )
1683
1684 for (UInt_t uMsIdx = 0; uMsIdx < fuMaxNbMicroslices; ++uMsIdx) {
1685 fvdMsTime[uMsIdx] = 0.0;
1686 } // for( UInt_t uMsIdx = 0; uMsIdx < fuMaxNbMicroslices; ++uMsIdx )
1687
1688 if (0 == ts.index() % 1000) {
1689 for (UInt_t uDpb = 0; uDpb < fuNrOfDpbs; ++uDpb) {
1690 Double_t dTsMsbTime =
1691 static_cast<ULong64_t>(stsxyter::kuHitNbTsBins) * static_cast<ULong64_t>(fvulCurrentTsMsb[fuCurrDpbIdx])
1692 + static_cast<ULong64_t>(stsxyter::kulTsCycleNbBins)
1693 * static_cast<ULong64_t>(fvuCurrentTsMsbCycle[fuCurrDpbIdx]);
1694 dTsMsbTime *= stsxyter::kdClockCycleNs * 1e-9;
1695
1696 LOG(info) << "End of TS " << std::setw(7) << ts.index() << " eDPB " << std::setw(2) << uDpb
1697 << " current TS MSB counter is " << std::setw(12) << fvulCurrentTsMsb[uDpb]
1698 << " current TS MSB cycle counter is " << std::setw(12) << fvuCurrentTsMsbCycle[uDpb]
1699 << " current TS MSB time is " << std::setw(12) << dTsMsbTime << " s";
1700 }
1701 } // if( 0 == ts.index() % 1000 )
1702
1703 if (0 == ts.index() % 10000) SaveAllHistos("data/PulserPeriodicHistosSave.root");
1704
1705 return kTRUE;
1706}
1707
1708Bool_t CbmMcbm2018MonitorSts::ProcessStsMs(const fles::Timeslice& ts, size_t uMsComp, UInt_t uMsIdx)
1709{
1710 auto msDescriptor = ts.descriptor(uMsComp, uMsIdx);
1711 fuCurrentEquipmentId = msDescriptor.eq_id;
1712 const uint8_t* msContent = reinterpret_cast<const uint8_t*>(ts.content(uMsComp, uMsIdx));
1713
1714 fulCurrentTsIdx = ts.index();
1715 if (0 == fvbMaskedComponents.size()) fvbMaskedComponents.resize(ts.num_components(), kFALSE);
1716
1717 if (0 == fulCurrentTsIdx && 0 == uMsIdx) {
1718 /*
1719 LOG(info) << "hi hv eqid flag si sv idx/start crc size offset";
1720 LOG(info) << Form( "%02x %02x %04x %04x %02x %02x %016llx %08x %08x %016llx",
1721 static_cast<unsigned int>(msDescriptor.hdr_id),
1722 static_cast<unsigned int>(msDescriptor.hdr_ver), msDescriptor.eq_id, msDescriptor.flags,
1723 static_cast<unsigned int>(msDescriptor.sys_id),
1724 static_cast<unsigned int>(msDescriptor.sys_ver), msDescriptor.idx, msDescriptor.crc,
1725 msDescriptor.size, msDescriptor.offset );
1726*/
1727 LOG(info) << FormatMsHeaderPrintout(msDescriptor);
1728 } // if( 0 == fulCurrentTsIndex && 0 == uMsIdx )
1729 if (kFALSE == fvbMaskedComponents[uMsComp] && 0 == uMsIdx) {
1730 auto it = fDpbIdIndexMap.find(fuCurrentEquipmentId);
1731 if (fDpbIdIndexMap.end() == it) {
1732 LOG(warning) << "Could not find the sDPB index for AFCK id 0x" << std::hex << fuCurrentEquipmentId << std::dec
1733 << " component " << uMsComp << "\n"
1734 << "If valid this index has to be added in the TOF parameter file in "
1735 "the RocIdArray field"
1736 << "\n"
1737 << "For now we remove it from the list of components analyzed";
1738 fvbMaskedComponents[uMsComp] = kTRUE;
1739 } // if( fDpbIdIndexMap.end() == it )
1740
1741 } // if( kFALSE == fvbMaskedComponents[ uMsComp ] && 0 == uMsIdx )
1742
1743 if (kTRUE == fvbMaskedComponents[uMsComp]) return kTRUE;
1744
1745 uint32_t uSize = msDescriptor.size;
1746 fulCurrentMsIdx = msDescriptor.idx;
1747 Double_t dMsTime = (1e-9) * static_cast<double>(fulCurrentMsIdx);
1748 LOG(debug) << "Microslice: " << fulCurrentMsIdx << " from EqId " << std::hex << fuCurrentEquipmentId << std::dec
1749 << " has size: " << uSize;
1750
1751 fuCurrDpbId = static_cast<uint32_t>(fuCurrentEquipmentId & 0xFFFF);
1753
1754 if (uMsComp < kiMaxNbFlibLinks) {
1755 if (fdStartTimeMsSz < 0) fdStartTimeMsSz = dMsTime;
1756 fhMsSz[uMsComp]->Fill(uSize);
1757 fhMsSzTime[uMsComp]->Fill(dMsTime - fdStartTimeMsSz, uSize);
1758 } // if( uMsComp < kiMaxNbFlibLinks )
1759
1761 if (static_cast<Int_t>(fvdPrevMsTime[uMsComp]) < static_cast<Int_t>(dMsTime)) {
1763 UInt_t uFebIdxOffset = fUnpackParSts->GetNbFebsPerDpb() * fuCurrDpbIdx;
1764 for (UInt_t uFebIdx = 0; uFebIdx < fUnpackParSts->GetNbFebsPerDpb(); ++uFebIdx) {
1765 UInt_t uFebIdxInSyst = uFebIdxOffset + uFebIdx;
1766
1768 if (fvdFebTimeSecLastRateUpdate[uFebIdxInSyst] < 0.0) {
1769 fvdFebTimeSecLastRateUpdate[uFebIdxInSyst] = dMsTime;
1770 fviFebCountsSinceLastRateUpdate[uFebIdxInSyst] = 0;
1771 for (UInt_t uChan = 0; uChan < fUnpackParSts->GetNbChanPerFeb(); ++uChan)
1772 fvdFebChanCountsSinceLastRateUpdate[uFebIdxInSyst][uChan] = 0.0;
1773 continue;
1774 } // if( fvdFebTimeSecLastRateUpdate[uFebIdxInSyst] < 0.0 )
1775
1776 Double_t dTimeInt = dMsTime - fvdFebTimeSecLastRateUpdate[uFebIdxInSyst];
1777 if (fiTimeIntervalRateUpdate <= dTimeInt) {
1779 if (0 == fviFebCountsSinceLastRateUpdate[uFebIdxInSyst]) {
1780 fvdFebTimeSecLastRateUpdate[uFebIdxInSyst] = dMsTime;
1781 continue;
1782 } // if( 0 == fviFebCountsSinceLastRateUpdate[uFebIdxInSyst] )
1783
1784 for (UInt_t uChan = 0; uChan < fUnpackParSts->GetNbChanPerFeb(); ++uChan) {
1785 fhStsFebChanHitRateProf[uFebIdxInSyst]->Fill(uChan, fvdFebChanCountsSinceLastRateUpdate[uFebIdxInSyst][uChan]
1786 / dTimeInt);
1787 fvdFebChanCountsSinceLastRateUpdate[uFebIdxInSyst][uChan] = 0.0;
1788 } // for( UInt_t uChan = 0; uChan < fUnpackParSts->GetNbChanPerFeb(); ++uChan )
1789
1790 fvdFebTimeSecLastRateUpdate[uFebIdxInSyst] = dMsTime;
1791 fviFebCountsSinceLastRateUpdate[uFebIdxInSyst] = 0;
1792 } // if( fiTimeIntervalRateUpdate <= dTimeInt )
1793 } // for( UInt_t uFebIdx = 0; uFebIdx < fUnpackParSts->GetNbFebsPerDpb(); ++uFebIdx )
1794 } // if( static_cast<Int_t>( fvdMsTime[ uMsCompIdx ] ) < static_cast<Int_t>( dMsTime ) )
1795
1796 // Store MS time for coincidence plots
1797 fvdPrevMsTime[uMsComp] = dMsTime;
1798 /*
1800// if( 0 == fuCurrDpbIdx && fulCurrentTsIdx < 100 )
1801 if( fulCurrentTsIdx < 100 )
1802 LOG(info) << "TS " << std::setw( 12 ) << fulCurrentTsIdx
1803 << " MS " << std::setw( 12 ) << fulCurrentMsIdx
1804 << " MsInTs " << std::setw( 3 ) << uMsIdx
1805 << " DPB " << fuCurrDpbIdx;
1806*/
1807
1809 uint16_t uMsHeaderFlags = msDescriptor.flags;
1810 for (UInt_t uBit = 0; uBit < 16; ++uBit)
1811 fhStsMsStatusFieldType->Fill(uBit, (uMsHeaderFlags >> uBit) & 0x1);
1812
1814 UInt_t uTsMsbCycleHeader = std::floor(fulCurrentMsIdx / (stsxyter::kulTsCycleNbBins * stsxyter::kdClockCycleNs))
1816
1818 if (kTRUE == fbBinningFw)
1821
1822 if (kFALSE == fvuInitialHeaderDone[fuCurrDpbIdx]) {
1823 fvuInitialTsMsbCycleHeader[fuCurrDpbIdx] = uTsMsbCycleHeader;
1826 } // if( kFALSE == fvuInitialHeaderDone[ fuCurrDpbIdx ] )
1827 else if (0 == uMsIdx) {
1828 if (uTsMsbCycleHeader != fvuCurrentTsMsbCycle[fuCurrDpbIdx])
1829 LOG(info) << " TS " << std::setw(12) << fulCurrentTsIdx << " MS " << std::setw(12) << fulCurrentMsIdx
1830 << " MS Idx " << std::setw(4) << uMsIdx << " Msg Idx " << std::setw(5) << 0 << " DPB " << std::setw(2)
1831 << fuCurrDpbIdx << " Old TsMsb " << std::setw(5) << fvulCurrentTsMsb[fuCurrDpbIdx] << " Old MsbCy "
1832 << std::setw(5) << fvuCurrentTsMsbCycle[fuCurrDpbIdx] << " New MsbCy " << uTsMsbCycleHeader;
1833 fvuCurrentTsMsbCycle[fuCurrDpbIdx] = uTsMsbCycleHeader;
1835 } // if( 0 == uMsIdx )
1836 else if (uTsMsbCycleHeader != fvuCurrentTsMsbCycle[fuCurrDpbIdx] && 4194303 != fvulCurrentTsMsb[fuCurrDpbIdx]) {
1837 LOG(warning) << "TS MSB cycle from MS header does not match current cycle from data "
1838 << "for TS " << std::setw(12) << fulCurrentTsIdx << " MS " << std::setw(12) << fulCurrentMsIdx
1839 << " MsInTs " << std::setw(3) << uMsIdx << " ====> " << fvuCurrentTsMsbCycle[fuCurrDpbIdx] << " VS "
1840 << uTsMsbCycleHeader;
1841 fvuCurrentTsMsbCycle[fuCurrDpbIdx] = uTsMsbCycleHeader;
1842 }
1843 else if (uTsMsbCycleHeader != fvuCurrentTsMsbCycle[fuCurrDpbIdx]) {
1845 if (4194303 == fvulCurrentTsMsb[fuCurrDpbIdx]) {
1846 LOG(info) << " TS " << std::setw(12) << fulCurrentTsIdx << " MS " << std::setw(12) << fulCurrentMsIdx
1847 << " MS Idx " << std::setw(4) << uMsIdx << " Msg Idx " << std::setw(5) << 0 << " DPB " << std::setw(2)
1848 << fuCurrDpbIdx << " Old TsMsb " << std::setw(5) << fvulCurrentTsMsb[fuCurrDpbIdx] << " Old MsbCy "
1849 << std::setw(5) << fvuCurrentTsMsbCycle[fuCurrDpbIdx] << " New MsbCy " << uTsMsbCycleHeader;
1850 }
1851 else {
1852 LOG(warning) << "TS MSB cycle from MS header does not match current cycle from data "
1853 << "for TS " << std::setw(12) << fulCurrentTsIdx << " MS " << std::setw(12) << fulCurrentMsIdx
1854 << " MsInTs " << std::setw(3) << uMsIdx << " ====> " << fvuCurrentTsMsbCycle[fuCurrDpbIdx]
1855 << " (cnt) VS " << uTsMsbCycleHeader << " (header)";
1856 } // else of if( 4194303 == fvulCurrentTsMsb[fuCurrDpbIdx] )
1857 fvuCurrentTsMsbCycle[fuCurrDpbIdx] = uTsMsbCycleHeader;
1858 } // else if( uTsMsbCycleHeader != fvuCurrentTsMsbCycle[ fuCurrDpbIdx ] )
1859
1860 // If not integer number of message in input buffer, print warning/error
1861 if (0 != (uSize % kuBytesPerMessage))
1862 LOG(error) << "The input microslice buffer does NOT "
1863 << "contain only complete nDPB messages!";
1864
1865 // Compute the number of complete messages in the input microslice buffer
1866 uint32_t uNbMessages = (uSize - (uSize % kuBytesPerMessage)) / kuBytesPerMessage;
1867
1868 // Prepare variables for the loop on contents
1869 const uint32_t* pInBuff = reinterpret_cast<const uint32_t*>(msContent);
1870
1871 for (uint32_t uIdx = 0; uIdx < uNbMessages; ++uIdx) {
1872 // Fill message
1873 uint32_t ulData = static_cast<uint32_t>(pInBuff[uIdx]);
1874
1875 stsxyter::Message mess(static_cast<uint32_t>(ulData & 0xFFFFFFFF));
1876
1877 // Print message if requested
1878 if (fbPrintMessages) mess.PrintMess(std::cout, fPrintMessCtrl);
1879
1880 /*
1881 if( 1000 == fulCurrentTsIdx )
1882 {
1883 mess.PrintMess( std::cout, fPrintMessCtrl );
1884 } // if( 0 == fulCurrentTsIdx )
1885*/
1886 stsxyter::MessType typeMess = mess.GetMessType();
1887 fmMsgCounter[typeMess]++;
1888 fhStsMessType->Fill(static_cast<uint16_t>(typeMess));
1889 fhStsMessTypePerDpb->Fill(fuCurrDpbIdx, static_cast<uint16_t>(typeMess));
1890
1891 switch (typeMess) {
1893 // Extract the eLink and Asic indices => Should GO IN the fill method now that obly hits are link/asic specific!
1894 UShort_t usElinkIdx = mess.GetLinkIndex();
1896 if (kTRUE == fbBinningFw) usElinkIdx = mess.GetLinkIndexHitBinning();
1897
1898 fhStsMessTypePerElink->Fill(fuCurrDpbIdx * fUnpackParSts->GetNbElinkPerDpb() + usElinkIdx,
1899 static_cast<uint16_t>(typeMess));
1900 fhStsHitsElinkPerDpb->Fill(fuCurrDpbIdx, usElinkIdx);
1901 /*
1903 if( ( 0 == fuCurrDpbIdx || 13 < usElinkIdx ) && fulCurrentTsIdx < 100 )
1904 mess.PrintMess( std::cout, stsxyter::MessagePrintMask::msg_print_Hex | stsxyter::MessagePrintMask::msg_print_Human );
1906*/
1907 UInt_t uCrobIdx = usElinkIdx / fUnpackParSts->GetNbElinkPerCrob();
1908 Int_t uFebIdx = fUnpackParSts->ElinkIdxToFebIdx(usElinkIdx);
1909
1910 if (-1 == uFebIdx) {
1911 LOG(warning) << "CbmMcbm2018MonitorSts::DoUnpack => "
1912 << "Wrong elink Idx! Elink raw " << Form("%d remap %d", usElinkIdx, uFebIdx);
1913 continue;
1914 } // if( -1 == uFebIdx )
1915
1916 UInt_t uAsicIdx =
1917 (fuCurrDpbIdx * fUnpackParSts->GetNbCrobsPerDpb() + uCrobIdx) * fUnpackParSts->GetNbAsicsPerCrob()
1918 + fUnpackParSts->ElinkIdxToAsicIdx(1 == fviFebType[fuCurrDpbIdx][uCrobIdx][uFebIdx], usElinkIdx);
1919
1920 FillHitInfo(mess, usElinkIdx, uAsicIdx, uMsIdx);
1921 break;
1922 } // case stsxyter::MessType::Hit :
1924 fhStsMessTypePerElink->Fill(fuCurrDpbIdx * fUnpackParSts->GetNbElinkPerDpb(), static_cast<uint16_t>(typeMess));
1925
1926 FillTsMsbInfo(mess, uIdx, uMsIdx);
1927 break;
1928 } // case stsxyter::MessType::TsMsb :
1930 fhStsMessTypePerElink->Fill(fuCurrDpbIdx * fUnpackParSts->GetNbElinkPerDpb(), static_cast<uint16_t>(typeMess));
1931
1932 // The first message in the TS is a special ones: EPOCH
1933 FillEpochInfo(mess);
1934
1935 if (0 < uIdx)
1936 LOG(info) << "CbmMcbm2018MonitorSts::DoUnpack => "
1937 << "EPOCH message at unexpected position in MS: message " << uIdx << " VS message 0 expected!";
1938 break;
1939 } // case stsxyter::MessType::TsMsb :
1941 UShort_t usElinkIdx = mess.GetStatusLink();
1942 fhStsMessTypePerElink->Fill(fuCurrDpbIdx * fUnpackParSts->GetNbElinkPerDpb() + usElinkIdx,
1943 static_cast<uint16_t>(typeMess));
1944 /*
1946 if( ( 0 == fuCurrDpbIdx || 13 < usElinkIdx ) && fulCurrentTsIdx < 100 )
1947 mess.PrintMess( std::cout, stsxyter::MessagePrintMask::msg_print_Hex | stsxyter::MessagePrintMask::msg_print_Human );
1949*/
1950 UInt_t uCrobIdx = usElinkIdx / fUnpackParSts->GetNbElinkPerCrob();
1951 Int_t uFebIdx = fUnpackParSts->ElinkIdxToFebIdx(usElinkIdx);
1952 UInt_t uAsicIdx =
1953 (fuCurrDpbIdx * fUnpackParSts->GetNbCrobsPerDpb() + uCrobIdx) * fUnpackParSts->GetNbAsicsPerCrob()
1954 + fUnpackParSts->ElinkIdxToAsicIdx(1 == fviFebType[fuCurrDpbIdx][uCrobIdx][uFebIdx], usElinkIdx);
1955
1956 UShort_t usStatusField = mess.GetStatusStatus();
1957
1958 fhStsStatusMessType->Fill(uAsicIdx, usStatusField);
1960 if (fbPrintMessages) {
1961 std::cout << Form("DPB %2u TS %12llu MS %12llu mess %5u ", fuCurrDpbIdx, fulCurrentTsIdx, fulCurrentMsIdx,
1962 uIdx);
1963 mess.PrintMess(std::cout, fPrintMessCtrl);
1964 } // if( fbPrintMessages )
1965 // FillTsMsbInfo( mess );
1966 break;
1967 } // case stsxyter::MessType::Status
1969 fhStsMessTypePerElink->Fill(fuCurrDpbIdx * fUnpackParSts->GetNbElinkPerDpb(), static_cast<uint16_t>(typeMess));
1970 // FillTsMsbInfo( mess );
1971 break;
1972 } // case stsxyter::MessType::Empty :
1974 // ++vNbMessType[5];
1975 // sMessPatt += " En";
1976 // bError = pMess[uIdx].IsMsErrorFlagOn();
1977
1978 // fhStsMessTypePerElink->Fill( fuCurrDpbIdx * fUnpackPar->GetNbElinkPerDpb(), static_cast< uint16_t > (typeMess) );
1979 // FillTsMsbInfo( pMess[uIdx] );
1980 if (mess.IsMsErrorFlagOn()) {
1981 fhMsErrorsEvo->Fill(1e-9 * fulCurrentMsIdx, mess.GetMsErrorType());
1982 } // if( pMess[uIdx].IsMsErrorFlagOn() )
1983 break;
1984 } // case stsxyter::MessType::EndOfMs :
1986 fhStsMessTypePerElink->Fill(fuCurrDpbIdx * fUnpackParSts->GetNbElinkPerDpb(), static_cast<uint16_t>(typeMess));
1987 break;
1988 } // case stsxyter::MessType::Dummy / ReadDataAck / Ack :
1989 default: {
1990 LOG(fatal) << "CbmMcbm2018MonitorSts::DoUnpack => "
1991 << "Unknown message type, should never happen, stopping "
1992 "here! Type found was: "
1993 << static_cast<int>(typeMess);
1994 }
1995 } // switch( mess.GetMessType() )
1996 } // for( uint32_t uIdx = 0; uIdx < uNbMessages; ++uIdx )
1997
1998 return kTRUE;
1999}
2000
2001void CbmMcbm2018MonitorSts::FillHitInfo(stsxyter::Message mess, const UShort_t& /*usElinkIdx*/, const UInt_t& uAsicIdx,
2002 const UInt_t& uMsIdx)
2003{
2004 UShort_t usChan = mess.GetHitChannel();
2005 UShort_t usRawAdc = mess.GetHitAdc();
2006 // UShort_t usFullTs = mess.GetHitTimeFull();
2007 // UShort_t usTsOver = mess.GetHitTimeOver();
2008 UShort_t usRawTs = mess.GetHitTime();
2009
2011 if (kTRUE == fbBinningFw) usRawTs = mess.GetHitTimeBinning();
2012
2014 // usChan = 127 - usChan;
2015
2016 /*
2017 fhStsChanCntRaw[ uAsicIdx ]->Fill( usChan );
2018 fhStsChanAdcRaw[ uAsicIdx ]->Fill( usChan, usRawAdc );
2019 fhStsChanAdcRawProf[ uAsicIdx ]->Fill( usChan, usRawAdc );
2020 fhStsChanRawTs[ uAsicIdx ]->Fill( usChan, usRawTs );
2021 fhStsChanMissEvt[ uAsicIdx ]->Fill( usChan, mess.IsHitMissedEvts() );
2022*/
2023 // UInt_t uCrobIdx = usElinkIdx / fUnpackParSts->GetNbElinkPerCrob();
2024 UInt_t uFebIdx = uAsicIdx / fUnpackParSts->GetNbAsicsPerFeb();
2025 UInt_t uAsicInFeb = uAsicIdx % fUnpackParSts->GetNbAsicsPerFeb();
2026 UInt_t uChanInFeb = usChan + fUnpackParSts->GetNbChanPerAsic() * (uAsicIdx % fUnpackParSts->GetNbAsicsPerFeb());
2027
2028 // Double_t dCalAdc = fvdFebAdcOffs[ fuCurrDpbIdx ][ uCrobIdx ][ uFebIdx ]
2029 // + (usRawAdc - 1)* fvdFebAdcGain[ fuCurrDpbIdx ][ uCrobIdx ][ uFebIdx ];
2030
2031 fhStsFebChanCntRaw[uFebIdx]->Fill(uChanInFeb);
2032 fhStsFebChanAdcRaw[uFebIdx]->Fill(uChanInFeb, usRawAdc);
2033 fhStsFebChanAdcRawProf[uFebIdx]->Fill(uChanInFeb, usRawAdc);
2034 // fhStsFebChanAdcCal[ uFebIdx ]->Fill( uChanInFeb, dCalAdc );
2035 // fhStsFebChanAdcCalProf[ uFebIdx ]->Fill( uChanInFeb, dCalAdc );
2036 fhStsFebChanRawTs[uFebIdx]->Fill(uChanInFeb, usRawTs);
2037 fhStsFebChanMissEvt[uFebIdx]->Fill(uChanInFeb, mess.IsHitMissedEvts());
2038
2039 // Compute the Full time stamp
2040 // Long64_t ulOldHitTime = fvulChanLastHitTime[ uAsicIdx ][ usChan ];
2041 // Double_t dOldHitTime = fvdChanLastHitTime[ uAsicIdx ][ usChan ];
2042
2043 // Use TS w/o overlap bits as they will anyway come from the TS_MSB
2044 fvulChanLastHitTime[uAsicIdx][usChan] = usRawTs;
2045
2046 fvulChanLastHitTime[uAsicIdx][usChan] +=
2047 static_cast<ULong64_t>(stsxyter::kuHitNbTsBins) * static_cast<ULong64_t>(fvulCurrentTsMsb[fuCurrDpbIdx])
2048 + static_cast<ULong64_t>(stsxyter::kulTsCycleNbBins) * static_cast<ULong64_t>(fvuCurrentTsMsbCycle[fuCurrDpbIdx]);
2049
2051 if (kTRUE == fbBinningFw)
2052 fvulChanLastHitTime[uAsicIdx][usChan] =
2053 usRawTs
2054 + static_cast<ULong64_t>(stsxyter::kuHitNbTsBinsBinning) * static_cast<ULong64_t>(fvulCurrentTsMsb[fuCurrDpbIdx])
2055 + static_cast<ULong64_t>(stsxyter::kulTsCycleNbBinsBinning)
2056 * static_cast<ULong64_t>(fvuCurrentTsMsbCycle[fuCurrDpbIdx]);
2057
2058 // fvuElinkLastTsHit[fuCurrDpbIdx] = usRawTs;
2059
2060 // Convert the Hit time in bins to Hit time in ns
2061 Double_t dHitTimeNs = fvulChanLastHitTime[uAsicIdx][usChan] * stsxyter::kdClockCycleNs;
2062 /*
2063 // If needed fill the hit interval plots
2064 if( fbChanHitDtEna )
2065 {
2066 Double_t dDeltaT = dHitTimeNs - fvdChanLastHitTime[ uAsicIdx ][ usChan ];
2067 if( 0 == dDeltaT )
2068 fhStsChanHitDtNeg[ uAsicIdx ]->Fill( 1, usChan );
2069 else if( 0 < dDeltaT )
2070 fhStsChanHitDt[ uAsicIdx ]->Fill( dDeltaT, usChan );
2071 else fhStsChanHitDtNeg[ uAsicIdx ]->Fill( -dDeltaT, usChan );
2072 } // if( fbChanHitDtEna )
2073*/
2075 Double_t dDeltaT = dHitTimeNs - fvdChanLastHitTime[uAsicIdx][usChan];
2076 fhStsFebChanCloseHitsCounts[uFebIdx]->Fill(uChanInFeb, 0);
2077 fhStsFebChanCloseHitsRatio[uFebIdx]->Fill(uChanInFeb, 0);
2078 if (dDeltaT < 80.0) {
2079 fhStsFebChanCloseHitsCounts[uFebIdx]->Fill(uChanInFeb, 1);
2080 fhStsFebChanCloseHitsRatio[uFebIdx]->Fill(uChanInFeb, 1);
2081 // LOG(info) << Form( "Potential duplicate (bad Vref2) in FEB %2u ASIC %u channel %4u (dt = %f) TS %7llu MS %3u",
2082 // uFebIdx, uAsicInFeb, uChanInFeb, dDeltaT, fulCurrentTsIdx, uMsIdx );
2083 } //
2084
2085 // Store new value of Hit time in ns
2086 fvdChanLastHitTime[uAsicIdx][usChan] = fvulChanLastHitTime[uAsicIdx][usChan] * stsxyter::kdClockCycleNs;
2087 /*
2088 LOG(info) << " Asic " << std::setw( 2 ) << uAsicIdx
2089 << " Channel " << std::setw( 3 ) << usChan
2090 << " Diff to last hit " << std::setw( 12 ) << ( fvulChanLastHitTime[ uAsicIdx ][ usChan ] - ulOldHitTime)
2091 << " in s " << std::setw( 12 ) << ( fvdChanLastHitTime[ uAsicIdx ][ usChan ] - dOldHitTime) * 1e-9;
2092*/
2093 // Pulser and MS
2094 fvuChanNbHitsInMs[uAsicIdx][usChan][uMsIdx]++;
2095 fvdChanLastHitTimeInMs[uAsicIdx][usChan][uMsIdx] = dHitTimeNs;
2096 fvusChanLastHitAdcInMs[uAsicIdx][usChan][uMsIdx] = usRawAdc;
2097 /*
2098 fvmChanHitsInTs[ uAsicIdx ][ usChan ].insert( stsxyter::FinalHit( fvulChanLastHitTime[ uAsicIdx ][ usChan ],
2099 usRawAdc, uAsicIdx, usChan ) );
2100*/
2101 fvmHitsInMs.push_back(stsxyter::FinalHit(fvulChanLastHitTime[uAsicIdx][usChan], usRawAdc, uAsicIdx, usChan));
2102
2103 /*
2104 if( ( 214514 < fulCurrentTsIdx && fulCurrentTsIdx < 214517 ) ||
2105 ( 260113 < fulCurrentTsIdx && fulCurrentTsIdx < 260116 ) ||
2106 ( 388109 < fulCurrentTsIdx && fulCurrentTsIdx < 388114 ) ||
2107 ( 573361 < fulCurrentTsIdx && fulCurrentTsIdx < 573364 ) ||
2108 ( 661731 < fulCurrentTsIdx && fulCurrentTsIdx < 661734 ) ||
2109 ( 712982 < fulCurrentTsIdx && fulCurrentTsIdx < 712985 ) ||
2110 ( 713857 < fulCurrentTsIdx && fulCurrentTsIdx < 713860 ) ||
2111 ( 739365 < fulCurrentTsIdx && fulCurrentTsIdx < 739368 ))
2112 if( 0 < fulCurrentTsIdx && fulCurrentTsIdx < 3 )
2113 {
2114 LOG(info) << " TS " << std::setw( 12 ) << fulCurrentTsIdx
2115 << " MS " << std::setw( 12 ) << fulCurrentMsIdx
2116 << " MsInTs " << std::setw( 3 ) << uMsIdx
2117 << " Asic " << std::setw( 2 ) << uAsicIdx
2118 << " Channel " << std::setw( 3 ) << usChan
2119 << " ADC " << std::setw( 3 ) << usRawAdc
2120 << " TS " << std::setw( 3 ) << usRawTs // 9 bits TS
2121 << " SX TsMsb " << std::setw( 2 ) << ( fvulCurrentTsMsb[fuCurrDpbIdx] & 0x1F ) // Total StsXyter TS = 14 bits => 9b Hit TS + lower 5b TS_MSB after DPB
2122 << " DPB TsMsb " << std::setw( 6 ) << ( fvulCurrentTsMsb[fuCurrDpbIdx] >> 5 ) // Total StsXyter TS = 14 bits => 9b Hit TS + lower 5b of TS_MSB after DPB
2123 << " TsMsb " << std::setw( 7 ) << fvulCurrentTsMsb[fuCurrDpbIdx]
2124 << " TS(b) " << std::bitset<9>(usRawTs)
2125 << " TSM(b) " << std::bitset<24>(fvulCurrentTsMsb[fuCurrDpbIdx])
2126 << " MsbCy " << std::setw( 4 ) << fvuCurrentTsMsbCycle[fuCurrDpbIdx]
2127 << " Time " << std::setw ( 12 ) << fvulChanLastHitTime[ uAsicIdx ][ usChan ];
2128 }
2129*/
2130 // Check Starting point of histos with time as X axis
2131 if (-1 == fdStartTime) fdStartTime = fvdChanLastHitTime[uAsicIdx][usChan];
2132
2133 // Fill histos with time as X axis
2134 Double_t dTimeSinceStartSec = (fvdChanLastHitTime[uAsicIdx][usChan] - fdStartTime) * 1e-9;
2135 Double_t dTimeSinceStartMin = dTimeSinceStartSec / 60.0;
2136
2138 fvdFebChanCountsSinceLastRateUpdate[uFebIdx][uChanInFeb] += 1;
2139
2140 fhStsFebAsicHitCounts->Fill(uFebIdx, uAsicInFeb);
2141 fhStsAllFebsHitRateEvo->Fill(dTimeSinceStartSec, uFebIdx);
2142 fhStsAllAsicsHitRateEvo->Fill(dTimeSinceStartSec, uAsicIdx);
2143 fhStsFebChanHitRateEvo[uFebIdx]->Fill(dTimeSinceStartSec, uChanInFeb);
2144 fhStsFebAsicHitRateEvo[uFebIdx]->Fill(dTimeSinceStartSec, uAsicInFeb);
2145 fhStsFebHitRateEvo[uFebIdx]->Fill(dTimeSinceStartSec);
2146 fhStsFebChanHitRateEvoLong[uFebIdx]->Fill(dTimeSinceStartMin, uChanInFeb, 1.0 / 60.0);
2147 fhStsFebAsicHitRateEvoLong[uFebIdx]->Fill(dTimeSinceStartMin, uAsicInFeb, 1.0 / 60.0);
2148 fhStsFebHitRateEvoLong[uFebIdx]->Fill(dTimeSinceStartMin, 1.0 / 60.0);
2149 if (mess.IsHitMissedEvts()) {
2150 fhStsFebChanMissEvtEvo[uFebIdx]->Fill(dTimeSinceStartSec, uChanInFeb);
2151 fhStsFebAsicMissEvtEvo[uFebIdx]->Fill(dTimeSinceStartSec, uAsicInFeb);
2152 fhStsFebMissEvtEvo[uFebIdx]->Fill(dTimeSinceStartSec);
2153 } // if( mess.IsHitMissedEvts() )
2154 /*
2155 if( kTRUE == fbLongHistoEnable )
2156 {
2157 std::chrono::steady_clock::time_point tNow = std::chrono::steady_clock::now();
2158 Double_t dUnixTimeInRun = std::chrono::duration_cast< std::chrono::seconds >(tNow - ftStartTimeUnix).count();
2159 fhFebRateEvoLong[ uAsicIdx ]->Fill( dUnixTimeInRun , 1.0 / fuLongHistoBinSizeSec );
2160 fhFebChRateEvoLong[ uAsicIdx ]->Fill( dUnixTimeInRun , usChan, 1.0 / fuLongHistoBinSizeSec );
2161 } // if( kTRUE == fbLongHistoEnable )
2162*/
2163}
2164
2165void CbmMcbm2018MonitorSts::FillTsMsbInfo(stsxyter::Message mess, UInt_t uMessIdx, UInt_t uMsIdx)
2166{
2167 UInt_t uVal = mess.GetTsMsbVal();
2169 if (kTRUE == fbBinningFw) uVal = mess.GetTsMsbValBinning();
2170
2171 /*
2172 if( ( 419369 < fulCurrentTsIdx && fulCurrentTsIdx < 419371 ) )
2173 LOG(info) << " TS " << std::setw( 12 ) << fulCurrentTsIdx
2174 << " MS " << std::setw( 12 ) << fulCurrentMsIdx
2175 << " MsInTs " << std::setw( 3 ) << uMsIdx
2176 << " DPB " << std::setw( 2 ) << fuCurrDpbIdx
2177 << " Mess " << std::setw( 5 ) << uMessIdx
2178 << " TsMsb " << std::setw( 5 ) << uVal;
2179*/
2180 /*
2181 if( (uVal != fvulCurrentTsMsb[fuCurrDpbIdx] + 1) && 0 < uVal &&
2182 !( 1 == uMessIdx && usVal == fvulCurrentTsMsb[fuCurrDpbIdx] ) ) // 1st TS_MSB in MS is always a repeat of the last one in previous MS!
2183 {
2184 LOG(info) << "TS MSB not increasing by 1! TS " << std::setw( 12 ) << fulCurrentTsIdx
2185 << " MS " << std::setw( 12 ) << fulCurrentMsIdx
2186 << " MsInTs " << std::setw( 3 ) << uMsIdx
2187 << " DPB " << std::setw( 2 ) << fuCurrDpbIdx
2188 << " Mess " << std::setw( 5 ) << uMessIdx
2189 << " Old TsMsb " << std::setw( 5 ) << fvulCurrentTsMsb[fuCurrDpbIdx]
2190 << " new TsMsb " << std::setw( 5 ) << uVal
2191 << " Diff " << std::setw( 5 ) << uVal - fvulCurrentTsMsb[fuCurrDpbIdx]
2192 << " Old MsbCy " << std::setw( 5 ) << fvuCurrentTsMsbCycle[fuCurrDpbIdx];
2193 } // if( (uVal != fvulCurrentTsMsb[fuCurrDpbIdx] + 1) && 0 < uVal )
2194*/
2195
2196 // Update Status counters
2197 if (uVal < fvulCurrentTsMsb[fuCurrDpbIdx]) {
2198
2199 LOG(info) << " TS " << std::setw(12) << fulCurrentTsIdx << " MS " << std::setw(12) << fulCurrentMsIdx << " MS Idx "
2200 << std::setw(4) << uMsIdx << " Msg Idx " << std::setw(5) << uMessIdx << " DPB " << std::setw(2)
2201 << fuCurrDpbIdx << " Old TsMsb " << std::setw(5) << fvulCurrentTsMsb[fuCurrDpbIdx] << " Old MsbCy "
2202 << std::setw(5) << fvuCurrentTsMsbCycle[fuCurrDpbIdx] << " new TsMsb " << std::setw(5) << uVal;
2203
2205 } // if( uVal < fvulCurrentTsMsb[fuCurrDpbIdx] )
2206 if (
2207 uVal != fvulCurrentTsMsb[fuCurrDpbIdx] + 1 && 0 != uVal && 4194303 != fvulCurrentTsMsb[fuCurrDpbIdx]
2208 && 1 != uMessIdx &&
2209 !(0 == uVal && 0 == fvulCurrentTsMsb[fuCurrDpbIdx] && 2 == uMessIdx) &&
2210 !(uVal == fvulCurrentTsMsb[fuCurrDpbIdx] && 2 == uMessIdx)
2211 &&
2212 uVal < fvulCurrentTsMsb
2213 [fuCurrDpbIdx]
2214 ) {
2215 LOG(info) << "TS MSb Jump in "
2216 << " TS " << std::setw(12) << fulCurrentTsIdx << " MS " << std::setw(12) << fulCurrentMsIdx << " MS Idx "
2217 << std::setw(4) << uMsIdx << " Msg Idx " << std::setw(5) << uMessIdx << " DPB " << std::setw(2)
2218 << fuCurrDpbIdx << " => Old TsMsb " << std::setw(5) << fvulCurrentTsMsb[fuCurrDpbIdx] << " new TsMsb "
2219 << std::setw(5) << uVal;
2220 } // if( uVal + 1 != fvulCurrentTsMsb[fuCurrDpbIdx] && 4194303 != uVal && 0 != fvulCurrentTsMsb[fuCurrDpbIdx] )
2222 /*
2223 if( 1 < uMessIdx )
2224 {
2225 fhStsDpbRawTsMsb->Fill( fuCurrDpbIdx, fvulCurrentTsMsb[fuCurrDpbIdx] );
2226 fhStsDpbRawTsMsbSx->Fill( fuCurrDpbIdx, ( fvulCurrentTsMsb[fuCurrDpbIdx] & 0x1F ) );
2227 fhStsDpbRawTsMsbDpb->Fill( fuCurrDpbIdx, ( fvulCurrentTsMsb[fuCurrDpbIdx] >> 5 ) );
2228 } // if( 0 < uMessIdx )
2229*/
2230 // fhStsAsicTsMsb->Fill( fvulCurrentTsMsb[fuCurrDpbIdx], uAsicIdx );
2231 /*
2232 ULong64_t ulNewTsMsbTime = static_cast< ULong64_t >( stsxyter::kuHitNbTsBins )
2233 * static_cast< ULong64_t >( fvulCurrentTsMsb[fuCurrDpbIdx])
2234 + static_cast< ULong64_t >( stsxyter::kulTsCycleNbBins )
2235 * static_cast< ULong64_t >( fvuCurrentTsMsbCycle[fuCurrDpbIdx] );
2236*/
2237}
2238
2240{
2241 // UInt_t uVal = mess.GetEpochVal();
2242 // UInt_t uCurrentCycle = uVal % stsxyter::kulTsCycleNbBins;
2243
2244 /*
2245 // Update Status counters
2246 if( usVal < fvulCurrentTsMsb[fuCurrDpbIdx] )
2247 fvuCurrentTsMsbCycle[fuCurrDpbIdx] ++;
2248 fvulCurrentTsMsb[fuCurrDpbIdx] = usVal;
2249
2250// fhStsAsicTsMsb->Fill( fvulCurrentTsMsb[fuCurrDpbIdx], uAsicIdx );
2251*/
2252}
2253
2255
2257{
2258
2259 LOG(info) << "-------------------------------------";
2260 LOG(info) << "CbmMcbm2018MonitorSts statistics are ";
2261 LOG(info) << " Hit messages: " << fmMsgCounter[stsxyter::MessType::Hit] << "\n"
2262 << " Ts MSB messages: " << fmMsgCounter[stsxyter::MessType::TsMsb] << "\n"
2263 << " Dummy messages: " << fmMsgCounter[stsxyter::MessType::Dummy] << "\n"
2264 << " Epoch messages: " << fmMsgCounter[stsxyter::MessType::Epoch] << "\n"
2265 << " Empty messages: " << fmMsgCounter[stsxyter::MessType::Empty];
2266
2267 LOG(info) << "-------------------------------------";
2268
2270 // SaveAllHistos();
2271}
2272
2273
2275{
2277 TFile* oldFile = gFile;
2278 TDirectory* oldDir = gDirectory;
2279
2280 TFile* histoFile = NULL;
2281 if ("" != sFileName) {
2282 // open separate histo file in recreate mode
2283 histoFile = new TFile(sFileName, "RECREATE");
2284 histoFile->cd();
2285 } // if( "" != sFileName )
2286
2287 /***************************/
2288 gDirectory->mkdir("Sts_Raw");
2289 gDirectory->cd("Sts_Raw");
2290
2291 fhStsMessType->Write();
2292 fhStsSysMessType->Write();
2293 fhStsMessTypePerDpb->Write();
2294 fhStsSysMessTypePerDpb->Write();
2295 fhStsStatusMessType->Write();
2296 fhStsMsStatusFieldType->Write();
2297 fhStsMessTypePerElink->Write();
2298 fhStsHitsElinkPerDpb->Write();
2299 fhStsAllFebsHitRateEvo->Write();
2300 fhStsAllAsicsHitRateEvo->Write();
2301 fhStsFebAsicHitCounts->Write();
2302 gDirectory->cd("..");
2303 /***************************/
2304
2305 /***************************/
2306 gDirectory->mkdir("Sts_Feb");
2307 gDirectory->cd("Sts_Feb");
2308 for (UInt_t uFebIdx = 0; uFebIdx < fuNbFebs; ++uFebIdx) {
2309 if (kTRUE == fUnpackParSts->IsFebActive(uFebIdx)) {
2310 fhStsFebChanCntRaw[uFebIdx]->Write();
2311 fhStsFebChanCntRawGood[uFebIdx]->Write();
2312 fhStsFebChanAdcRaw[uFebIdx]->Write();
2313 fhStsFebChanAdcRawProf[uFebIdx]->Write();
2314 // fhStsFebChanAdcCal[ uFebIdx ]->Write();
2315 // fhStsFebChanAdcCalProf[ uFebIdx ]->Write();
2316 fhStsFebChanRawTs[uFebIdx]->Write();
2317 fhStsFebChanMissEvt[uFebIdx]->Write();
2318 fhStsFebChanMissEvtEvo[uFebIdx]->Write();
2319 fhStsFebAsicMissEvtEvo[uFebIdx]->Write();
2320 fhStsFebMissEvtEvo[uFebIdx]->Write();
2321 fhStsFebChanHitRateEvo[uFebIdx]->Write();
2322 fhStsFebChanHitRateProf[uFebIdx]->Write();
2323 fhStsFebAsicHitRateEvo[uFebIdx]->Write();
2324 fhStsFebHitRateEvo[uFebIdx]->Write();
2325 fhStsFebChanHitRateEvoLong[uFebIdx]->Write();
2326 fhStsFebAsicHitRateEvoLong[uFebIdx]->Write();
2327 fhStsFebHitRateEvoLong[uFebIdx]->Write();
2328 fhStsFebChanDistT[uFebIdx]->Write();
2329 fhStsFebChanCloseHitsCounts[uFebIdx]->Write();
2330 fhStsFebChanCloseHitsRatio[uFebIdx]->Write();
2331 /*
2332 for( UInt_t uFebIdxB = uFebIdx; uFebIdxB < fuNbFebs; ++uFebIdxB )
2333 {
2334 fhStsFebChanDtCoinc[ uFebIdx ][ uFebIdxB ]->Write();
2335 fhStsFebChanCoinc[ uFebIdx ][ uFebIdxB ]->Write();
2336 } // for( UInt_t uFebIdxB = uFebIdx; uFebIdxB < fuNbFebs; ++uFebIdxB )
2337*/
2338 } // if( kTRUE == fUnpackParSts->IsFebActive( uFebIdx ) )
2339 } // for( UInt_t uFebIdx = 0; uFebIdx < fuNbFebs; ++uFebIdx )
2340 gDirectory->cd("..");
2341 /***************************/
2342
2343 /***************************/
2344 /*
2345 gDirectory->mkdir("Sts_Module");
2346 gDirectory->cd("Sts_Module");
2347 for( UInt_t uModIdx = 0; uModIdx < fuNbModules; ++ uModIdx )
2348 {
2349 fhStsModulePNCoincDt[ uModIdx ]->Write();
2350 fhStsModulePNCoincDtAsicP[ uModIdx ]->Write();
2351 fhStsModulePNCoincDtAsicN[ uModIdx ]->Write();
2352 fhStsModulePNCoincChan[ uModIdx ]->Write();
2353 fhStsModulePNCoincAdc[ uModIdx ]->Write();
2354 fhStsModuleCoincAdcChanP[ uModIdx ]->Write();
2355 fhStsModuleCoincAdcChanN[ uModIdx ]->Write();
2356 fhStsModuleCoincMap[ uModIdx ]->Write();
2357 } // for( UInt_t uModIdx = 0; uModIdx < fuNbModules; ++ uModIdx )
2358 gDirectory->cd("..");
2359*/
2360 /***************************/
2361
2363 if (kTRUE == fbEnableCheckBugSmx20) {
2364 gDirectory->mkdir("Sts_SmxErr");
2365 gDirectory->cd("Sts_SmxErr");
2366 for (UInt_t uFebIdx = 0; uFebIdx < fuNbFebs; ++uFebIdx) {
2367 if (kTRUE == fUnpackParSts->IsFebActive(uFebIdx)) {
2368 fhStsFebSmxErrRatioEvo[uFebIdx]->Write();
2369 fhStsFebSmxErrRatioEvoAsic[uFebIdx]->Write();
2370 fhStsFebSmxErrRatioCopyEvo[uFebIdx]->Write();
2371 fhStsFebSmxErrRatioCopyEvoAsic[uFebIdx]->Write();
2372 fhStsFebSmxErrRatioCopySameAdcEvo[uFebIdx]->Write();
2373 fhStsFebSmxErrRatioCopySameAdcEvoAsic[uFebIdx]->Write();
2374 } // if( kTRUE == fUnpackParSts->IsFebActive( uFebIdx ) )
2375 } // for( UInt_t uFebIdx = 0; uFebIdx < fuNbFebs; ++uFebIdx )
2376 gDirectory->cd("..");
2377 } // if( kTRUE == fbEnableCheckBugSmx20 )
2379
2380 /***************************/
2381 // Flib Histos
2382 gDirectory->mkdir("Flib_Raw");
2383 gDirectory->cd("Flib_Raw");
2384 for (UInt_t uLinks = 0; uLinks < kiMaxNbFlibLinks; uLinks++) {
2385 TString sMsSzName = Form("MsSz_link_%02u", uLinks);
2386 if (fHM->Exists(sMsSzName.Data())) fHM->H1(sMsSzName.Data())->Write();
2387
2388 sMsSzName = Form("MsSzTime_link_%02u", uLinks);
2389 if (fHM->Exists(sMsSzName.Data())) fHM->P1(sMsSzName.Data())->Write();
2390 } // for( UInt_t uLinks = 0; uLinks < 16; uLinks ++)
2391
2392 TH1* pMissedTsH1 = dynamic_cast<TH1*>(gROOT->FindObjectAny("Missed_TS"));
2393 if (NULL != pMissedTsH1) pMissedTsH1->Write();
2394
2395 TProfile* pMissedTsEvoP = dynamic_cast<TProfile*>(gROOT->FindObjectAny("Missed_TS_Evo"));
2396 if (NULL != pMissedTsEvoP) pMissedTsEvoP->Write();
2397
2398 gDirectory->cd("..");
2399 /***************************/
2400 for (UInt_t uFebIdx = 0; uFebIdx < fuNbFebs; ++uFebIdx) {
2401 if (kTRUE == fUnpackParSts->IsFebActive(uFebIdx)) {
2402 fvcStsSumm[uFebIdx]->Write();
2403 if (kTRUE == fbEnableCheckBugSmx20) fvcStsSmxErr[uFebIdx]->Write();
2404 } // if( kTRUE == fUnpackParSts->IsFebActive( uFebIdx ) )
2405 } // for( UInt_t uFebIdx = 0; uFebIdx < fuNbFebs; ++uFebIdx )
2406 /***************************/
2407
2408 if ("" != sFileName) {
2409 // Restore original directory position
2410 histoFile->Close();
2411 } // if( "" != sFileName )
2412
2414 gFile = oldFile;
2415 gDirectory = oldDir;
2416}
2418{
2419 LOG(info) << "Reseting all STS histograms.";
2420
2421 fhStsMessType->Reset();
2422 fhStsSysMessType->Reset();
2423 fhStsMessTypePerDpb->Reset();
2424 fhStsSysMessTypePerDpb->Reset();
2425 fhStsStatusMessType->Reset();
2426 fhStsMsStatusFieldType->Reset();
2427 fhStsMessTypePerElink->Reset();
2428 fhStsHitsElinkPerDpb->Reset();
2429 fhStsAllFebsHitRateEvo->Reset();
2430 fhStsAllAsicsHitRateEvo->Reset();
2431 fhStsFebAsicHitCounts->Reset();
2432
2433 for (UInt_t uFebIdx = 0; uFebIdx < fuNbFebs; ++uFebIdx) {
2434 if (kTRUE == fUnpackParSts->IsFebActive(uFebIdx)) {
2435 fhStsFebChanCntRaw[uFebIdx]->Reset();
2436 fhStsFebChanCntRawGood[uFebIdx]->Reset();
2437 fhStsFebChanAdcRaw[uFebIdx]->Reset();
2438 fhStsFebChanAdcRawProf[uFebIdx]->Reset();
2439 // fhStsFebChanAdcCal[ uFebIdx ]->Reset();
2440 // fhStsFebChanAdcCalProf[ uFebIdx ]->Reset();
2441 fhStsFebChanRawTs[uFebIdx]->Reset();
2442 fhStsFebChanMissEvt[uFebIdx]->Reset();
2443 fhStsFebChanMissEvtEvo[uFebIdx]->Reset();
2444 fhStsFebAsicMissEvtEvo[uFebIdx]->Reset();
2445 fhStsFebMissEvtEvo[uFebIdx]->Reset();
2446 fhStsFebChanHitRateEvo[uFebIdx]->Reset();
2447 fhStsFebChanHitRateProf[uFebIdx]->Reset();
2448 fhStsFebAsicHitRateEvo[uFebIdx]->Reset();
2449 fhStsFebHitRateEvo[uFebIdx]->Reset();
2450 fhStsFebChanHitRateEvoLong[uFebIdx]->Reset();
2451 fhStsFebAsicHitRateEvoLong[uFebIdx]->Reset();
2452 fhStsFebHitRateEvoLong[uFebIdx]->Reset();
2453 fhStsFebChanDistT[uFebIdx]->Reset();
2454 fhStsFebChanCloseHitsCounts[uFebIdx]->Reset();
2455 fhStsFebChanCloseHitsRatio[uFebIdx]->Reset();
2456 /*
2457 for( UInt_t uFebIdxB = uFebIdx; uFebIdxB < fuNbFebs; ++uFebIdxB )
2458 {
2459 fhStsFebChanDtCoinc[ uFebIdx ][ uFebIdxB ]->Reset();
2460 fhStsFebChanCoinc[ uFebIdx ][ uFebIdxB ]->Reset();
2461 } // for( UInt_t uFebIdxB = uFebIdx; uFebIdxB < fuNbFebs; ++uFebIdxB )
2462*/
2464 if (kTRUE == fbEnableCheckBugSmx20) {
2465 fhStsFebSmxErrRatioEvo[uFebIdx]->Reset();
2466 fhStsFebSmxErrRatioEvoAsic[uFebIdx]->Reset();
2467 fhStsFebSmxErrRatioCopyEvo[uFebIdx]->Reset();
2468 fhStsFebSmxErrRatioCopyEvoAsic[uFebIdx]->Reset();
2469 fhStsFebSmxErrRatioCopySameAdcEvo[uFebIdx]->Reset();
2470 fhStsFebSmxErrRatioCopySameAdcEvoAsic[uFebIdx]->Reset();
2471 } // if( kTRUE == fbEnableCheckBugSmx20 )
2473 } // if( kTRUE == fUnpackParSts->IsFebActive( uFebIdx ) )
2474 } // for( UInt_t uFebIdx = 0; uFebIdx < fuNbFebs; ++uFebIdx )
2475 /*
2476 for( UInt_t uModIdx = 0; uModIdx < fuNbModules; ++ uModIdx )
2477 {
2478 fhStsModulePNCoincDt[ uModIdx ]->Reset();
2479 fhStsModulePNCoincDtAsicP[ uModIdx ]->Reset();
2480 fhStsModulePNCoincDtAsicN[ uModIdx ]->Reset();
2481 fhStsModulePNCoincChan[ uModIdx ]->Reset();
2482 fhStsModulePNCoincAdc[ uModIdx ]->Reset();
2483 fhStsModuleCoincAdcChanP[ uModIdx ]->Reset();
2484 fhStsModuleCoincAdcChanN[ uModIdx ]->Reset();
2485 fhStsModuleCoincMap[ uModIdx ]->Reset();
2486 } // for( UInt_t uModIdx = 0; uModIdx < fuNbModules; ++ uModIdx )
2487*/
2488 for (UInt_t uLinks = 0; uLinks < kiMaxNbFlibLinks; ++uLinks) {
2489 TString sMsSzName = Form("MsSz_link_%02u", uLinks);
2490 if (fHM->Exists(sMsSzName.Data())) fHM->H1(sMsSzName.Data())->Reset();
2491
2492 sMsSzName = Form("MsSzTime_link_%02u", uLinks);
2493 if (fHM->Exists(sMsSzName.Data())) fHM->P1(sMsSzName.Data())->Reset();
2494 } // for( UInt_t uLinks = 0; uLinks < kiMaxNbFlibLinks; ++uLinks )
2495
2496 fdStartTime = -1;
2497 fdStartTimeMsSz = -1;
2498}
2499/*
2500void CbmMcbm2018MonitorSts::UpdatePairMeanValue( UInt_t uAsicA, UInt_t uAsicB, Double_t dNewValue )
2501{
2502 if( fvuLastTimeDiffSlotAsicPair.size() < uAsicA )
2503 {
2504 LOG(warning) << "CbmMcbm2018MonitorSts::UpdatePairMeanValue => wrong ASIC A value: " << uAsicA
2505 << " VS " << fvuLastTimeDiffSlotAsicPair.size();
2506 return;
2507 } // if( fvuLastTimeDiffSlotAsicPair.size() < uAsicA )
2508 if( fvuLastTimeDiffSlotAsicPair[ uAsicA ].size() < uAsicB )
2509 {
2510 LOG(warning) << "CbmMcbm2018MonitorSts::UpdatePairMeanValue => wrong ASIC B value: " << uAsicB
2511 << " VS " << fvuLastTimeDiffSlotAsicPair[ uAsicA ].size();
2512 return;
2513 } // if( fvuLastTimeDiffSlotAsicPair[ uAsicA ].size() < uAsicB )
2514
2515 if( kuNbValuesForTimeDiffMean == fvdLastTimeDiffValuesAsicPair[ uAsicA ][ uAsicB ].size() )
2516 {
2517 fvdLastTimeDiffValuesAsicPair[ uAsicA ][ uAsicB ][ fvuLastTimeDiffSlotAsicPair[ uAsicA ][ uAsicB ] ] = dNewValue;
2518 fvuLastTimeDiffSlotAsicPair[ uAsicA ][ uAsicB ] = ( fvuLastTimeDiffSlotAsicPair[ uAsicA ][ uAsicB ] + 1 )
2519 % kuNbValuesForTimeDiffMean;
2520 } // if( kuNbValuesForTimeDiffMean == fvdLastTimeDiffValuesAsicPair[ uAsicA ][ uAsicB ].size() )
2521 else fvdLastTimeDiffValuesAsicPair[ uAsicA ][ uAsicB ].push_back( dNewValue );
2522
2523 Double_t dNewMean = 0.0;
2524 UInt_t uNbVal = fvdLastTimeDiffValuesAsicPair[ uAsicA ][ uAsicB ].size();
2525 for( UInt_t uIdx = 0; uIdx < uNbVal; ++uIdx )
2526 dNewMean += fvdLastTimeDiffValuesAsicPair[ uAsicA ][ uAsicB ][ uIdx ];
2527 dNewMean /= uNbVal;
2528
2529 fvdMeanTimeDiffAsicPair[ uAsicA ][ uAsicB ] = dNewMean;
2530}
2531*/
2533{
2534 TDatime* fRunStartDateTime = new TDatime(dateIn, timeIn);
2535 fiRunStartDateTimeSec = fRunStartDateTime->Convert();
2536 fiBinSizeDatePlots = iBinSize;
2537
2538 LOG(info) << "Assigned new MUCH Run Start Date-Time: " << fRunStartDateTime->AsString();
2539}
2540
2541void CbmMcbm2018MonitorSts::SetLongDurationLimits(UInt_t uDurationSeconds, UInt_t uBinSize)
2542{
2543 fbLongHistoEnable = kTRUE;
2544 fuLongHistoNbSeconds = uDurationSeconds;
2545 fuLongHistoBinSizeSec = uBinSize;
2546}
2547
2551{
2552
2553 if (kTRUE == fbSmx2ErrorUseNoiseLevels) {
2559 } // if( kTRUE == fbSmx2ErrorUseNoiseLevels )
2560
2568 for (UInt_t uFebIdx = 0; uFebIdx < fuNbFebs; ++uFebIdx) {
2569 fvdSmxErrTimeLastHits[uFebIdx].resize(fUnpackParSts->GetNbAsicsPerFeb());
2570 fvuSmxErrIdxFirstHitM07[uFebIdx].resize(fUnpackParSts->GetNbAsicsPerFeb());
2571 fvuSmxErrIdxFirstHitM08[uFebIdx].resize(fUnpackParSts->GetNbAsicsPerFeb());
2572 fvuSmxErrIdxFirstHitM09[uFebIdx].resize(fUnpackParSts->GetNbAsicsPerFeb());
2573 fvuSmxErrIdxFirstHitM10[uFebIdx].resize(fUnpackParSts->GetNbAsicsPerFeb());
2574 fvuSmxErrIdxFirstHitM11[uFebIdx].resize(fUnpackParSts->GetNbAsicsPerFeb());
2575 fvuSmxErrIdxLastHit[uFebIdx].resize(fUnpackParSts->GetNbAsicsPerFeb());
2576 for (UInt_t uAsicIdx = 0; uAsicIdx < fUnpackParSts->GetNbAsicsPerFeb(); ++uAsicIdx) {
2577 //fvulLastHitTs[ uXyterIdx ].resize( fuNbChanPerAsic, 0 );
2578 fvdSmxErrTimeLastHits[uFebIdx][uAsicIdx].resize(kuSmxErrCoincWinNbHits, -1.0);
2584 fvuSmxErrIdxLastHit[uFebIdx][uAsicIdx] = kuSmxErrCoincWinNbHits;
2585 } // for( UInt_t uAsicIdx = 0; uAsicIdx < fUnpackParSts->GetNbAsicsPerFeb(); ++uAsicIdx )
2586 } // for( UInt_t uFebIdx = 0; uFebIdx < fuNbFebs; ++uFebIdx )
2587}
2588
2589Bool_t CbmMcbm2018MonitorSts::SmxErrCheckCoinc(UInt_t uFebIdx, UInt_t uAsicIdx, Double_t dNewHitTime)
2590{
2591 if (kuSmxErrCoincWinNbHits == fvuSmxErrIdxFirstHitM07[uFebIdx][uAsicIdx]
2592 && kuSmxErrCoincWinNbHits == fvuSmxErrIdxFirstHitM08[uFebIdx][uAsicIdx]
2593 && kuSmxErrCoincWinNbHits == fvuSmxErrIdxFirstHitM09[uFebIdx][uAsicIdx]
2594 && kuSmxErrCoincWinNbHits == fvuSmxErrIdxFirstHitM10[uFebIdx][uAsicIdx]
2595 && kuSmxErrCoincWinNbHits == fvuSmxErrIdxFirstHitM11[uFebIdx][uAsicIdx]
2596 && kuSmxErrCoincWinNbHits == fvuSmxErrIdxLastHit[uFebIdx][uAsicIdx]) {
2598 fvuSmxErrIdxFirstHitM07[uFebIdx][uAsicIdx] = 4;
2599 fvuSmxErrIdxFirstHitM08[uFebIdx][uAsicIdx] = 3;
2600 fvuSmxErrIdxFirstHitM09[uFebIdx][uAsicIdx] = 2;
2601 fvuSmxErrIdxFirstHitM10[uFebIdx][uAsicIdx] = 1;
2602 fvuSmxErrIdxFirstHitM11[uFebIdx][uAsicIdx] = 0;
2603 fvuSmxErrIdxLastHit[uFebIdx][uAsicIdx] = 0;
2604 fvdSmxErrTimeLastHits[uFebIdx][uAsicIdx][fvuSmxErrIdxLastHit[uFebIdx][uAsicIdx]] = dNewHitTime;
2605
2607 return kFALSE;
2608 }
2609 else if (kuSmxErrCoincWinNbHits - 1
2610 == fvuSmxErrIdxLastHit[uFebIdx][uAsicIdx] - fvuSmxErrIdxFirstHitM11[uFebIdx][uAsicIdx]
2611 || fvuSmxErrIdxFirstHitM11[uFebIdx][uAsicIdx] - 1 == fvuSmxErrIdxLastHit[uFebIdx][uAsicIdx]) {
2613 fvuSmxErrIdxFirstHitM07[uFebIdx][uAsicIdx] =
2614 (fvuSmxErrIdxFirstHitM07[uFebIdx][uAsicIdx] + 1) % kuSmxErrCoincWinNbHits;
2615 fvuSmxErrIdxFirstHitM08[uFebIdx][uAsicIdx] =
2616 (fvuSmxErrIdxFirstHitM08[uFebIdx][uAsicIdx] + 1) % kuSmxErrCoincWinNbHits;
2617 fvuSmxErrIdxFirstHitM09[uFebIdx][uAsicIdx] =
2618 (fvuSmxErrIdxFirstHitM09[uFebIdx][uAsicIdx] + 1) % kuSmxErrCoincWinNbHits;
2619 fvuSmxErrIdxFirstHitM10[uFebIdx][uAsicIdx] =
2620 (fvuSmxErrIdxFirstHitM10[uFebIdx][uAsicIdx] + 1) % kuSmxErrCoincWinNbHits;
2621 fvuSmxErrIdxFirstHitM11[uFebIdx][uAsicIdx] =
2622 (fvuSmxErrIdxFirstHitM11[uFebIdx][uAsicIdx] + 1) % kuSmxErrCoincWinNbHits;
2623 fvuSmxErrIdxLastHit[uFebIdx][uAsicIdx] = (fvuSmxErrIdxLastHit[uFebIdx][uAsicIdx] + 1) % kuSmxErrCoincWinNbHits;
2624 fvdSmxErrTimeLastHits[uFebIdx][uAsicIdx][fvuSmxErrIdxLastHit[uFebIdx][uAsicIdx]] = dNewHitTime;
2625 }
2626 else {
2628 fvuSmxErrIdxLastHit[uFebIdx][uAsicIdx] = fvuSmxErrIdxLastHit[uFebIdx][uAsicIdx] + 1;
2629 fvdSmxErrTimeLastHits[uFebIdx][uAsicIdx][fvuSmxErrIdxLastHit[uFebIdx][uAsicIdx]] = dNewHitTime;
2630
2632 return kFALSE;
2633 }
2634
2636 Double_t dTimeDiff07 = fvdSmxErrTimeLastHits[uFebIdx][uAsicIdx][fvuSmxErrIdxLastHit[uFebIdx][uAsicIdx]]
2637 - fvdSmxErrTimeLastHits[uFebIdx][uAsicIdx][fvuSmxErrIdxFirstHitM07[uFebIdx][uAsicIdx]];
2638 Double_t dTimeDiff08 = fvdSmxErrTimeLastHits[uFebIdx][uAsicIdx][fvuSmxErrIdxLastHit[uFebIdx][uAsicIdx]]
2639 - fvdSmxErrTimeLastHits[uFebIdx][uAsicIdx][fvuSmxErrIdxFirstHitM08[uFebIdx][uAsicIdx]];
2640 Double_t dTimeDiff09 = fvdSmxErrTimeLastHits[uFebIdx][uAsicIdx][fvuSmxErrIdxLastHit[uFebIdx][uAsicIdx]]
2641 - fvdSmxErrTimeLastHits[uFebIdx][uAsicIdx][fvuSmxErrIdxFirstHitM09[uFebIdx][uAsicIdx]];
2642 Double_t dTimeDiff10 = fvdSmxErrTimeLastHits[uFebIdx][uAsicIdx][fvuSmxErrIdxLastHit[uFebIdx][uAsicIdx]]
2643 - fvdSmxErrTimeLastHits[uFebIdx][uAsicIdx][fvuSmxErrIdxFirstHitM10[uFebIdx][uAsicIdx]];
2644 Double_t dTimeDiff11 = fvdSmxErrTimeLastHits[uFebIdx][uAsicIdx][fvuSmxErrIdxLastHit[uFebIdx][uAsicIdx]]
2645 - fvdSmxErrTimeLastHits[uFebIdx][uAsicIdx][fvuSmxErrIdxFirstHitM11[uFebIdx][uAsicIdx]];
2646
2648 /*
2649 if( 3200 <= dTimeDiff07 && dTimeDiff07 <= 4500 )
2650 dTimeDiff7 -= 3200;
2651 if( 3200 <= dTimeDiff08 && dTimeDiff08 <= 4500 )
2652 dTimeDiff8 -= 3200;
2653 if( 3200 <= dTimeDiff09 && dTimeDiff09 <= 4500 )
2654 dTimeDiff9 -= 3200;
2655 if( 3200 <= dTimeDiff10 && dTimeDiff10 <= 4500 )
2656 dTimeDiff10 -= 3200;
2657 if( 3200 <= dTimeDiff11 && dTimeDiff11 <= 4500 )
2658 dTimeDiff11 -= 3200;
2659*/
2660
2662 if ((kdSmxErrCoincWinBeg <= dTimeDiff07 && dTimeDiff07 <= fdSmxErrCoincWinM07)
2663 || (kdSmxErrCoincWinBeg <= dTimeDiff08 && dTimeDiff08 <= fdSmxErrCoincWinM08)
2664 || (kdSmxErrCoincWinBeg <= dTimeDiff09 && dTimeDiff09 <= fdSmxErrCoincWinM09)
2665 || (kdSmxErrCoincWinBeg <= dTimeDiff10 && dTimeDiff10 <= fdSmxErrCoincWinM10)
2666 || (kdSmxErrCoincWinBeg <= dTimeDiff11 && dTimeDiff11 <= fdSmxErrCoincWinM11)) {
2667 return kTRUE;
2668 }
2669 else {
2670 return kFALSE;
2671 }
2672}
2673
2675{
2676 for (UInt_t uFebIdx = 0; uFebIdx < fuNbFebs; ++uFebIdx) {
2677 if (kTRUE == fUnpackParSts->IsFebActive(uFebIdx)) {
2678 LOG(info) << Form(" ------------------ Noisy channels scan for FEB %2d ------------", uFebIdx);
2679 for (UInt_t uAsicIdx = 0; uAsicIdx < fUnpackParSts->GetNbAsicsPerFeb(); ++uAsicIdx)
2680 for (UInt_t uChanIdx = 0; uChanIdx < fUnpackParSts->GetNbChanPerAsic(); ++uChanIdx) {
2681 UInt_t uChanInFeb = uAsicIdx * fUnpackParSts->GetNbChanPerAsic() + uChanIdx;
2682 if (dNoiseThreshold < fhStsFebChanHitRateProf[uFebIdx]->GetBinContent(1 + uChanInFeb))
2683 LOG(info) << Form("Noisy Channel ASIC %d channel %3d (%4d) level %6.0f", uAsicIdx, uChanIdx, uChanInFeb,
2684 fhStsFebChanHitRateProf[uFebIdx]->GetBinContent(1 + uChanInFeb));
2685
2686 } // for( UInt_t uChanIdx = 0; uChanIdx < fUnpackParSts->GetNbChanPerAsic(); ++uChanIdx )
2687
2688 LOG(info) << " ---------------------------------------------------------------";
2689 } // if( kTRUE == fUnpackParSts->IsFebActive( uFebIdx ) )
2690 } // for( UInt_t uFebIdx = 0; uFebIdx < fuNbFebs; ++uFebIdx )
2691 return kTRUE;
2692}
2693
2694
ClassImp(CbmConverterManager)
std::string FormatMsHeaderPrintout(const fles::MicrosliceDescriptor &msDescriptor)
Histogram manager.
Bool_t bMcbm2018ScanNoisySts
Bool_t bMcbm2018WriteSts
Bool_t bMcbm2018ResetSts
int Int_t
bool Bool_t
Histogram manager.
std::vector< std::vector< std::vector< Int_t > > > fviFebModuleIdx
Number of StsXyter ASICs.
std::vector< TProfile2D * > fhStsFebSmxErrRatioCopySameAdcEvoAsic
Int_t fiTimeIntervalRateUpdate
Mean Rate per channel plots.
std::vector< TH2 * > fhStsFebAsicHitRateEvo
std::vector< std::vector< ULong64_t > > fvulChanLastHitTime
std::vector< UInt_t > fvuInitialTsMsbCycleHeader
Flag set after seeing MS header in 1st MS for DPB.
std::vector< std::vector< Double_t > > fvdChanLastHitTime
Last hit time in bins for each Channel.
std::vector< std::vector< std::vector< Int_t > > > fviFebModuleSide
Idx of the STS module for each FEB, [ NbDpb ][ NbCrobPerDpb ][ NbFebsPerCrob ], -1 if inactive.
std::vector< TH2 * > fhStsFebChanDistT
std::vector< ULong64_t > fvulCurrentTsMsb
Bin size in s for the plots with date as X axis.
Bool_t ScanForNoisyChannels(Double_t dNoiseThreshold=1e3)
---------------------------------------------------------------—///
std::vector< std::vector< stsxyter::FinalHit > > fvmAsicHitsInMs
All hits (time in bins, ADC in bins, asic, channel) in last TS, sorted with "<" operator.
std::map< stsxyter::MessType, UInt_t > fmMsgCounter
Current data properties.
std::vector< TProfile * > fhStsFebChanAdcRawProf
std::vector< TH1 * > fhStsFebChanCntRaw
std::vector< Double_t > fvdFebTimeSecLastRateUpdate
static constexpr const Double_t kdSmxErrCoincWinNoiseM11
TH2 * fhMsErrorsEvo
Binning FW error flag.
std::vector< TCanvas * > fvcStsSumm
Canvases.
Bool_t fbLongHistoEnable
Rate evolution histos.
std::vector< TProfile * > fhStsFebChanCloseHitsRatio
TH1 * fhMsSz[kiMaxNbFlibLinks]
std::vector< TH2 * > fhStsFebChanCloseHitsCounts
void SaveAllHistos(TString sFileName="")
std::vector< std::vector< std::vector< Int_t > > > fviFebType
STS module side for each FEB, [ NbDpb ][ NbCrobPerDpb ][ NbFebsPerCrob ], 0 = P, 1 = N,...
Double_t fdCoincCenter
Coincidences in sorted hits.
std::vector< TCanvas * > fvcStsSmxErr
Bool_t fbSmx2ErrorUseNoiseLevels
SXM 2.0 logic error detection and tagging, 1 eLink case.
std::vector< TH1 * > fhStsFebChanCntRawGood
std::vector< UInt_t > fvuInitialHeaderDone
Current TS MSB cycle for DPB.
std::map< UInt_t, UInt_t > fDpbIdIndexMap
Total number of STS DPBs in system.
stsxyter::MessagePrintMask fPrintMessCtrl
virtual void SetNbMsInTs(size_t uCoreMsNb, size_t uOverlapMsNb)
UInt_t fuCurrDpbId
Current equipment ID, tells from which DPB the current MS is originating.
Int_t fiBinSizeDatePlots
Start of run time since "epoch" in s, for the plots with date as X axis.
std::vector< Int_t > fviFebCountsSinceLastRateUpdate
std::vector< std::vector< UInt_t > > fvuSmxErrIdxLastHit
[ NbFebs ][ NbSmxPerFeb ]
std::vector< std::vector< Bool_t > > fvbCrobActiveFlag
Map of DPB Identifier to DPB index.
std::vector< Double_t > fvdPrevMsTime
Last hit time in ns for each Channel.
UInt_t fuNrOfDpbs
STS address for the first strip of each module.
std::vector< TProfile * > fhStsFebSmxErrRatioCopySameAdcEvo
std::vector< TH2 * > fhStsFebAsicHitRateEvoLong
std::vector< std::vector< Double_t > > fdStsFebChanLastTimeForDist
static constexpr const Double_t kdSmxErrCoincWinNoiseM10
std::vector< TProfile2D * > fhStsFebSmxErrRatioCopyEvoAsic
std::vector< std::vector< UInt_t > > fvuSmxErrIdxFirstHitM07
[ NbFebs ][ NbSmxPerFeb ][ kuSmxErrCoincWinNbHits ]
std::vector< std::vector< stsxyter::FinalHit > > fvmFebHitsInMs
All hits (time in bins, ADC in bins, asic, channel) in last TS, per ASIC, sorted with "<" operator.
std::vector< TProfile2D * > fhStsFebSmxErrRatioEvoAsic
static constexpr const Double_t kdSmxErrCoincWinMainM11
std::vector< TH2 * > fhStsFebChanRawTs
Double_t fdStartTime
Last hit ADC in bins in each MS for each Channel.
std::vector< UInt_t > fvuElinkLastTsHit
TS MSB cycle from MS header in 1st MS for DPB.
static constexpr const Double_t kdSmxErrCoincWinNoiseM09
static constexpr const Double_t kdSmxErrCoincWinMainM08
static constexpr const Double_t kdSmxErrCoincWinMainM09
std::vector< TH2 * > fhStsFebChanMissEvt
std::vector< TH1 * > fhStsFebHitRateEvo
UInt_t fuCurrDpbIdx
Temp holder until Current equipment ID is properly filled in MS.
std::vector< Int_t > fviModAddress
Type of each module: 0 for connectors on the right, 1 for connectors on the left.
static constexpr const Double_t kdSmxErrCoincWinMainM10
std::vector< std::vector< UInt_t > > fvuSmxErrIdxFirstHitM08
[ NbFebs ][ NbSmxPerFeb ]
std::vector< TH2 * > fhStsFebChanHitRateEvo
virtual Bool_t DoUnpack(const fles::Timeslice &ts, size_t component)
static constexpr const Double_t kdSmxErrCoincWinBeg
std::vector< std::vector< std::vector< Double_t > > > fvdChanLastHitTimeInMs
Number of hits in each MS for each Channel.
std::vector< std::vector< std::vector< Double_t > > > fvdFebAdcOffs
ADC gain in e-/b, [ NbDpb ][ NbCrobPerDpb ][ NbFebsPerCrob ].
std::vector< std::vector< std::vector< Double_t > > > fvdFebAdcGain
FEB type, [ NbDpb ][ NbCrobPerDpb ][ NbFebsPerCrob ], 0 = A, 1 = B, -1 if inactive.
static constexpr const Double_t kdSmxErrCoincWinNoiseM07
Coincidence windows 99.9% tagging (up to 0.1% of corruption not detected)
static const UInt_t kuBytesPerMessage
TH1 * fhStsMessType
Histogram manager.
std::vector< TProfile * > fhStsFebSmxErrRatioCopyEvo
static constexpr const Double_t kdSmxErrCoincWinNoiseM08
Bool_t ProcessStsMs(const fles::Timeslice &ts, size_t uMsComp, UInt_t uMsIdx)
UInt_t fuNbFebs
Array to hold the active flag for all CROBs, [ NbDpb ][ NbCrobPerDpb ].
std::vector< Int_t > fviModuleType
Total number of STS modules in the setup.
Double_t fdFebChanCoincidenceLimit
Plots per FEB-8.
std::vector< TH2 * > fhStsFebChanMissEvtEvo
std::vector< std::vector< UInt_t > > fvuSmxErrIdxFirstHitM10
[ NbFebs ][ NbSmxPerFeb ]
std::vector< TH2 * > fhStsFebChanAdcRaw
std::vector< TProfile * > fhStsFebChanHitRateProf
ULong64_t fulCurrentTsIdx
TS/MS info.
void SetRunStart(Int_t dateIn, Int_t timeIn, Int_t iBinSize=5)
std::vector< Bool_t > fvbMaskedComponents
Double_t fdSmxErrCoincWinM07
Tagging variables.
std::vector< std::vector< std::vector< Double_t > > > fvdSmxErrTimeLastHits
std::vector< TProfile * > fhStsFebSmxErrRatioEvo
Histograms.
std::vector< std::vector< std::vector< UInt_t > > > fvuChanNbHitsInMs
Header time of each MS.
std::vector< std::vector< Double_t > > fvdFebChanCountsSinceLastRateUpdate
std::vector< TH2 * > fhStsFebChanHitRateEvoLong
void SetLongDurationLimits(UInt_t uDurationSeconds=3600, UInt_t uBinSize=1)
CbmMcbm2018StsPar * fUnpackParSts
/‍** Ignore Overlap Ms: all fuOverlapMsNb MS at the end of timeslice **‍/
std::vector< TH1 * > fhStsFebHitRateEvoLong
Bool_t fbBinningFw
=> Quick and dirty hack for binning FW!!!
static constexpr const Double_t kdSmxErrCoincWinMainM07
Coincidence windows for 99.0% tagging (up to 1% of corruption not detected)
Bool_t SmxErrCheckCoinc(UInt_t uFebIdx, UInt_t uAsicIdx, Double_t dNewHitTime)
std::vector< TH1 * > fhStsFebMissEvtEvo
std::vector< std::vector< std::vector< UShort_t > > > fvusChanLastHitAdcInMs
Last hit time in bins in each MS for each Channel.
UInt_t fuNbStsXyters
Number of StsXyter ASICs.
std::chrono::steady_clock::time_point ftStartTimeUnix
std::vector< stsxyter::FinalHit > fvmHitsInMs
Hits time-sorting.
std::vector< TH2 * > fhStsFebAsicMissEvtEvo
std::vector< UInt_t > fvuCurrentTsMsbCycle
Current TS MSB for each DPB.
TProfile * fhMsSzTime[kiMaxNbFlibLinks]
std::vector< size_t > fvMsComponentsList
virtual void AddMsComponentToList(size_t component, UShort_t usDetectorId)
Bool_t fbPrintMessages
Task configuration values.
std::vector< std::vector< UInt_t > > fvuSmxErrIdxFirstHitM09
[ NbFebs ][ NbSmxPerFeb ]
void FillHitInfo(stsxyter::Message mess, const UShort_t &usElinkIdx, const UInt_t &uAsicIdx, const UInt_t &uMsIdx)
Int_t fiRunStartDateTimeSec
Index of the DPB from which the MS currently unpacked is coming.
void FillEpochInfo(stsxyter::Message mess)
static const Int_t kiMaxNbFlibLinks
ADC offset in e-, [ NbDpb ][ NbCrobPerDpb ][ NbFebsPerCrob ].
static const UInt_t kuSmxErrCoincWinNbHits
std::vector< std::vector< UInt_t > > fvuSmxErrIdxFirstHitM11
[ NbFebs ][ NbSmxPerFeb ]
void FillTsMsbInfo(stsxyter::Message mess, UInt_t uMessIdx=0, UInt_t uMsIdx=0)
std::vector< Double_t > fvdMsTime
Header time of previous MS per link.
CbmHistManager * fHM
Histograms.
Data class with information on a STS local track.
XPU_D bool IsHitMissedEvts() const
For Hit data: Returns Missed event flag (1 bit field)
XPU_D uint32_t GetTsMsbValBinning() const
For TS MSB data: Returns the TS MSB 29 bit field)
XPU_D uint16_t GetHitAdc() const
For Hit data: Returns ADC value (5 bit field)
XPU_D uint16_t GetHitChannel() const
For Hit data: Returns StsXYTER channel number (7 bit field)
bool PrintMess(std::ostream &os, MessagePrintMask ctrl=MessagePrintMask::msg_print_Human, bool bBinning=true) const
XPU_D uint16_t GetLinkIndex() const
For all data: Returns the (global) index of the eLink on which the message was received (n bit field)
XPU_D uint16_t GetStatusStatus() const
For Status data: Returns the Status field from ACK frame (4 bit field)
XPU_D MessType GetMessType() const
Returns the message type, see enum MessType.
XPU_D uint16_t GetLinkIndexHitBinning() const
XPU_D uint16_t GetHitTime() const
For Hit data: Returns timestamp (8 bit field, 2 MSB bits overlap removed)
XPU_D uint32_t GetTsMsbVal() const
For TS MSB data: Returns the TS MSB 22 bit field)
XPU_D uint16_t GetHitTimeBinning() const
XPU_D uint16_t GetStatusLink() const
For Status data: Returns the Link Inedx (9 bit field)
XPU_D uint16_t GetMsErrorType() const
For End of MS data: Returns the MS error type field (2 bit field)
XPU_D bool IsMsErrorFlagOn() const
For End of MS data: Returns the MS error flag (1 bit field)
Hash for CbmL1LinkKey.
static constexpr uint32_t kuHitNbTsBinsBinning
static constexpr uint64_t kulTsCycleNbBins
static constexpr uint64_t kulTsCycleNbBinsBinning
static constexpr uint32_t kuHitNbTsBins
static constexpr uint32_t kuHitNbAdcBins
Status/properties constants.
static constexpr double kdClockCycleNs
MessType
Message types.