CbmRoot
Loading...
Searching...
No Matches
CbmMcbm2018UnpackerAlgoSts.cxx
Go to the documentation of this file.
1/* Copyright (C) 2019-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// ----- CbmMcbm2018UnpackerAlgoSts -----
8// ----- Created 26.01.2019 by P.-A. Loizeau -----
9// ----- -----
10// -----------------------------------------------------------------------------
11
13
15#include "CbmMcbm2018StsPar.h"
16#include "TCanvas.h"
17#include "TH1.h"
18#include "TH2.h"
19#include "TList.h"
20#include "TProfile.h"
21#include "TROOT.h"
22#include "TString.h"
23
24#include <Logger.h>
25
26#include <cstdint>
27#include <fstream>
28#include <iomanip>
29#include <iostream>
30
31// -------------------------------------------------------------------------
34 ,
36 fbMonitorMode(kFALSE)
37 , fbDebugMonitorMode(kFALSE)
39 , fUnpackPar(nullptr)
40 , fuNbModules(0)
43 , fuNrOfDpbs(0)
46 , fuNbFebs(0)
47 , fuNbStsXyters(0)
50 , fviFebType()
52 , fviFebSide()
55 , fdTimeOffsetNs(0.0)
57 , fbUseChannelMask(kFALSE)
59 , fdAdcCut(0)
62 , fdTsStartTime(-1.0)
63 , fdTsStopTimeCore(-1.0)
64 , fdMsTime(-1.0)
65 , fuMsIndex(0)
66 , fmMsgCounter()
68 , fuCurrDpbId(0)
69 , fuCurrDpbIdx(0)
73 , fdStartTime(0.0)
74 , fdStartTimeMsSz(0.0)
75 , ftStartTimeUnix(std::chrono::steady_clock::now())
76 , fvmHitsInMs()
81 , fhDigisTimeInRun(nullptr)
82/*
83 fvhHitsTimeToTriggerRaw(),
84 fvhMessDistributionInMs(),
85 fhEventNbPerTs( nullptr ),
86 fcTimeToTrigRaw( nullptr )
87*/
88{
89}
91{
93 fvmHitsInMs.clear();
94 /*
95 for( UInt_t uDpb = 0; uDpb < fuNrOfDpbs; ++uDpb )
96 {
97 fvmHitsInMs[ uDpb ].clear();
98 } // for( UInt_t uDpb = 0; uDpb < fuNrOfDpbs; ++uDpb )
99*/
100 if (nullptr != fParCList) delete fParCList;
101 if (nullptr != fUnpackPar) delete fUnpackPar;
102}
103
104// -------------------------------------------------------------------------
106{
107 LOG(info) << "Initializing mCBM STS 2019 unpacker algo";
108
109 return kTRUE;
110}
118
119// -------------------------------------------------------------------------
121{
122 LOG(info) << "Init parameter containers for CbmMcbm2018UnpackerAlgoSts";
123 Bool_t initOK = ReInitContainers();
124
125 return initOK;
126}
128{
129 LOG(info) << "**********************************************";
130 LOG(info) << "ReInit parameter containers for CbmMcbm2018UnpackerAlgoSts";
131
132 fUnpackPar = (CbmMcbm2018StsPar*) fParCList->FindObject("CbmMcbm2018StsPar");
133 if (nullptr == fUnpackPar) return kFALSE;
134
135 Bool_t initOK = InitParameters();
136
137 return initOK;
138}
140{
141 if (nullptr == fParCList) fParCList = new TList();
142 fUnpackPar = new CbmMcbm2018StsPar("CbmMcbm2018StsPar");
143 fParCList->Add(fUnpackPar);
144
145 return fParCList;
146}
148{
149 fuNbModules = fUnpackPar->GetNbOfModules();
150 LOG(info) << "Nr. of STS Modules: " << fuNbModules;
151
154 for (UInt_t uModIdx = 0; uModIdx < fuNbModules; ++uModIdx) {
155 fviModuleType[uModIdx] = fUnpackPar->GetModuleType(uModIdx);
156 fviModAddress[uModIdx] = fUnpackPar->GetModuleAddress(uModIdx);
157 LOG(info) << "Module #" << std::setw(2) << uModIdx << " Type " << std::setw(4) << fviModuleType[uModIdx]
158 << " Address 0x" << std::setw(8) << std::hex << fviModAddress[uModIdx] << std::dec;
159 } // for( UInt_t uModIdx = 0; uModIdx < fuNbModules; ++uModIdx)
160
161 fuNrOfDpbs = fUnpackPar->GetNrOfDpbs();
162 LOG(info) << "Nr. of STS DPBs: " << fuNrOfDpbs;
163
164 fDpbIdIndexMap.clear();
165 for (UInt_t uDpb = 0; uDpb < fuNrOfDpbs; ++uDpb) {
166 fDpbIdIndexMap[fUnpackPar->GetDpbId(uDpb)] = uDpb;
167 LOG(info) << "Eq. ID for DPB #" << std::setw(2) << uDpb << " = 0x" << std::setw(4) << std::hex
168 << fUnpackPar->GetDpbId(uDpb) << std::dec << " => " << fDpbIdIndexMap[fUnpackPar->GetDpbId(uDpb)];
169 } // for( UInt_t uDpb = 0; uDpb < fuNrOfDpbs; ++uDpb )
170
171 fuNbFebs = fUnpackPar->GetNrOfFebs();
172 LOG(info) << "Nr. of FEBs: " << fuNbFebs;
173
174 fuNbStsXyters = fUnpackPar->GetNrOfAsics();
175 LOG(info) << "Nr. of StsXyter ASICs: " << fuNbStsXyters;
176
177 if (fvdTimeOffsetNsAsics.size() < fuNbStsXyters) {
179 } // if( fvdTimeOffsetNsAsics.size() < fuNbStsXyters )
180
184 fviFebType.resize(fuNrOfDpbs);
185 for (UInt_t uDpb = 0; uDpb < fuNrOfDpbs; ++uDpb) {
186 fvbCrobActiveFlag[uDpb].resize(fUnpackPar->GetNbCrobsPerDpb());
187 fviFebModuleIdx[uDpb].resize(fUnpackPar->GetNbCrobsPerDpb());
188 fviFebModuleSide[uDpb].resize(fUnpackPar->GetNbCrobsPerDpb());
189 fviFebType[uDpb].resize(fUnpackPar->GetNbCrobsPerDpb());
190 for (UInt_t uCrobIdx = 0; uCrobIdx < fUnpackPar->GetNbCrobsPerDpb(); ++uCrobIdx) {
191 fvbCrobActiveFlag[uDpb][uCrobIdx] = fUnpackPar->IsCrobActive(uDpb, uCrobIdx);
192
193 fviFebModuleIdx[uDpb][uCrobIdx].resize(fUnpackPar->GetNbFebsPerCrob());
194 fviFebModuleSide[uDpb][uCrobIdx].resize(fUnpackPar->GetNbFebsPerCrob());
195 fviFebType[uDpb][uCrobIdx].resize(fUnpackPar->GetNbFebsPerCrob(), -1);
196 for (UInt_t uFebIdx = 0; uFebIdx < fUnpackPar->GetNbFebsPerCrob(); ++uFebIdx) {
197 fviFebModuleIdx[uDpb][uCrobIdx][uFebIdx] = fUnpackPar->GetFebModuleIdx(uDpb, uCrobIdx, uFebIdx);
198 fviFebModuleSide[uDpb][uCrobIdx][uFebIdx] = fUnpackPar->GetFebModuleSide(uDpb, uCrobIdx, uFebIdx);
199
200 fvbFebPulser.push_back(fUnpackPar->IsFebPulser(uDpb, uCrobIdx, uFebIdx));
201 fvdFebAdcGain.push_back(fUnpackPar->GetFebAdcGain(uDpb, uCrobIdx, uFebIdx));
202 fvdFebAdcOffs.push_back(fUnpackPar->GetFebAdcOffset(uDpb, uCrobIdx, uFebIdx));
203
204 if (0 <= fviFebModuleIdx[uDpb][uCrobIdx][uFebIdx]
205 && static_cast<UInt_t>(fviFebModuleIdx[uDpb][uCrobIdx][uFebIdx]) < fuNbModules
206 && 0 <= fviFebModuleSide[uDpb][uCrobIdx][uFebIdx] && fviFebModuleSide[uDpb][uCrobIdx][uFebIdx] < 2) {
207 switch (fviModuleType[fviFebModuleIdx[uDpb][uCrobIdx][uFebIdx]]) {
208 case 0: // FEB-8-1 with ZIF connector on the right
209 {
210 // P side (0) has type A (0)
211 // N side (1) has type B (1)
212 fviFebType[uDpb][uCrobIdx][uFebIdx] = fviFebModuleSide[uDpb][uCrobIdx][uFebIdx];
213
220 fviFebAddress.push_back(fviModAddress[fviFebModuleIdx[uDpb][uCrobIdx][uFebIdx]]
221 + (fviFebModuleSide[uDpb][uCrobIdx][uFebIdx] << 25));
222 fviFebSide.push_back(fviFebModuleSide[uDpb][uCrobIdx][uFebIdx]);
223 break;
224 } // case 0: // FEB-8-1 with ZIF connector on the right
225 case 1: // FEB-8-1 with ZIF connector on the left
226 {
227 // P side (0) has type B (1)
228 // N side (1) has type A (0)
229 fviFebType[uDpb][uCrobIdx][uFebIdx] = !(fviFebModuleSide[uDpb][uCrobIdx][uFebIdx]);
230
237 fviFebAddress.push_back(fviModAddress[fviFebModuleIdx[uDpb][uCrobIdx][uFebIdx]]
238 + ((!fviFebModuleSide[uDpb][uCrobIdx][uFebIdx]) << 25));
239 fviFebSide.push_back(fviFebModuleSide[uDpb][uCrobIdx][uFebIdx]);
240 break;
241 } // case 1: // FEB-8-1 with ZIF connector on the left
242 default:
243 LOG(fatal) << Form("Bad module type for DPB #%02u CROB #%u FEB %02u: %d", uDpb, uCrobIdx, uFebIdx,
244 fviModuleType[fviFebModuleIdx[uDpb][uCrobIdx][uFebIdx]]);
245 break;
246 } // switch( fviModuleType[ fviFebModuleIdx[ uDpb ][ uCrobIdx ][ uFebIdx ] ] )
247 } // FEB active and module index OK
248 else if (-1 == fviFebModuleIdx[uDpb][uCrobIdx][uFebIdx] || -1 == fviFebModuleSide[uDpb][uCrobIdx][uFebIdx]) {
249 fviFebAddress.push_back(0);
250 fviFebSide.push_back(-1);
251 } // Module index or type is set to inactive
252 else {
253 LOG(fatal) << Form("Bad module Index and/or Side for DPB #%02u CROB "
254 "#%u FEB %02u: %d %d",
255 uDpb, uCrobIdx, uFebIdx, fviFebModuleIdx[uDpb][uCrobIdx][uFebIdx],
256 fviFebModuleSide[uDpb][uCrobIdx][uFebIdx]);
257 } // Bad module index or type for this FEB
258 } // for( UInt_t uFebIdx = 0; uFebIdx < fUnpackPar->GetNbFebsPerCrob(); ++ uFebIdx )
259 } // for( UInt_t uCrobIdx = 0; uCrobIdx < fUnpackPar->GetNbCrobsPerDpb(); ++uCrobIdx )
260 } // for( UInt_t uDpb = 0; uDpb < fuNrOfDpbs; ++uDpb )
261
262 for (UInt_t uDpb = 0; uDpb < fuNrOfDpbs; ++uDpb) {
263 TString sPrintoutLine = Form("DPB #%02u CROB Active ?: ", uDpb);
264 for (UInt_t uCrobIdx = 0; uCrobIdx < fUnpackPar->GetNbCrobsPerDpb(); ++uCrobIdx) {
265 sPrintoutLine += Form("%1u", (fvbCrobActiveFlag[uDpb][uCrobIdx] == kTRUE));
266 } // for( UInt_t uCrobIdx = 0; uCrobIdx < fUnpackPar->GetNbCrobsPerDpb(); ++uCrobIdx )
267 LOG(info) << sPrintoutLine;
268 } // for( UInt_t uDpb = 0; uDpb < fuNrOfDpbs; ++uDpb )
269
270 UInt_t uGlobalFebIdx = 0;
271 for (UInt_t uDpb = 0; uDpb < fuNrOfDpbs; ++uDpb) {
272 for (UInt_t uCrobIdx = 0; uCrobIdx < fUnpackPar->GetNbCrobsPerDpb(); ++uCrobIdx) {
273 LOG(info) << Form("DPB #%02u CROB #%u: ", uDpb, uCrobIdx);
274 for (UInt_t uFebIdx = 0; uFebIdx < fUnpackPar->GetNbFebsPerCrob(); ++uFebIdx) {
275 if (0 <= fviFebModuleIdx[uDpb][uCrobIdx][uFebIdx])
276 LOG(info) << Form(" FEB #%02u (%02u): Mod. Idx = %03d Side %c (%2d) Type %c "
277 "(%2d) (Addr. 0x%08x) ADC gain %4.0f e- ADC Offs %5.0f e-",
278 uFebIdx, uGlobalFebIdx, fviFebModuleIdx[uDpb][uCrobIdx][uFebIdx],
279 1 == fviFebModuleSide[uDpb][uCrobIdx][uFebIdx] ? 'N' : 'P',
280 fviFebModuleSide[uDpb][uCrobIdx][uFebIdx],
281 1 == fviFebType[uDpb][uCrobIdx][uFebIdx] ? 'B' : 'A', fviFebType[uDpb][uCrobIdx][uFebIdx],
282 fviFebAddress[uGlobalFebIdx], fvdFebAdcGain[uGlobalFebIdx], fvdFebAdcOffs[uGlobalFebIdx]);
283 else
284 LOG(info) << Form("Disabled FEB #%02u (%02u): Mod. Idx = %03d Side %c (%2d) Type %c "
285 "(%2d) (Addr. 0x%08x) ADC gain %4.0f e- ADC Offs %5.0f e-",
286 uFebIdx, uGlobalFebIdx, fviFebModuleIdx[uDpb][uCrobIdx][uFebIdx],
287 1 == fviFebModuleSide[uDpb][uCrobIdx][uFebIdx] ? 'N' : 'P',
288 fviFebModuleSide[uDpb][uCrobIdx][uFebIdx],
289 1 == fviFebType[uDpb][uCrobIdx][uFebIdx] ? 'B' : 'A', fviFebType[uDpb][uCrobIdx][uFebIdx],
290 fviFebAddress[uGlobalFebIdx], fvdFebAdcGain[uGlobalFebIdx], fvdFebAdcOffs[uGlobalFebIdx]);
291 uGlobalFebIdx++;
292 } // for( UInt_t uFebIdx = 0; uFebIdx < fUnpackPar->GetNbFebsPerCrob(); ++ uFebIdx )
293 } // for( UInt_t uCrobIdx = 0; uCrobIdx < fUnpackPar->GetNbCrobsPerDpb(); ++uCrobIdx )
294 } // for( UInt_t uDpb = 0; uDpb < fuNrOfDpbs; ++uDpb )
295
296 if (fbBinningFw) LOG(info) << "Unpacking data in bin sorter FW mode";
297 else
298 LOG(info) << "Unpacking data in full time sorter FW mode (legacy)";
299
300 // Internal status initialization
303 for (UInt_t uDpb = 0; uDpb < fuNrOfDpbs; ++uDpb) {
304 fvulCurrentTsMsb[uDpb] = 0;
305 fvuCurrentTsMsbCycle[uDpb] = 0;
306 } // for( UInt_t uDpb = 0; uDpb < fuNrOfDpbs; ++uDpb )
307
312 for (UInt_t uAsicIdx = 0; uAsicIdx < fuNbStsXyters; ++uAsicIdx) {
313 fvvusLastTsChan[uAsicIdx].resize(fUnpackPar->GetNbChanPerAsic(), 0);
314 fvvusLastAdcChan[uAsicIdx].resize(fUnpackPar->GetNbChanPerAsic(), 0);
315 fvvusLastTsMsbChan[uAsicIdx].resize(fUnpackPar->GetNbChanPerAsic(), 0);
316 fvvusLastTsMsbCycleChan[uAsicIdx].resize(fUnpackPar->GetNbChanPerAsic(), 0);
317 } // for( UInt_t uAsicIdx = 0; uAsicIdx < fuNbStsXyters; ++ uAsicIdx )
318
319 return kTRUE;
320}
321// -------------------------------------------------------------------------
322
323void CbmMcbm2018UnpackerAlgoSts::AddMsComponentToList(size_t component, UShort_t usDetectorId)
324{
326 for (UInt_t uCompIdx = 0; uCompIdx < fvMsComponentsList.size(); ++uCompIdx)
327 if (component == fvMsComponentsList[uCompIdx]) return;
328
330 fvMsComponentsList.push_back(component);
331
332 LOG(info) << "CbmMcbm2018UnpackerAlgoSts::AddMsComponentToList => Component " << component << " with detector ID 0x"
333 << std::hex << usDetectorId << std::dec << " added to list";
334}
335// -------------------------------------------------------------------------
336
338{
339 fulCurrentTsIdx = ts.index();
340 fdTsStartTime = static_cast<Double_t>(ts.descriptor(0, 0).idx);
341
343 if (0 == fulCurrentTsIdx) { return kTRUE; } // if( 0 == fulCurrentTsIdx )
344
346 if (-1.0 == fdTsCoreSizeInNs) {
347 fuNbCoreMsPerTs = ts.num_core_microslices();
348 fuNbOverMsPerTs = ts.num_microslices(0) - ts.num_core_microslices();
351 LOG(info) << "Timeslice parameters: each TS has " << fuNbCoreMsPerTs << " Core MS and " << fuNbOverMsPerTs
352 << " Overlap MS, for a core duration of " << fdTsCoreSizeInNs << " ns and a full duration of "
353 << fdTsFullSizeInNs << " ns";
354
358 LOG(info) << "In each TS " << fuNbMsLoop << " MS will be looped over";
359 } // if( -1.0 == fdTsCoreSizeInNs )
360
363 // LOG(info) << Form( "TS %5d Start %12f Stop %12f", fulCurrentTsIdx, fdTsStartTime, fdTsStopTimeCore );
364
366 for (fuMsIndex = 0; fuMsIndex < fuNbMsLoop; fuMsIndex++) {
368 for (UInt_t uMsCompIdx = 0; uMsCompIdx < fvMsComponentsList.size(); ++uMsCompIdx) {
369 UInt_t uMsComp = fvMsComponentsList[uMsCompIdx];
370
371 if (kFALSE == ProcessMs(ts, uMsComp, fuMsIndex)) {
372 LOG(error) << "Failed to process ts " << fulCurrentTsIdx << " MS " << fuMsIndex << " for component " << uMsComp;
373 return kFALSE;
374 } // if( kFALSE == ProcessMs( ts, uMsCompIdx, fuMsIndex ) )
375 } // for( UInt_t uMsCompIdx = 0; uMsCompIdx < fvMsComponentsList.size(); ++uMsCompIdx )
376
379 // std::sort( fvmHitsInMs.begin(), fvmHitsInMs.end() );
380
382 for (auto itHitIn = fvmHitsInMs.begin(); itHitIn < fvmHitsInMs.end(); ++itHitIn) {
383 /*
384 UInt_t uFebIdx = itHitIn->GetAsic() / fUnpackPar->GetNbAsicsPerFeb()
385 + ( itHitIn->GetDpb() * fUnpackPar->GetNbCrobsPerDpb() + itHitIn->GetCrob() )
386 * fUnpackPar->GetNbFebsPerCrob();
387*/
388 UInt_t uAsicIdx = itHitIn->GetAsic();
389 UInt_t uFebIdx = itHitIn->GetAsic() / fUnpackPar->GetNbAsicsPerFeb();
390 UInt_t uChanInMod =
391 itHitIn->GetChan() + fUnpackPar->GetNbChanPerAsic() * (itHitIn->GetAsic() % fUnpackPar->GetNbAsicsPerFeb());
395 if (0 == fviFebSide[uFebIdx])
396 uChanInMod = fUnpackPar->GetNbChanPerFeb() - uChanInMod - 1 // Invert channel order
397 + fUnpackPar->GetNbChanPerFeb(); // Offset for P (back) side
398
399 Double_t dTimeInNs = itHitIn->GetTs() * stsxyter::kdClockCycleNs - fdTimeOffsetNs;
400 if (uAsicIdx < fvdTimeOffsetNsAsics.size()) dTimeInNs -= fvdTimeOffsetNsAsics[uAsicIdx];
401 ULong64_t ulTimeInNs = static_cast<ULong64_t>(dTimeInNs);
402
403 Double_t dCalAdc = fvdFebAdcOffs[uFebIdx] + (itHitIn->GetAdc() - 1) * fvdFebAdcGain[uFebIdx];
404
405 if (0 == fviFebAddress[uFebIdx] || -1 == fviFebSide[uFebIdx]) {
406 LOG(error) << Form("Digi on disabled FEB %02u has address 0x%08x and side %d", uFebIdx, fviFebAddress[uFebIdx],
407 fviFebSide[uFebIdx]);
408 } // if( 0 == fviFebAddress[ uFebIdx ] || -1 == fviFebSide[ uFebIdx ] )
409
410
412 if (fbPulserOutput && fvbFebPulser[uFebIdx]) {
413 fPulserDigiVect.emplace_back(CbmStsDigi(fviFebAddress[uFebIdx], uChanInMod, ulTimeInNs, dCalAdc));
414 } // if (fvbFebPulser[uFebIdx])
415 else {
416 fDigiVect.emplace_back(fviFebAddress[uFebIdx], uChanInMod, ulTimeInNs, dCalAdc);
417 } // else of if (fvbFebPulser[uFebIdx])
418 } // for( auto itHitIn = fvmHitsInMs.begin(); itHitIn < fvmHitsInMs.end(); ++itHitIn )
419
421 fvmHitsInMs.clear();
422 } // for( fuMsIndex = 0; fuMsIndex < uNbMsLoop; fuMsIndex ++ )
423
425 fvmHitsInMs.clear();
426 /*
427 for( UInt_t uDpb = 0; uDpb < fuNrOfDpbs; ++uDpb )
428 {
429 fvmHitsInMs[ uDpb ].clear();
430 } // for( UInt_t uDpb = 0; uDpb < fuNrOfDpbs; ++uDpb )
431*/
432
434 std::sort(fDigiVect.begin(), fDigiVect.end(),
435 [](const CbmStsDigi& a, const CbmStsDigi& b) -> bool { return a.GetTime() < b.GetTime(); });
436
438 if (fbPulserOutput) {
439 std::sort(fPulserDigiVect.begin(), fPulserDigiVect.end(),
440 [](const CbmStsDigi& a, const CbmStsDigi& b) -> bool { return a.GetTime() < b.GetTime(); });
441 }
442
444 if (fbMonitorMode) {
445 if (kFALSE == FillHistograms()) {
446 LOG(error) << "Failed to fill histos in ts " << fulCurrentTsIdx;
447 return kFALSE;
448 } // if( kFALSE == FillHistograms() )
449
451 fhVectorCapacity->Fill(fulCurrentTsIdx, fDigiVect.capacity());
452 } // if( fbMonitorMode )
453
454 if (fuTsMaxVectorSize < fDigiVect.size()) {
456 fDigiVect.shrink_to_fit();
458 } // if( fuTsMaxVectorSize < fDigiVect.size() )
459
460 return kTRUE;
461}
462
463Bool_t CbmMcbm2018UnpackerAlgoSts::ProcessMs(const fles::Timeslice& ts, size_t uMsCompIdx, size_t uMsIdx)
464{
465 auto msDescriptor = ts.descriptor(uMsCompIdx, uMsIdx);
466 fuCurrentEquipmentId = msDescriptor.eq_id;
467 const uint8_t* msContent = reinterpret_cast<const uint8_t*>(ts.content(uMsCompIdx, uMsIdx));
468
469 uint32_t uSize = msDescriptor.size;
470 fulCurrentMsIdx = msDescriptor.idx;
471 // Double_t dMsTime = (1e-9) * static_cast<double>(fulCurrentMsIdx);
472 LOG(debug) << "Microslice: " << fulCurrentMsIdx << " from EqId " << std::hex << fuCurrentEquipmentId << std::dec
473 << " has size: " << uSize;
474
475 if (0 == fvbMaskedComponents.size()) fvbMaskedComponents.resize(ts.num_components(), kFALSE);
476
477 fuCurrDpbId = static_cast<uint32_t>(fuCurrentEquipmentId & 0xFFFF);
478 // fuCurrDpbIdx = fDpbIdIndexMap[ fuCurrDpbId ];
479
481 auto it = fDpbIdIndexMap.find(fuCurrDpbId);
482 if (it == fDpbIdIndexMap.end()) {
483 if (kFALSE == fvbMaskedComponents[uMsCompIdx]) {
484 LOG(info) << "---------------------------------------------------------------";
485 /*
486 LOG(info) << "hi hv eqid flag si sv idx/start crc size offset";
487 LOG(info) << Form( "%02x %02x %04x %04x %02x %02x %016llx %08x %08x %016llx",
488 static_cast<unsigned int>(msDescriptor.hdr_id),
489 static_cast<unsigned int>(msDescriptor.hdr_ver), msDescriptor.eq_id, msDescriptor.flags,
490 static_cast<unsigned int>(msDescriptor.sys_id),
491 static_cast<unsigned int>(msDescriptor.sys_ver), msDescriptor.idx, msDescriptor.crc,
492 msDescriptor.size, msDescriptor.offset );
493*/
494 LOG(info) << FormatMsHeaderPrintout(msDescriptor);
495 LOG(warning) << "Could not find the sDPB index for AFCK id 0x" << std::hex << fuCurrDpbId << std::dec
496 << " in timeslice " << fulCurrentTsIdx << " in microslice " << uMsIdx << " component " << uMsCompIdx
497 << "\n"
498 << "If valid this index has to be added in the STS "
499 "parameter file in the DbpIdArray field";
500 fvbMaskedComponents[uMsCompIdx] = kTRUE;
501
504 if (1 == fulCurrentTsIdx) return kTRUE;
505 } // if( kFALSE == fvbMaskedComponents[ uMsComp ] )
506 else
507 return kTRUE;
508
509 return kFALSE;
510 } // if( it == fDpbIdIndexMap.end() )
511 else
513
515
517 UInt_t uTsMsbCycleHeader = std::floor(fulCurrentMsIdx / (stsxyter::kulTsCycleNbBins * stsxyter::kdClockCycleNs));
518
520 if (kTRUE == fbBinningFw)
522
523 if (0 == uMsIdx) {
524 if (uTsMsbCycleHeader != fvuCurrentTsMsbCycle[fuCurrDpbIdx])
525 LOG(info) << " TS " << std::setw(12) << fulCurrentTsIdx << " MS " << std::setw(12) << fulCurrentMsIdx
526 << " MS Idx " << std::setw(4) << uMsIdx << " Msg Idx " << std::setw(5) << 0 << " DPB " << std::setw(2)
527 << fuCurrDpbIdx << " Old TsMsb " << std::setw(5) << fvulCurrentTsMsb[fuCurrDpbIdx] << " Old MsbCy "
528 << std::setw(5) << fvuCurrentTsMsbCycle[fuCurrDpbIdx] << " New MsbCy " << uTsMsbCycleHeader;
529 fvuCurrentTsMsbCycle[fuCurrDpbIdx] = uTsMsbCycleHeader;
531 } // if( 0 == uMsIdx )
532 else if (uTsMsbCycleHeader != fvuCurrentTsMsbCycle[fuCurrDpbIdx]) {
533 if (4194303 == fvulCurrentTsMsb[fuCurrDpbIdx]) {
534 LOG(info) << " TS " << std::setw(12) << fulCurrentTsIdx << " MS " << std::setw(12) << fulCurrentMsIdx
535 << " MS Idx " << std::setw(4) << uMsIdx << " Msg Idx " << std::setw(5) << 0 << " DPB " << std::setw(2)
536 << fuCurrDpbIdx << " Old TsMsb " << std::setw(5) << fvulCurrentTsMsb[fuCurrDpbIdx] << " Old MsbCy "
537 << std::setw(5) << fvuCurrentTsMsbCycle[fuCurrDpbIdx] << " New MsbCy " << uTsMsbCycleHeader;
538 }
539 else {
540 LOG(warning) << "TS MSB cycle from MS header does not match current cycle from data "
541 << "for TS " << std::setw(12) << fulCurrentTsIdx << " MS " << std::setw(12) << fulCurrentMsIdx
542 << " MsInTs " << std::setw(3) << uMsIdx << " ====> " << fvuCurrentTsMsbCycle[fuCurrDpbIdx]
543 << " (cnt) VS " << uTsMsbCycleHeader << " (header)";
544 } // else of if( 4194303 == fvulCurrentTsMsb[fuCurrDpbIdx] )
545 fvuCurrentTsMsbCycle[fuCurrDpbIdx] = uTsMsbCycleHeader;
546 } // else if( uTsMsbCycleHeader != fvuCurrentTsMsbCycle[ fuCurrDpbIdx ] )
547
548 // If not integer number of message in input buffer, print warning/error
549 if (0 != (uSize % sizeof(stsxyter::Message)))
550 LOG(error) << "The input microslice buffer does NOT "
551 << "contain only complete sDPB messages!";
552
553 // Compute the number of complete messages in the input microslice buffer
554 uint32_t uNbMessages = (uSize - (uSize % sizeof(stsxyter::Message))) / sizeof(stsxyter::Message);
555
556 // std::vector< uint32_t > vNbMessType( 7, 0 );
557 // std::string sMessPatt = "";
558 // bool bError = false;
559
560 // Prepare variables for the loop on contents
561 const stsxyter::Message* pMess = reinterpret_cast<const stsxyter::Message*>(msContent);
562 for (uint32_t uIdx = 0; uIdx < uNbMessages; uIdx++) {
564 stsxyter::MessType typeMess = pMess[uIdx].GetMessType();
565 fmMsgCounter[typeMess]++;
566 // fhStsMessType->Fill( static_cast< uint16_t > (typeMess) );
567 // fhStsMessTypePerDpb->Fill( fuCurrDpbIdx, static_cast< uint16_t > (typeMess) );
568
569 switch (typeMess) {
571 // ++vNbMessType[0];
572 // sMessPatt += " H ";
573 // sMessPatt += ".";
574
575 // Extract the eLink and Asic indices => Should GO IN the fill method now that obly hits are link/asic specific!
576 UShort_t usElinkIdx = pMess[uIdx].GetLinkIndex();
578 if (kTRUE == fbBinningFw) usElinkIdx = pMess[uIdx].GetLinkIndexHitBinning();
579 // fhStsMessTypePerElink->Fill( usElinkIdx, static_cast< uint16_t > (typeMess) );
580 // fhStsHitsElinkPerDpb->Fill( fuCurrDpbIdx, usElinkIdx );
581
582 UInt_t uCrobIdx = usElinkIdx / fUnpackPar->GetNbElinkPerCrob();
583 Int_t uFebIdx = fUnpackPar->ElinkIdxToFebIdx(usElinkIdx);
584
585 if (-1 == uFebIdx) {
586 LOG(warning) << "CbmMcbm2018UnpackerAlgoSts::DoUnpack => "
587 << "Wrong elink Idx! Elink raw " << Form("%d remap %d", usElinkIdx, uFebIdx);
588 continue;
589 } // if( -1 == uFebIdx )
590
591 UInt_t uAsicIdx = (fuCurrDpbIdx * fUnpackPar->GetNbCrobsPerDpb() + uCrobIdx) * fUnpackPar->GetNbAsicsPerCrob()
592 + fUnpackPar->ElinkIdxToAsicIdx(1 == fviFebType[fuCurrDpbIdx][uCrobIdx][uFebIdx], usElinkIdx);
593
594 ProcessHitInfo(pMess[uIdx], usElinkIdx, uAsicIdx, uMsIdx);
595 break;
596 } // case stsxyter::MessType::Hit :
598 // ++vNbMessType[1];
599 // sMessPatt += " T ";
600
601 // fhStsMessTypePerElink->Fill( fuCurrDpbIdx * fUnpackPar->GetNbElinkPerDpb(), static_cast< uint16_t > (typeMess) );
602
603 ProcessTsMsbInfo(pMess[uIdx], uIdx, uMsIdx);
604 break;
605 } // case stsxyter::MessType::TsMsb :
607 // ++vNbMessType[2];
608 // sMessPatt += " E ";
609
610 // fhStsMessTypePerElink->Fill( fuCurrDpbIdx * fUnpackPar->GetNbElinkPerDpb(), static_cast< uint16_t > (typeMess) );
611
612 // The first message in the TS is a special ones: EPOCH
613 ProcessEpochInfo(pMess[uIdx]);
614
615 if (0 < uIdx)
616 LOG(info) << "CbmMcbm2018UnpackerAlgoSts::DoUnpack => "
617 << "EPOCH message at unexpected position in MS: message " << uIdx << " VS message 0 expected!";
618 break;
619 } // case stsxyter::MessType::TsMsb :
621 // ++vNbMessType[3];
622 // sMessPatt += " S ";
623
624 // Extract the eLink and Asic indices => Should GO IN the fill method now that obly hits are link/asic specific!
625 UShort_t usElinkIdx = pMess[uIdx].GetStatusLink();
626
627 UInt_t uCrobIdx = usElinkIdx / fUnpackPar->GetNbElinkPerCrob();
628 Int_t uFebIdx = fUnpackPar->ElinkIdxToFebIdx(usElinkIdx);
629
630 if (-1 == uFebIdx) {
631 LOG(warning) << "CbmMcbm2018UnpackerAlgoSts::DoUnpack => "
632 << "Wrong elink Idx! Elink raw " << Form("%d remap %d", usElinkIdx, uFebIdx);
633 continue;
634 } // if( -1 == uFebIdx )
635 UInt_t uAsicIdx = (fuCurrDpbIdx * fUnpackPar->GetNbCrobsPerDpb() + uCrobIdx) * fUnpackPar->GetNbAsicsPerCrob()
636 + fUnpackPar->ElinkIdxToAsicIdx(1 == fviFebType[fuCurrDpbIdx][uCrobIdx][uFebIdx], usElinkIdx);
637
638 ProcessStatusInfo(pMess[uIdx], uAsicIdx);
639 break;
640 } // case stsxyter::MessType::Status
642 // ++vNbMessType[4];
643 // sMessPatt += " Em";
644
645 // fhStsMessTypePerElink->Fill( fuCurrDpbIdx * fUnpackPar->GetNbElinkPerDpb(), static_cast< uint16_t > (typeMess) );
646 // FillTsMsbInfo( pMess[uIdx] );
647 break;
648 } // case stsxyter::MessType::Empty :
650 // ++vNbMessType[5];
651 // sMessPatt += " En";
652 // bError = pMess[uIdx].IsMsErrorFlagOn();
653
654 // fhStsMessTypePerElink->Fill( fuCurrDpbIdx * fUnpackPar->GetNbElinkPerDpb(), static_cast< uint16_t > (typeMess) );
655 // FillTsMsbInfo( pMess[uIdx] );
656 if (pMess[uIdx].IsMsErrorFlagOn()) {
658 fhMsErrorsEvo->Fill(1e-9 * fulCurrentMsIdx, pMess[uIdx].GetMsErrorType());
659 fErrVect.push_back(
660 CbmErrorMessage(ECbmModuleId::kSts, fulCurrentMsIdx, fuCurrDpbIdx, 0x20, pMess[uIdx].GetMsErrorType()));
661 } // if( pMess[uIdx].IsMsErrorFlagOn() )
662 break;
663 } // case stsxyter::MessType::EndOfMs :
665 // ++vNbMessType[6];
666 // sMessPatt += " D ";
667
668 // fhStsMessTypePerElink->Fill( fuCurrDpbIdx * fUnpackPar->GetNbElinkPerDpb(), static_cast< uint16_t > (typeMess) );
669 break;
670 } // case stsxyter::MessType::Dummy / ReadDataAck / Ack :
671 default: {
672 LOG(fatal) << "CbmMcbm2018UnpackerAlgoSts::DoUnpack => "
673 << "Unknown message type, should never happen, stopping "
674 "here! Type found was: "
675 << static_cast<int>(typeMess);
676 }
677 } // switch( typeMess )
678 } // for (uint32_t uIdx = 0; uIdx < uNbMessages; uIdx ++)
679
680 // if( 18967040000 == fulCurrentMsIdx || 18968320000 == fulCurrentMsIdx )
681 // LOG( info ) << sMessPatt;
682
683 /*
684 if( static_cast< uint16_t >( fles::MicrosliceFlags::CrcValid ) != msDescriptor.flags )
685 {
686 LOG(info) << "STS unp "
687 << " TS " << std::setw( 12 ) << fulCurrentTsIdx
688 << " MS " << std::setw( 12 ) << fulCurrentMsIdx
689 << " MS flags 0x" << std::setw( 4 ) << std::hex << msDescriptor.flags << std::dec
690 << " Size " << std::setw( 8 ) << uSize << " bytes "
691 << " H " << std::setw( 5 ) << vNbMessType[0]
692 << " T " << std::setw( 5 ) << vNbMessType[1]
693 << " E " << std::setw( 5 ) << vNbMessType[2]
694 << " S " << std::setw( 5 ) << vNbMessType[3]
695 << " Em " << std::setw( 5 ) << vNbMessType[4]
696 << " En " << std::setw( 5 ) << vNbMessType[5]
697 << " D " << std::setw( 5 ) << vNbMessType[6]
698 << " Err " << bError
699 << " Undet. bad " << ( !bError && 400 != vNbMessType[1] );
700 } // if( MicrosliceFlags::CrcValid != msDescriptor.flags )
701*/
702 return kTRUE;
703}
704
705// -------------------------------------------------------------------------
706void CbmMcbm2018UnpackerAlgoSts::ProcessHitInfo(const stsxyter::Message& mess, const UShort_t& usElinkIdx,
707 const UInt_t& uAsicIdx, const UInt_t& /*uMsIdx*/)
708{
709 UShort_t usChan = mess.GetHitChannel();
710 UShort_t usRawAdc = mess.GetHitAdc();
711 // UShort_t usFullTs = mess.GetHitTimeFull();
712 // UShort_t usTsOver = mess.GetHitTimeOver();
713 UShort_t usRawTs = mess.GetHitTime();
714
716 if (kTRUE == fbBinningFw) usRawTs = mess.GetHitTimeBinning();
717
719 // usChan = 127 - usChan;
720
721 /*
722 fhStsChanCntRaw[ uAsicIdx ]->Fill( usChan );
723 fhStsChanAdcRaw[ uAsicIdx ]->Fill( usChan, usRawAdc );
724 fhStsChanAdcRawProf[ uAsicIdx ]->Fill( usChan, usRawAdc );
725 fhStsChanRawTs[ uAsicIdx ]->Fill( usChan, usRawTs );
726 fhStsChanMissEvt[ uAsicIdx ]->Fill( usChan, mess.IsHitMissedEvts() );
727*/
728 UInt_t uCrobIdx = usElinkIdx / fUnpackPar->GetNbElinkPerCrob();
729 UInt_t uFebIdx = uAsicIdx / fUnpackPar->GetNbAsicsPerFeb();
730 // UInt_t uAsicInFeb = uAsicIdx % fUnpackPar->GetNbAsicsPerFeb();
731 UInt_t uChanInFeb = usChan + fUnpackPar->GetNbChanPerAsic() * (uAsicIdx % fUnpackPar->GetNbAsicsPerFeb());
732
734 if (usRawTs == fvvusLastTsChan[uAsicIdx][usChan] &&
735 // usRawAdc == fvvusLastAdcChan[ uAsicIdx ][ usChan ] &&
739 return;
740 } // if SMX 2.0 DPB and same TS, ADC, TS MSB, TS MSB cycle!
741 fvvusLastTsChan[uAsicIdx][usChan] = usRawTs;
742 fvvusLastAdcChan[uAsicIdx][usChan] = usRawAdc;
745
746 /*
747 fhStsFebChanCntRaw[ uFebIdx ]->Fill( uChanInFeb );
748 fhStsFebChanAdcRaw[ uFebIdx ]->Fill( uChanInFeb, usRawAdc );
749 fhStsFebChanAdcRawProf[ uFebIdx ]->Fill( uChanInFeb, usRawAdc );
750 fhStsFebChanAdcCal[ uFebIdx ]->Fill( uChanInFeb, dCalAdc );
751 fhStsFebChanAdcCalProf[ uFebIdx ]->Fill( uChanInFeb, dCalAdc );
752 fhStsFebChanRawTs[ uFebIdx ]->Fill( usChan, usRawTs );
753 fhStsFebChanMissEvt[ uFebIdx ]->Fill( usChan, mess.IsHitMissedEvts() );
754*/
755
756 // Compute the Full time stamp
757 // Use TS w/o overlap bits as they will anyway come from the TS_MSB
758 Long64_t ulHitTime = usRawTs;
759 ulHitTime +=
760 static_cast<ULong64_t>(stsxyter::kuHitNbTsBins) * static_cast<ULong64_t>(fvulCurrentTsMsb[fuCurrDpbIdx])
761 + static_cast<ULong64_t>(stsxyter::kulTsCycleNbBins) * static_cast<ULong64_t>(fvuCurrentTsMsbCycle[fuCurrDpbIdx]);
762
764 if (kTRUE == fbBinningFw)
765 ulHitTime =
766 usRawTs
767 + static_cast<ULong64_t>(stsxyter::kuHitNbTsBinsBinning) * static_cast<ULong64_t>(fvulCurrentTsMsb[fuCurrDpbIdx])
768 + static_cast<ULong64_t>(stsxyter::kulTsCycleNbBinsBinning)
769 * static_cast<ULong64_t>(fvuCurrentTsMsbCycle[fuCurrDpbIdx]);
770
771 // Convert the Hit time in bins to Hit time in ns
772 Double_t dHitTimeNs = ulHitTime * stsxyter::kdClockCycleNs;
773 /*
774 if( 3 == uAsicIdx && 32 == usChan && 10 < usRawAdc )
775 {
776
777 LOG(info) << "; TS; " << std::setw( 12 ) << fulCurrentTsIdx
778 << "; MS; " << std::setw( 12 ) << fulCurrentMsIdx
779 << "; time; " << std::setw( 12 ) << dHitTimeNs
780 << "; time to MS; " << std::setw( 12 ) << (dHitTimeNs / fulCurrentMsIdx)
781 << "; time; " << std::setw( 12 ) << ulHitTime
782 << "; TS; " << std::setw( 4 ) << usRawTs
783 << "; MSB; " << std::setw( 12 ) << fvulCurrentTsMsb[fuCurrDpbIdx]
784 << "; Cy; " << std::setw( 4 ) << fvuCurrentTsMsbCycle[fuCurrDpbIdx]
785 << "; dt; " << std::setw( 12 ) << dHitTimeNs - fdTimeA - 20971600;
786
787 if( 0.0 < fdTimeA )
788 fhPulserVsTsAB->Fill( fusTsA, usRawTs, dHitTimeNs - fdTimeA - 20971600 );
789
790 fusTsA = usRawTs;
791 fdTimeA = dHitTimeNs;
792 } // if( 3 == uAsicIdx && 32 == usChan && 10 < usRawAdc )
793*/
795 if (0 != fviFebAddress[uFebIdx] && fdAdcCut < usRawAdc) {
798 if (kFALSE == fbUseChannelMask || false == fvvbMaskedChannels[uFebIdx][uChanInFeb])
799 fvmHitsInMs.push_back(stsxyter::FinalHit(ulHitTime, usRawAdc, uAsicIdx, usChan, fuCurrDpbIdx, uCrobIdx));
800 } // if( 0 != fviFebAddress[ uFebIdx ] )
801
803 if (mess.IsHitMissedEvts())
804 fErrVect.push_back(
805 CbmErrorMessage(ECbmModuleId::kSts, dHitTimeNs, uAsicIdx, 1 << stsxyter::kusLenStatStatus, usChan));
806
807 // Check Starting point of histos with time as X axis
808 if (-1 == fdStartTime) fdStartTime = dHitTimeNs;
809 /*
810 // Fill histos with time as X axis
811 Double_t dTimeSinceStartSec = (dHitTimeNs - fdStartTime)* 1e-9;
812 Double_t dTimeSinceStartMin = dTimeSinceStartSec / 60.0;
813
814 fviFebCountsSinceLastRateUpdate[uFebIdx]++;
815 fvdFebChanCountsSinceLastRateUpdate[uFebIdx][uChanInFeb] += 1;
816
817 fhStsFebChanHitRateEvo[ uFebIdx ]->Fill( dTimeSinceStartSec, uChanInFeb );
818 fhStsFebAsicHitRateEvo[ uFebIdx ]->Fill( dTimeSinceStartSec, uAsicInFeb );
819 fhStsFebHitRateEvo[ uFebIdx ]->Fill( dTimeSinceStartSec );
820 fhStsFebChanHitRateEvoLong[ uFebIdx ]->Fill( dTimeSinceStartMin, uChanInFeb, 1.0/60.0 );
821 fhStsFebAsicHitRateEvoLong[ uFebIdx ]->Fill( dTimeSinceStartMin, uAsicInFeb, 1.0/60.0 );
822 fhStsFebHitRateEvoLong[ uFebIdx ]->Fill( dTimeSinceStartMin, 1.0/60.0 );
823 if( mess.IsHitMissedEvts() )
824 {
825 fhStsFebChanMissEvtEvo[ uFebIdx ]->Fill( dTimeSinceStartSec, uChanInFeb );
826 fhStsFebAsicMissEvtEvo[ uFebIdx ]->Fill( dTimeSinceStartSec, uAsicInFeb );
827 fhStsFebMissEvtEvo[ uFebIdx ]->Fill( dTimeSinceStartSec );
828 } // if( mess.IsHitMissedEvts() )
829*/
830 /*
831 if( kTRUE == fbLongHistoEnable )
832 {
833 std::chrono::steady_clock::time_point tNow = std::chrono::steady_clock::now();
834 Double_t dUnixTimeInRun = std::chrono::duration_cast< std::chrono::seconds >(tNow - ftStartTimeUnix).count();
835 fhFebRateEvoLong[ uAsicIdx ]->Fill( dUnixTimeInRun , 1.0 / fuLongHistoBinSizeSec );
836 fhFebChRateEvoLong[ uAsicIdx ]->Fill( dUnixTimeInRun , usChan, 1.0 / fuLongHistoBinSizeSec );
837 } // if( kTRUE == fbLongHistoEnable )
838*/
839}
840
841void CbmMcbm2018UnpackerAlgoSts::ProcessTsMsbInfo(const stsxyter::Message& mess, UInt_t uMessIdx, UInt_t uMsIdx)
842{
843 UInt_t uVal = mess.GetTsMsbVal();
844
846 if (kTRUE == fbBinningFw) uVal = mess.GetTsMsbValBinning();
847
848 /*
849 if( (uVal != fvulCurrentTsMsb[fuCurrDpbIdx] + 1) && 0 < uVal &&
850 !( 1 == uMessIdx && usVal == fvulCurrentTsMsb[fuCurrDpbIdx] ) ) // 1st TS_MSB in MS is always a repeat of the last one in previous MS!
851 {
852 LOG(info) << "TS MSB not increasing by 1! TS " << std::setw( 12 ) << fulCurrentTsIdx
853 << " MS " << std::setw( 12 ) << fulCurrentMsIdx
854 << " MsInTs " << std::setw( 3 ) << uMsIdx
855 << " DPB " << std::setw( 2 ) << fuCurrDpbIdx
856 << " Mess " << std::setw( 5 ) << uMessIdx
857 << " Old TsMsb " << std::setw( 5 ) << fvulCurrentTsMsb[fuCurrDpbIdx]
858 << " new TsMsb " << std::setw( 5 ) << uVal
859 << " Diff " << std::setw( 5 ) << uVal - fvulCurrentTsMsb[fuCurrDpbIdx]
860 << " Old MsbCy " << std::setw( 5 ) << fvuCurrentTsMsbCycle[fuCurrDpbIdx];
861 } // if( (uVal != fvulCurrentTsMsb[fuCurrDpbIdx] + 1) && 0 < uVal )
862*/
863
864 // Update Status counters
865 if (uVal < fvulCurrentTsMsb[fuCurrDpbIdx]) {
866
867 LOG(info) << " TS " << std::setw(12) << fulCurrentTsIdx << " MS " << std::setw(12) << fulCurrentMsIdx << " MS Idx "
868 << std::setw(4) << uMsIdx << " Msg Idx " << std::setw(5) << uMessIdx << " DPB " << std::setw(2)
869 << fuCurrDpbIdx << " Old TsMsb " << std::setw(5) << fvulCurrentTsMsb[fuCurrDpbIdx] << " Old MsbCy "
870 << std::setw(5) << fvuCurrentTsMsbCycle[fuCurrDpbIdx] << " new TsMsb " << std::setw(5) << uVal;
871
873 } // if( uVal < fvulCurrentTsMsb[fuCurrDpbIdx] )
874 if (
875 uVal != fvulCurrentTsMsb[fuCurrDpbIdx] + 1 && !(0 == uVal && 4194303 == fvulCurrentTsMsb[fuCurrDpbIdx])
876 &&
877 1 != uMessIdx &&
878 !(0 == uVal && 0 == fvulCurrentTsMsb[fuCurrDpbIdx] && 2 == uMessIdx) &&
879 !(uVal == fvulCurrentTsMsb[fuCurrDpbIdx] && 2 == uMessIdx)
880 &&
881 uVal < fvulCurrentTsMsb
882 [fuCurrDpbIdx]
883 ) {
884 LOG(info) << "TS MSb Jump in "
885 << " TS " << std::setw(12) << fulCurrentTsIdx << " MS " << std::setw(12) << fulCurrentMsIdx << " MS Idx "
886 << std::setw(4) << uMsIdx << " Msg Idx " << std::setw(5) << uMessIdx << " DPB " << std::setw(2)
887 << fuCurrDpbIdx << " => Old TsMsb " << std::setw(5) << fvulCurrentTsMsb[fuCurrDpbIdx] << " new TsMsb "
888 << std::setw(5) << uVal;
889 } // if( uVal + 1 != fvulCurrentTsMsb[fuCurrDpbIdx] && 4194303 != uVal && 0 != fvulCurrentTsMsb[fuCurrDpbIdx] )
890
893 if (4194303 == uVal && 1 == uMessIdx) fvulCurrentTsMsb[fuCurrDpbIdx] = 0;
894 else
896 /*
897 if( 1 < uMessIdx )
898 {
899 fhStsDpbRawTsMsb->Fill( fuCurrDpbIdx, fvulCurrentTsMsb[fuCurrDpbIdx] );
900 fhStsDpbRawTsMsbSx->Fill( fuCurrDpbIdx, ( fvulCurrentTsMsb[fuCurrDpbIdx] & 0x1F ) );
901 fhStsDpbRawTsMsbDpb->Fill( fuCurrDpbIdx, ( fvulCurrentTsMsb[fuCurrDpbIdx] >> 5 ) );
902 } // if( 0 < uMessIdx )
903*/
904 // fhStsAsicTsMsb->Fill( fvulCurrentTsMsb[fuCurrDpbIdx], uAsicIdx );
905 /*
906 ULong64_t ulNewTsMsbTime = static_cast< ULong64_t >( stsxyter::kuHitNbTsBins )
907 * static_cast< ULong64_t >( fvulCurrentTsMsb[fuCurrDpbIdx])
908 + static_cast< ULong64_t >( stsxyter::kulTsCycleNbBins )
909 * static_cast< ULong64_t >( fvuCurrentTsMsbCycle[fuCurrDpbIdx] );
910*/
911}
912
914{
915 // UInt_t uVal = mess.GetEpochVal();
916 // UInt_t uCurrentCycle = uVal % stsxyter::kulTsCycleNbBins;
917
918 /*
919 // Update Status counters
920 if( usVal < fvulCurrentTsMsb[fuCurrDpbIdx] )
921 fvuCurrentTsMsbCycle[fuCurrDpbIdx] ++;
922 fvulCurrentTsMsb[fuCurrDpbIdx] = usVal;
923
924// fhStsAsicTsMsb->Fill( fvulCurrentTsMsb[fuCurrDpbIdx], uAsicIdx );
925*/
926}
927
929{
930 /*
931 UShort_t usStatusField = mess.GetStatusStatus();
932
933 fhPulserStatusMessType->Fill( uAsicIdx, usStatusField );
935 if( fbPrintMessages )
936 {
937 std::cout << Form("DPB %2u TS %12u MS %12u mess %5u ", fuCurrDpbIdx, fulCurrentTsIdx, fulCurrentMsIdx, uIdx );
938 mess.PrintMess( std::cout, fPrintMessCtrl );
939 } // if( fbPrintMessages )
940*/
943 Long64_t ulTime =
944 static_cast<ULong64_t>(stsxyter::kuHitNbTsBins) * static_cast<ULong64_t>(fvulCurrentTsMsb[fuCurrDpbIdx])
945 + static_cast<ULong64_t>(stsxyter::kulTsCycleNbBins) * static_cast<ULong64_t>(fvuCurrentTsMsbCycle[fuCurrDpbIdx]);
946
948 if (kTRUE == fbBinningFw)
949 ulTime =
950 static_cast<ULong64_t>(stsxyter::kuHitNbTsBinsBinning) * static_cast<ULong64_t>(fvulCurrentTsMsb[fuCurrDpbIdx])
951 + static_cast<ULong64_t>(stsxyter::kulTsCycleNbBinsBinning)
952 * static_cast<ULong64_t>(fvuCurrentTsMsbCycle[fuCurrDpbIdx]);
953
955 Double_t dTimeNs = ulTime * stsxyter::kdClockCycleNs;
956
957 fErrVect.push_back(CbmErrorMessage(ECbmModuleId::kSts, dTimeNs, uAsicIdx, mess.GetStatusStatus(), mess.GetData()));
958}
959
960// -------------------------------------------------------------------------
961
962// -------------------------------------------------------------------------
963
965{
968 new TH1I("hStsDigisTimeInRun", "Digis Nb vs Time in Run; Time in run [s]; Digis Nb []", 36000, 0, 3600);
970
971 fhVectorSize = new TH1I("fhVectorSize", "Size of the vector VS TS index; TS index; Size [bytes]", 10000, 0., 10000.);
973 new TH1I("fhVectorCapacity", "Size of the vector VS TS index; TS index; Size [bytes]", 10000, 0., 10000.);
976
977
978 fhMsCntEvo = new TH1I("fhMsCntEvo", "; MS index [s]; Counts []", 600, 0.0, 600.0);
980
981 fhMsErrorsEvo = new TH2I("fhMsErrorsEvo", "; MS index [s]; Error type []; Counts []", 600, 0.0, 600.0, 4, -0.5, 3.5);
983
984 /*
985 fhPulserVsTsAB = new TProfile2D( "fhPulserVsTsAB", "; TS A [bin]; TS B [bin]; dT B - A [ns]",
986 1024, -0.5, 1023.5,
987 1024, -0.5, 1023.5 );
988 AddHistoToVector( fhPulserVsTsAB, "" );
989*/
990 /*
992 for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
993 {
994 UInt_t uSector = fUnpackPar->GetGdpbToSectorOffset() + uGdpb;
995 std::string sFolder = Form( "sector%2u", uSector);
996
997 LOG(info) << "gDPB " << uGdpb << " is " << sFolder;
998
999 fvhHitsTimeToTriggerRaw.push_back( new TH1D(
1000 Form( "hHitsTimeToTriggerRawSect%2u", uSector ),
1001 Form( "Time to trigger for all neighboring hits in sector %2u; t - Ttrigg [ns]; Hits []", uSector ),
1002 2000, -5000, 5000 ) );
1003
1004 UInt_t uNbBinsDtSel = fdStarTriggerWinSize[ uGdpb ];
1005 Double_t dMaxDtSel = fdStarTriggerDelay[ uGdpb ] + fdStarTriggerWinSize[ uGdpb ];
1006 fvhHitsTimeToTriggerSel.push_back( new TH1D(
1007 Form( "hHitsTimeToTriggerSelSect%2u", uSector ),
1008 Form( "Time to trigger for all selected hits in sector %2u; t - Ttrigg [ns]; Hits []", uSector ),
1009 uNbBinsDtSel, fdStarTriggerDelay[ uGdpb ], dMaxDtSel ) );
1010
1012 AddHistoToVector( fvhHitsTimeToTriggerRaw[ uGdpb ], sFolder );
1013 AddHistoToVector( fvhHitsTimeToTriggerSel[ uGdpb ], sFolder );
1014
1015 if( kTRUE == fbDebugMonitorMode )
1016 {
1017 fvhHitsTimeToTriggerSelVsDaq.push_back( new TH2D(
1018 Form( "hHitsTimeToTriggerSelVsDaqSect%2u", uSector ),
1019 Form( "Time to trigger for all selected hits vs DAQ CMD in sector %2u; t - Ttrigg [ns]; DAQ CMD []; Hits []", uSector ),
1020 uNbBinsDtSel, fdStarTriggerDelay[ uGdpb ], dMaxDtSel,
1021 16, 0., 16. ) );
1022
1023 fvhHitsTimeToTriggerSelVsTrig.push_back( new TH2D(
1024 Form( "hHitsTimeToTriggerSelVsTrigSect%2u", uSector ),
1025 Form( "Time to trigger for all selected hits vs TRIG CMD in sector %2u; t - Ttrigg [ns]; TRIG CMD []; Hits []", uSector ),
1026 uNbBinsDtSel, fdStarTriggerDelay[ uGdpb ], dMaxDtSel,
1027 16, 0., 16. ) );
1028
1029 fvhTriggerDt.push_back( new TH1I(
1030 Form( "hTriggerDtSect%2u", uSector ),
1031 Form( "Trigger time difference between sector %2u and the first sector, full events only; Ttrigg%2u - TtriggRef [Clk]; events []",
1032 uSector, uSector ),
1033 200, -100, 100 ) );
1034
1037 UInt_t uNbBinsInTs = fdMsSizeInNs * 111 / 1000. / 10.;
1038 UInt_t uNbBinsInMs = fdMsSizeInNs * 20 / 1000. / 10.;
1039
1040 fvhTriggerDistributionInTs.push_back( new TH1I( Form( "hTriggerDistInTsSect%2u", uSector ),
1041 Form( "Trigger distribution inside TS in sector %2u; Time in TS [us]; Trigger [];", uSector ),
1042 uNbBinsInTs, -0.5 - fdMsSizeInNs * 10 / 1000., fdMsSizeInNs * 101 / 1000. - 0.5 ) );
1043
1044 fvhTriggerDistributionInMs.push_back( new TH1I( Form( "hTriggerDistInMsSect%2u", uSector ),
1045 Form( "Trigger distribution inside MS in sector %2u; Time in MS [us]; Trigger [];", uSector ),
1046 uNbBinsInMs, -0.5 - fdMsSizeInNs * 10 / 1000., fdMsSizeInNs * 10 / 1000. - 0.5 ) );
1047
1048 fvhMessDistributionInMs.push_back( new TH1I( Form( "hMessDistInMsSect%2u", uSector ),
1049 Form( "Messages distribution inside MS in sector %2u; Time in MS [us]; Trigger [];", uSector ),
1050 uNbBinsInMs, -0.5 - fdMsSizeInNs * 10 / 1000., fdMsSizeInNs * 10 / 1000. - 0.5 ) );
1051
1053 AddHistoToVector( fvhHitsTimeToTriggerSelVsDaq[ uGdpb ], sFolder );
1054 AddHistoToVector( fvhHitsTimeToTriggerSelVsTrig[ uGdpb ], sFolder );
1055 AddHistoToVector( fvhTriggerDt[ uGdpb ], sFolder );
1056 AddHistoToVector( fvhTriggerDistributionInTs[ uGdpb ], sFolder );
1057 AddHistoToVector( fvhTriggerDistributionInMs[ uGdpb ], sFolder );
1058 AddHistoToVector( fvhMessDistributionInMs[ uGdpb ], sFolder );
1059 } // if( kTRUE == fbDebugMonitorMode )
1060 } // for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
1061
1063 fhEventNbPerTs = new TH1I( "hEventNbPerTs",
1064 "Number of Events per TS; Events []; TS []",
1065 1000 , 0, 1000 );
1066
1067 fhEventSizeDistribution = new TH1I( "hEventSizeDistribution",
1068 "Event size distribution; Event size [byte]; Events []",
1069 CbmTofStarSubevent2019::GetMaxOutputSize()/8 , 0, CbmTofStarSubevent2019::GetMaxOutputSize() );
1070
1071 fhEventSizeEvolution = new TProfile( "hEventSizeEvolution",
1072 "Event size evolution; Time in run [min]; mean Event size [byte];",
1073 14400, 0, 14400 );
1074
1075 fhEventNbEvolution = new TH1I( "hEventNbEvolution",
1076 "Event number evolution; Time in run [min]; Events [];",
1077 14400, 0, 14400 );
1078
1080 AddHistoToVector( fhEventNbPerTs, "eventbuilder" );
1081 AddHistoToVector( fhEventSizeDistribution, "eventbuilder" );
1082 AddHistoToVector( fhEventSizeEvolution, "eventbuilder" );
1083 AddHistoToVector( fhEventNbEvolution, "eventbuilder" );
1084
1085 if( kTRUE == fbDebugMonitorMode )
1086 {
1089 UInt_t uNbBinsInTs = fdMsSizeInNs * 101 / 1000. / 10.;
1090
1091 fhEventNbDistributionInTs = new TH1I( "hEventNbDistributionInTs",
1092 "Event number distribution inside TS; Time in TS [us]; Events [];",
1093 uNbBinsInTs, -0.5, fdMsSizeInNs * 101 / 1000. - 0.5 );
1094
1095 fhEventSizeDistributionInTs = new TProfile( "hEventSizeDistributionInTs",
1096 "Event size distribution inside TS; Time in TS [us]; mean Event size [Byte];",
1097 uNbBinsInTs, -0.5, fdMsSizeInNs * 101 / 1000. - 0.5 );
1098
1099 fhRawTriggersStats = new TH2I(
1100 "hRawTriggersStats",
1101 "Raw triggers statistics per sector; ; Sector []; Messages []",
1102 5, 0, 5,
1103 12, 13, 25 );
1104 fhRawTriggersStats->GetXaxis()->SetBinLabel( 1, "A");
1105 fhRawTriggersStats->GetXaxis()->SetBinLabel( 2, "B");
1106 fhRawTriggersStats->GetXaxis()->SetBinLabel( 3, "C");
1107 fhRawTriggersStats->GetXaxis()->SetBinLabel( 4, "D");
1108 fhRawTriggersStats->GetXaxis()->SetBinLabel( 5, "F");
1109
1110 fhMissingTriggersEvolution = new TH2I(
1111 "hMissingTriggersEvolution",
1112 "Missing trigger counts per sector vs time in run; Time in run [min]; Sector []; Missing triggers []",
1113 14400, 0, 14400,
1114 12, 13, 25 );
1115
1117 AddHistoToVector( fhEventNbDistributionInTs, "eventbuilder" );
1118 AddHistoToVector( fhEventSizeDistributionInTs, "eventbuilder" );
1119 AddHistoToVector( fhRawTriggersStats, "eventbuilder" );
1120 AddHistoToVector( fhMissingTriggersEvolution, "eventbuilder" );
1121 } // if( kTRUE == fbDebugMonitorMode )
1122
1124 Double_t w = 10;
1125 Double_t h = 10;
1126
1128 fcTimeToTrigRaw = new TCanvas( "cTimeToTrigRaw", "Raw Time to trig for all sectors", w, h);
1129 fcTimeToTrigRaw->Divide( 2, fuNrOfGdpbs / 2 );
1130 for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
1131 {
1132 fcTimeToTrigRaw->cd( 1 + uGdpb );
1133 gPad->SetGridx();
1134 gPad->SetGridy();
1135 gPad->SetLogy();
1136 fvhHitsTimeToTriggerRaw[ uGdpb ]->Draw();
1137 } // for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
1138
1140 fcTimeToTrigSel = new TCanvas( "cTimeToTrigSel", "Selected Time to trig for all sectors", w, h);
1141 fcTimeToTrigSel->Divide( 2, fuNrOfGdpbs / 2 );
1142 for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
1143 {
1144 fcTimeToTrigSel->cd( 1 + uGdpb );
1145 gPad->SetGridx();
1146 gPad->SetGridy();
1147 gPad->SetLogy();
1148 fvhHitsTimeToTriggerSel[ uGdpb ]->Draw();
1149 } // for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
1150
1151 if( kTRUE == fbDebugMonitorMode )
1152 {
1154 fcTrigDistMs = new TCanvas( "cTrigDistMs", "Trigger time to MS start for all sectors", w, h);
1155 fcTrigDistMs->Divide( 2, fuNrOfGdpbs / 2 );
1156 for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
1157 {
1158 fcTrigDistMs->cd( 1 + uGdpb );
1159 gPad->SetGridx();
1160 gPad->SetGridy();
1161 gPad->SetLogy();
1162 fvhTriggerDistributionInMs[ uGdpb ]->Draw();
1163 } // for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
1164
1166 fcMessDistMs = new TCanvas( "cMessDistMs", "Message time to MS start for all sectors", w, h);
1167 fcMessDistMs->Divide( 2, fuNrOfGdpbs / 2 );
1168 for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
1169 {
1170 fcMessDistMs->cd( 1 + uGdpb );
1171 gPad->SetGridx();
1172 gPad->SetGridy();
1173 gPad->SetLogy();
1174 fvhMessDistributionInMs[ uGdpb ]->Draw();
1175 } // for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
1176 } // if( kTRUE == fbDebugMonitorMode )
1177
1179 fcEventBuildStats = new TCanvas( "cEvtBuildStats", "Event building statistics", w, h);
1180 if( kTRUE == fbDebugMonitorMode )
1181 fcEventBuildStats->Divide( 2, 3 );
1182 else fcEventBuildStats->Divide( 2, 2 );
1183
1184 fcEventBuildStats->cd( 1 );
1185 gPad->SetGridx();
1186 gPad->SetGridy();
1187 gPad->SetLogy();
1188 fhEventNbPerTs->Draw();
1189
1190 fcEventBuildStats->cd( 2 );
1191 gPad->SetGridx();
1192 gPad->SetGridy();
1193 gPad->SetLogy();
1194 fhEventSizeDistribution->Draw();
1195
1196 fcEventBuildStats->cd( 3 );
1197 gPad->SetGridx();
1198 gPad->SetGridy();
1199 gPad->SetLogy();
1200 fhEventSizeEvolution->Draw();
1201
1202 fcEventBuildStats->cd( 4 );
1203 gPad->SetGridx();
1204 gPad->SetGridy();
1205 gPad->SetLogy();
1206 fhEventNbEvolution->Draw();
1207
1208 if( kTRUE == fbDebugMonitorMode )
1209 {
1210 fcEventBuildStats->cd( 5 );
1211 gPad->SetGridx();
1212 gPad->SetGridy();
1213 gPad->SetLogy();
1214 fhEventNbDistributionInTs->Draw();
1215
1216 fcEventBuildStats->cd( 6 );
1217 gPad->SetGridx();
1218 gPad->SetGridy();
1219 gPad->SetLogy();
1220 fhEventSizeDistributionInTs->Draw();
1221 } // if( kTRUE == fbDebugMonitorMode )
1222
1223 AddCanvasToVector( fcEventBuildStats, "canvases" );
1224*/
1225 return kTRUE;
1226}
1228{
1229 for (auto itHit = fDigiVect.begin(); itHit != fDigiVect.end(); ++itHit) {
1230 fhDigisTimeInRun->Fill(itHit->GetTime() * 1e-9);
1231 } // for( auto itHit = fDigiVect.begin(); itHit != fDigiVect.end(); ++itHit)
1232 if (fbPulserOutput) {
1233 for (auto itHit = fPulserDigiVect.begin(); itHit != fPulserDigiVect.end(); ++itHit) {
1234 fhDigisTimeInRun->Fill(itHit->GetTime() * 1e-9);
1235 } // for( auto itHit = fPulserDigiVect.begin(); itHit != fPulserDigiVect.end(); ++itHit)
1236 }
1237 return kTRUE;
1238}
1240{
1241 fhDigisTimeInRun->Reset();
1242 /*
1243 for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
1244 {
1245 fvhHitsTimeToTriggerRaw[ uGdpb ]->Reset();
1246 fvhHitsTimeToTriggerSel[ uGdpb ]->Reset();
1247
1248 if( kTRUE == fbDebugMonitorMode )
1249 {
1250 fvhHitsTimeToTriggerSelVsDaq[ uGdpb ]->Reset();
1251 fvhHitsTimeToTriggerSelVsTrig[ uGdpb ]->Reset();
1252 fvhTriggerDt[ uGdpb ]->Reset();
1253 fvhTriggerDistributionInTs[ uGdpb ]->Reset();
1254 fvhTriggerDistributionInMs[ uGdpb ]->Reset();
1255 fvhMessDistributionInMs[ uGdpb ]->Reset();
1256 } // if( kTRUE == fbDebugMonitorMode )
1257 } // for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
1258
1260 fhEventNbPerTs->Reset();
1261 fhEventSizeDistribution->Reset();
1262 fhEventSizeEvolution->Reset();
1263 fhEventNbEvolution->Reset();
1264
1265 if( kTRUE == fbDebugMonitorMode )
1266 {
1267 fhEventNbDistributionInTs->Reset();
1268 fhEventSizeDistributionInTs->Reset();
1269 fhRawTriggersStats->Reset();
1270 fhMissingTriggersEvolution->Reset();
1271 } // if( kTRUE == fbDebugMonitorMode )
1272*/
1273 return kTRUE;
1274}
1275// -------------------------------------------------------------------------
1276
1277void CbmMcbm2018UnpackerAlgoSts::SetTimeOffsetNsAsic(UInt_t uAsicIdx, Double_t dOffsetIn)
1278{
1279 if (uAsicIdx >= fvdTimeOffsetNsAsics.size()) {
1280 fvdTimeOffsetNsAsics.resize(uAsicIdx + 1, 0.0);
1281 } // if( uAsicIdx < fvdTimeOffsetNsAsics.size() )
1282
1283 fvdTimeOffsetNsAsics[uAsicIdx] = dOffsetIn;
1284}
1285// -------------------------------------------------------------------------
1286void CbmMcbm2018UnpackerAlgoSts::MaskNoisyChannel(UInt_t uFeb, UInt_t uChan, Bool_t bMasked)
1287{
1288 if (kFALSE == fbUseChannelMask) {
1289 fbUseChannelMask = kTRUE;
1291 for (UInt_t uFebIdx = 0; uFebIdx < fuNbFebs; ++uFebIdx) {
1292 fvvbMaskedChannels[uFebIdx].resize(fUnpackPar->GetNbChanPerFeb(), false);
1293 } // for( UInt_t uFeb = 0; uFeb < fuNbFebs; ++uFeb )
1294 } // if( kFALSE == fbUseChannelMask )
1295
1296 if (uFeb < fuNbFebs && uChan < fUnpackPar->GetNbChanPerFeb()) fvvbMaskedChannels[uFeb][uChan] = bMasked;
1297 else
1298 LOG(fatal) << "CbmMcbm2018UnpackerAlgoSts::MaskNoisyChannel => Invalid FEB "
1299 "and/or CHAN index:"
1300 << Form(" %u vs %u and %u vs %u", uFeb, fuNbFebs, uChan, fUnpackPar->GetNbChanPerFeb());
1301}
1302// -------------------------------------------------------------------------
@ kSts
Silicon Tracking System.
Definition CbmDefs.h:48
std::string FormatMsHeaderPrintout(const fles::MicrosliceDescriptor &msDescriptor)
int Int_t
bool Bool_t
CbmRoot (+externals) headers.
void AddMsComponentToList(size_t component, UShort_t usDetectorId)
std::vector< stsxyter::FinalHit > fvmHitsInMs
Hits time-sorting.
std::vector< std::vector< UShort_t > > fvvusLastTsChan
std::vector< Bool_t > fvbMaskedComponents
Switch ON the filling of a additional set of histograms.
Bool_t fbBinningFw
=> Quick and dirty hack for binning FW!!!
std::vector< Double_t > fvdFebAdcOffs
ADC gain in e-/b, [ NbDpb * NbCrobPerDpb * NbFebsPerCrob ].
std::vector< Int_t > fviFebAddress
Pulser flag for each FEB, [ NbDpb * NbCrobPerDpb * NbFebsPerCrob ].
std::vector< CbmStsDigi > fPulserDigiVect
All hits (time in bins, ADC in bins, asic, channel) in last MS, sorted with "<" operator.
UInt_t fuNbFebs
Array to hold the active flag for all CROBs, [ NbDpb ][ NbCrobPerDpb ].
std::vector< Double_t > fvdFebAdcGain
Module side for each FEB, [ NbDpb * NbCrobPerDpb * NbFebsPerCrob ].
TH1 * fhDigisTimeInRun
TS MSB cycle of last hit message for each channel, [ AsicIdx ][ Chan ].
Int_t fiRunStartDateTimeSec
Index of the DPB from which the MS currently unpacked is coming.
std::vector< Int_t > fviModuleType
Total number of STS modules in the setup.
std::vector< UInt_t > fvuCurrentTsMsbCycle
Current TS MSB for each DPB.
void ProcessEpochInfo(const stsxyter::Message &mess)
UInt_t fuMsIndex
Start Time in ns of current MS from its index field in header.
std::vector< Int_t > fviFebSide
STS address for each FEB, [ NbDpb * NbCrobPerDpb * NbFebsPerCrob ].
std::vector< std::vector< std::vector< Int_t > > > fviFebModuleIdx
Number of StsXyter ASICs.
std::vector< Int_t > fviModAddress
Type of each module: 0 for connectors on the right, 1 for connectors on the left.
UInt_t fuNbStsXyters
Number of FEBs with StsXyter ASICs.
std::vector< std::vector< UShort_t > > fvvusLastTsMsbCycleChan
TS MSB of last hit message for each channel, [ AsicIdx ][ Chan ].
void ProcessHitInfo(const stsxyter::Message &mess, const UShort_t &usElinkIdx, const UInt_t &uAsicIdx, const UInt_t &uMsIdx)
UInt_t fuCurrDpbId
Current equipment ID, tells from which DPB the current MS is originating.
std::vector< std::vector< Bool_t > > fvbCrobActiveFlag
Map of DPB Identifier to DPB index.
Bool_t fbDebugMonitorMode
Switch ON the filling of a minimal set of histograms.
void ProcessStatusInfo(const stsxyter::Message &mess, const UInt_t &uAsicIdx)
void ProcessTsMsbInfo(const stsxyter::Message &mess, UInt_t uMessIdx=0, UInt_t uMsIdx=0)
std::vector< std::vector< UShort_t > > fvvusLastTsMsbChan
ADC of last hit message for each channel, [ AsicIdx ][ Chan ].
std::vector< std::vector< std::vector< Int_t > > > fviFebType
STS module side for each FEB, [ NbDpb ][ NbCrobPerDpb ][ NbFebsPerCrob ], 0 = P, 1 = N,...
std::vector< ULong64_t > fvulCurrentTsMsb
std::map< UInt_t, UInt_t > fDpbIdIndexMap
Total number of STS DPBs in system.
static const UInt_t kuMaxTsMsbDiffDuplicates
Duplicate hits suppression.
Double_t fdTsStopTimeCore
Time in ns of current TS from the index of the first MS first component.
std::chrono::steady_clock::time_point ftStartTimeUnix
std::vector< std::vector< bool > > fvvbMaskedChannels
UInt_t fuNrOfDpbs
STS address for the first strip of each module.
Bool_t ProcessMs(const fles::Timeslice &ts, size_t uMsCompIdx, size_t uMsIdx)
Double_t fdMsTime
End Time in ns of current TS Core from the index of the first MS first component.
CbmMcbm2018StsPar * fUnpackPar
If ON a separate output vector of digi is used for the pulser.
UInt_t fdAdcCut
Vector of channel masks, [ NbFeb ][ NbCHanInFeb ], used only if fbUseChannelMask is true.
std::vector< bool > fvbFebPulser
FEB type, [ NbDpb ][ NbCrobPerDpb ][ NbFebsPerCrob ], 0 = A, 1 = B, -1 if inactive.
std::vector< Double_t > fvdTimeOffsetNsAsics
Double_t fdTimeOffsetNs
ADC offset in e-, [ NbDpb * NbCrobPerDpb * NbFebsPerCrob ].
Int_t fiBinSizeDatePlots
Start of run time since "epoch" in s, for the plots with date as X axis.
void SetTimeOffsetNsAsic(UInt_t uAsicIdx, Double_t dOffsetIn=0.0)
UInt_t fuCurrDpbIdx
Temp holder until Current equipment ID is properly filled in MS.
std::vector< std::vector< UShort_t > > fvvusLastAdcChan
TS of last hit message for each channel, [ AsicIdx ][ Chan ].
Bool_t ProcessTs(const fles::Timeslice &ts)
std::vector< std::vector< std::vector< Int_t > > > fviFebModuleSide
Idx of the STS module for each FEB, [ NbDpb ][ NbCrobPerDpb ][ NbFebsPerCrob ], -1 if inactive.
void MaskNoisyChannel(UInt_t uFeb, UInt_t uChan, Bool_t bMasked=kTRUE)
std::map< stsxyter::MessType, UInt_t > fmMsgCounter
std::vector< CbmErrorMessage > fErrVect
void AddHistoToVector(TNamed *pointer, std::string sFolder="")
std::vector< size_t > fvMsComponentsList
std::vector< CbmStsDigi > fDigiVect
Data class for a single-channel message in the STS.
Definition CbmStsDigi.h:40
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)
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 uint32_t GetData() 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)
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 double kdClockCycleNs
MessType
Message types.
static constexpr uint16_t kusLenStatStatus