CbmRoot
Loading...
Searching...
No Matches
CbmMcbm2018UnpackerAlgoMuch.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// ----- CbmMcbm2018UnpackerAlgoMuch -----
8// ----- Created 02.02.2019 by P.-A. Loizeau -----
9// ----- -----
10// -----------------------------------------------------------------------------
11
13
15#include "CbmMcbm2018MuchPar.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 , fiFlag(0)
40 , fUnpackPar(nullptr)
41 , fuNrOfDpbs(0)
44 , fuNbFebs(0)
45 ,
46 //fviFebAddress(),
48 , fdTimeOffsetNs(0.0)
50 , fbUseChannelMask(kFALSE)
52 , fdAdcCut(0)
55 , fdTsStartTime(-1.0)
56 , fdTsStopTimeCore(-1.0)
57 , fdMsTime(-1.0)
58 , fuMsIndex(0)
59 , fmMsgCounter()
61 , fuCurrDpbId(0)
62 , fuCurrDpbIdx(0)
67 , fdStartTime(0.0)
68 , fdStartTimeMsSz(0.0)
69 , ftStartTimeUnix(std::chrono::steady_clock::now())
70 , fvmHitsInMs()
75 , fhDigisTimeInRun(nullptr)
76/*
77 fvhHitsTimeToTriggerRaw(),
78 fvhMessDistributionInMs(),
79 fhEventNbPerTs( nullptr ),
80 fcTimeToTrigRaw( nullptr )
81*/
82{
83}
85{
87 fvmHitsInMs.clear();
88 /*
89 for( UInt_t uDpb = 0; uDpb < fuNrOfDpbs; ++uDpb )
90 {
91 fvmHitsInMs[ uDpb ].clear();
92 } // for( UInt_t uDpb = 0; uDpb < fuNrOfDpbs; ++uDpb )
93*/
94
95 if (nullptr != fParCList) delete fParCList;
96 if (nullptr != fUnpackPar) delete fUnpackPar;
97}
98
99// -------------------------------------------------------------------------
101{
102 LOG(info) << "Initializing mCBM MUCH 2019 unpacker algo";
103
104 return kTRUE;
105}
113
114// -------------------------------------------------------------------------
116{
117 LOG(info) << "Init parameter containers for CbmMcbm2018UnpackerAlgoMuch";
118 Bool_t initOK = ReInitContainers();
119
120 return initOK;
121}
123{
124 LOG(info) << "**********************************************";
125 LOG(info) << "ReInit parameter containers for CbmMcbm2018UnpackerAlgoMuch";
126
127 fUnpackPar = (CbmMcbm2018MuchPar*) fParCList->FindObject("CbmMcbm2018MuchPar");
128 if (nullptr == fUnpackPar) return kFALSE;
129
130 Bool_t initOK = InitParameters();
131
132 return initOK;
133}
135{
136 if (nullptr == fParCList) fParCList = new TList();
137 fUnpackPar = new CbmMcbm2018MuchPar("CbmMcbm2018MuchPar");
138 fParCList->Add(fUnpackPar);
139
140 return fParCList;
141}
143{
144 fuNrOfDpbs = fUnpackPar->GetNrOfDpbs();
145 LOG(info) << "Nr. of MUCH DPBs: " << fuNrOfDpbs;
146
147 fDpbIdIndexMap.clear();
148 for (UInt_t uDpb = 0; uDpb < fuNrOfDpbs; ++uDpb) {
149 fDpbIdIndexMap[fUnpackPar->GetDpbId(uDpb)] = uDpb;
150 LOG(info) << "Eq. ID for DPB #" << std::setw(2) << uDpb << " = 0x" << std::setw(4) << std::hex
151 << fUnpackPar->GetDpbId(uDpb) << std::dec << " => " << fDpbIdIndexMap[fUnpackPar->GetDpbId(uDpb)];
152 } // for( UInt_t uDpb = 0; uDpb < fuNrOfDpbs; ++uDpb )
153
154 fuNbFebs = fUnpackPar->GetNrOfFebs();
155 LOG(info) << "Nr. of FEBs: " << fuNbFebs;
156
157 fuNbStsXyters = fUnpackPar->GetNrOfAsics();
158 LOG(info) << "Nr. of StsXyter ASICs: " << fuNbStsXyters;
159
160 if (fvdTimeOffsetNsAsics.size() < fuNbStsXyters) {
162 } // if( fvdTimeOffsetNsAsics.size() < fuNbStsXyters )
163
165 for (UInt_t uDpb = 0; uDpb < fuNrOfDpbs; ++uDpb) {
166 fvbCrobActiveFlag[uDpb].resize(fUnpackPar->GetNbCrobsPerDpb());
167 for (UInt_t uCrobIdx = 0; uCrobIdx < fUnpackPar->GetNbCrobsPerDpb(); ++uCrobIdx) {
168 fvbCrobActiveFlag[uDpb][uCrobIdx] = fUnpackPar->IsCrobActive(uDpb, uCrobIdx);
169 } // for( UInt_t uCrobIdx = 0; uCrobIdx < fUnpackPar->GetNbCrobsPerDpb(); ++uCrobIdx )
170 } // for( UInt_t uDpb = 0; uDpb < fuNrOfDpbs; ++uDpb )
171
172 for (UInt_t uDpb = 0; uDpb < fuNrOfDpbs; ++uDpb) {
173 TString sPrintoutLine = Form("DPB #%02u CROB Active ?: ", uDpb);
174 for (UInt_t uCrobIdx = 0; uCrobIdx < fUnpackPar->GetNbCrobsPerDpb(); ++uCrobIdx) {
175 sPrintoutLine += Form("%1u", (fvbCrobActiveFlag[uDpb][uCrobIdx] == kTRUE));
176 } // for( UInt_t uCrobIdx = 0; uCrobIdx < fUnpackPar->GetNbCrobsPerDpb(); ++uCrobIdx )
177 LOG(info) << sPrintoutLine;
178 } // for( UInt_t uDpb = 0; uDpb < fuNrOfDpbs; ++uDpb )
179 /*
180 for( UInt_t uDpb = 0; uDpb < fuNrOfDpbs; ++uDpb )
181 {
182 for( UInt_t uCrobIdx = 0; uCrobIdx < fUnpackPar->GetNbCrobsPerDpb(); ++uCrobIdx )
183 {
184 LOG(info) << Form( "DPB #%02u CROB #%u: ", uDpb, uCrobIdx);
185 for( UInt_t uFebIdx = 0; uFebIdx < fUnpackPar->GetNbFebsPerCrob(); ++ uFebIdx )
186 if( 0 <= fviFebModuleIdx[ uDpb ][ uCrobIdx ][ uFebIdx ] )
187 {
188 LOG(info) << Form( " FEB #%02u: Mod. Idx = %03d Side %c (%2d) Type %c (%2d) ADC gain %4.0f e- ADC Offs %5.0f e-",
189 uFebIdx,
190 fviFebModuleIdx[ uDpb ][ uCrobIdx ][ uFebIdx ],
191 1 == fviFebModuleSide[ uDpb ][ uCrobIdx ][ uFebIdx ] ? 'N': 'P',
192 fviFebModuleSide[ uDpb ][ uCrobIdx ][ uFebIdx ],
193 1 == fviFebType[ uDpb ][ uCrobIdx ][ uFebIdx ] ? 'B' : 'A',
194 fviFebType[ uDpb ][ uCrobIdx ][ uFebIdx ],
195 fvdFebAdcGain[ uDpb ][ uCrobIdx ][ uFebIdx ],
196 fvdFebAdcOffs[ uDpb ][ uCrobIdx ][ uFebIdx ]
197 );
198 } // for( UInt_t uFebIdx = 0; uFebIdx < fUnpackPar->GetNbFebsPerCrob(); ++ uFebIdx )
199 } // for( UInt_t uCrobIdx = 0; uCrobIdx < fUnpackPar->GetNbCrobsPerDpb(); ++uCrobIdx )
200 } // for( UInt_t uDpb = 0; uDpb < fuNrOfDpbs; ++uDpb )
201*/
202
203 if (fbBinningFw) LOG(info) << "Unpacking data in bin sorter FW mode";
204 else
205 LOG(info) << "Unpacking data in full time sorter FW mode (legacy)";
206
207 // Internal status initialization
210 for (UInt_t uDpb = 0; uDpb < fuNrOfDpbs; ++uDpb) {
211 fvulCurrentTsMsb[uDpb] = 0;
212 fvuCurrentTsMsbCycle[uDpb] = 0;
213 } // for( UInt_t uDpb = 0; uDpb < fuNrOfDpbs; ++uDpb )
214
219 for (UInt_t uAsicIdx = 0; uAsicIdx < fuNbStsXyters; ++uAsicIdx) {
220 fvvusLastTsChan[uAsicIdx].resize(fUnpackPar->GetNbChanPerAsic(), 0);
221 fvvusLastAdcChan[uAsicIdx].resize(fUnpackPar->GetNbChanPerAsic(), 0);
222 fvvusLastTsMsbChan[uAsicIdx].resize(fUnpackPar->GetNbChanPerAsic(), 0);
223 fvvusLastTsMsbCycleChan[uAsicIdx].resize(fUnpackPar->GetNbChanPerAsic(), 0);
224 } // for( UInt_t uAsicIdx = 0; uAsicIdx < fuNbStsXyters; ++ uAsicIdx )
225
226 return kTRUE;
227}
228// -------------------------------------------------------------------------
229
230void CbmMcbm2018UnpackerAlgoMuch::AddMsComponentToList(size_t component, UShort_t usDetectorId)
231{
233 for (UInt_t uCompIdx = 0; uCompIdx < fvMsComponentsList.size(); ++uCompIdx)
234 if (component == fvMsComponentsList[uCompIdx]) return;
235
237 fvMsComponentsList.push_back(component);
238
239 LOG(info) << "CbmMcbm2018UnpackerAlgoMuch::AddMsComponentToList => Component " << component << " with detector ID 0x"
240 << std::hex << usDetectorId << std::dec << " added to list";
241}
242// -------------------------------------------------------------------------
243
245{
246 fulCurrentTsIdx = ts.index();
247 fdTsStartTime = static_cast<Double_t>(ts.descriptor(0, 0).idx);
248
250 if (0 == fulCurrentTsIdx) { return kTRUE; } // if( 0 == fulCurrentTsIdx )
251
253 if (-1.0 == fdTsCoreSizeInNs) {
254 fuNbCoreMsPerTs = ts.num_core_microslices();
255 fuNbOverMsPerTs = ts.num_microslices(0) - ts.num_core_microslices();
258 LOG(info) << "Timeslice parameters: each TS has " << fuNbCoreMsPerTs << " Core MS and " << fuNbOverMsPerTs
259 << " Overlap MS, for a core duration of " << fdTsCoreSizeInNs << " ns and a full duration of "
260 << fdTsFullSizeInNs << " ns";
261
265 LOG(info) << "In each TS " << fuNbMsLoop << " MS will be looped over";
266 } // if( -1.0 == fdTsCoreSizeInNs )
267
270 // LOG(info) << Form( "TS %5d Start %12f Stop %12f", fulCurrentTsIdx, fdTsStartTime, fdTsStopTimeCore );
271
273 for (fuMsIndex = 0; fuMsIndex < fuNbMsLoop; fuMsIndex++) {
275 for (UInt_t uMsCompIdx = 0; uMsCompIdx < fvMsComponentsList.size(); ++uMsCompIdx) {
276 UInt_t uMsComp = fvMsComponentsList[uMsCompIdx];
277
278 if (kFALSE == ProcessMs(ts, uMsComp, fuMsIndex)) {
279 LOG(error) << "Failed to process ts " << fulCurrentTsIdx << " MS " << fuMsIndex << " for component " << uMsComp;
280 return kFALSE;
281 } // if( kFALSE == ProcessMs( ts, uMsCompIdx, fuMsIndex ) )
282 } // for( UInt_t uMsCompIdx = 0; uMsCompIdx < fvMsComponentsList.size(); ++uMsCompIdx )
283
286 // std::sort( fvmHitsInMs.begin(), fvmHitsInMs.end() );
287
289 for (auto itHitIn = fvmHitsInMs.begin(); itHitIn < fvmHitsInMs.end(); ++itHitIn) {
291
292 CbmMuchBeamTimeDigi* digi = CreateMuchDigi(&(*itHitIn));
293 if (digi == NULL) continue;
294 fDigiVect.push_back(*digi);
295 delete digi;
296 } // for( auto itHitIn = fvmHitsInMs.begin(); itHitIn < fvmHitsInMs.end(); ++itHitIn )
297
299 fvmHitsInMs.clear();
300 } // for( fuMsIndex = 0; fuMsIndex < uNbMsLoop; fuMsIndex ++ )
301
303 fvmHitsInMs.clear();
304 /*
305 for( UInt_t uDpb = 0; uDpb < fuNrOfDpbs; ++uDpb )
306 {
307 fvmHitsInMs[ uDpb ].clear();
308 } // for( UInt_t uDpb = 0; uDpb < fuNrOfDpbs; ++uDpb )
309*/
310
312 std::sort(fDigiVect.begin(), fDigiVect.end(), [](const CbmMuchBeamTimeDigi& a, const CbmMuchBeamTimeDigi& b) -> bool {
313 return a.GetTime() < b.GetTime();
314 });
315
317 if (fbMonitorMode) {
318 if (kFALSE == FillHistograms()) {
319 LOG(error) << "Failed to fill histos in ts " << fulCurrentTsIdx;
320 return kFALSE;
321 } // if( kFALSE == FillHistograms() )
323 fhVectorCapacity->Fill(fulCurrentTsIdx, fDigiVect.capacity());
324 } // if( fbMonitorMode )
325
326 if (fuTsMaxVectorSize < fDigiVect.size()) {
328 fDigiVect.shrink_to_fit();
330 } // if( fuTsMaxVectorSize < fDigiVect.size() )
331
332 return kTRUE;
333}
334
335// -------------------------------------------------------------------------
337{
338 /*
339 UInt_t uFebIdx = itHitIn->GetAsic() / fUnpackPar->GetNbAsicsPerFeb()
340 + ( itHitIn->GetDpb() * fUnpackPar->GetNbCrobsPerDpb() + itHitIn->GetCrob() )
341 * fUnpackPar->GetNbFebsPerCrob();
342 UInt_t uChanInFeb = itHitIn->GetChan()
343 + fUnpackPar->GetNbChanPerAsic() * (itHitIn->GetAsic() % fUnpackPar->GetNbAsicsPerFeb());
344 */
345 UInt_t uAsicIdx = itHitIn->GetAsic();
347 Int_t iFebId = fUnpackPar->GetFebId(itHitIn->GetAsic());
348
349 Short_t station = 0; // for mCBM setup only one station
350 Short_t layer = fUnpackPar->GetModule(itHitIn->GetAsic()); // 0 for Module GEM-A and 1 for Module GEM-B
351 Short_t layerside = 0; // 0 in mCBM
352 Short_t module = 0; // 0 in mCBM
353 Short_t sSector = 0; //channel values are from 0-96 therefore as per CbmMuchAddress it is sector
354 Short_t sChannel = 0; //sector values are from 0-22 therefore as per CbmMuchAddress it is channel
355
356 Int_t usChan = itHitIn->GetChan();
357
358 // ----------------------------- //
359 // Channel flip in stsXYTER v2.1 : 0<->1, 2<->3, 3<->4 and so on...
360 if (fiFlag == 1) {
361 if (usChan % 2 == 0) usChan = usChan + 1;
362 else
363 usChan = usChan - 1;
364 }
365 // ---------------------------- //
366
367 // ---------------------------- //
368 //Due to two FLEX cable connected to single FEB; First Flex Connector number 1 - 63 and second flex connector number 64 - 127
369 //in few FEB positioned these flex connectors are flipped so below correction applied.
370 if (fiFlag == 1 && layer == 0) { // First layer (GEM1) has old readout PCB
371
372 if (iFebId == 0 || iFebId == 1 || iFebId == 2 || iFebId == 3 || iFebId == 4 || iFebId == 8 || iFebId == 9
373 || iFebId == 10 || iFebId == 11 || iFebId == 17) {
374 sSector = fUnpackPar->GetPadXA(iFebId, 127 - usChan);
375 sChannel = fUnpackPar->GetPadYA(iFebId, 127 - usChan);
376 }
377 else {
378 sSector = fUnpackPar->GetPadXA(iFebId, usChan);
379 sChannel = fUnpackPar->GetPadYA(iFebId, usChan);
380 }
381 } // if( fiFlag == 1 && layer == 0 )
382
383 else if (fiFlag == 1 && layer == 1) { // second layer (GEM2) has new readout PCB
384
385 if (iFebId == 0 || iFebId == 1 || iFebId == 2 || iFebId == 3 || iFebId == 4 || iFebId == 8 || iFebId == 9
386 || iFebId == 10 || iFebId == 11 || iFebId == 17) {
387 sSector = fUnpackPar->GetPadXB(iFebId, 127 - usChan);
388 sChannel = fUnpackPar->GetPadYB(iFebId, 127 - usChan);
389 }
390 else {
391 sSector = fUnpackPar->GetPadXB(iFebId, usChan);
392 sChannel = fUnpackPar->GetPadYB(iFebId, usChan);
393 }
394 } // else if( fiFlag == 1 && layer == 1 )
395
396 else { // Both layer with same type of PCB
397 if (iFebId == 0 || iFebId == 1 || iFebId == 2 || iFebId == 3 || iFebId == 4 || iFebId == 8 || iFebId == 9
398 || iFebId == 10 || iFebId == 11 || iFebId == 17) {
399 sSector = fUnpackPar->GetPadXA(iFebId, 127 - usChan);
400 sChannel = fUnpackPar->GetPadYA(iFebId, 127 - usChan);
401 }
402 else {
403 sSector = fUnpackPar->GetPadXA(iFebId, itHitIn->GetChan());
404 sChannel = fUnpackPar->GetPadYA(iFebId, itHitIn->GetChan());
405 }
406 } // else{
407 // ---------------------------- //
408
409 // Checking for the not connected or misconnected pads
410 if (sSector < 0 || sChannel < 0) {
411 LOG(debug) << "Sector " << sSector << " channel " << sChannel
412 << " is not connected or misconnected to pad. Skipping this hit message.";
413 return NULL;
414 }
415 //Creating Unique address of the particular channel of GEM
416 UInt_t address = CbmMuchAddress::GetAddress(station, layer, layerside, module, sChannel, sSector);
417 Double_t dTimeInNs = itHitIn->GetTs() * stsxyter::kdClockCycleNs - fdTimeOffsetNs;
418 if (uAsicIdx < fvdTimeOffsetNsAsics.size()) dTimeInNs -= fvdTimeOffsetNsAsics[uAsicIdx];
419 ULong64_t ulTimeInNs = static_cast<ULong64_t>(dTimeInNs);
420
422 // fDigiVect.push_back( CbmMuchDigi( address, itHitIn->GetAdc(), ulTimeInNs ) );
423
425 CbmMuchBeamTimeDigi* digi = new CbmMuchBeamTimeDigi(address, itHitIn->GetAdc(), ulTimeInNs);
426 LOG(debug) << "Sector " << sSector << " channel " << sChannel << " layer " << layer << " Address " << address
427 << " Time " << ulTimeInNs;
428
429 digi->SetPadX(sSector);
430 digi->SetPadY(sChannel);
431 digi->SetRocId(itHitIn->GetDpb());
432 digi->SetNxId(itHitIn->GetAsic());
433 digi->SetNxCh(itHitIn->GetChan());
434 // Add Elink also in the stsxyter::FinalHit as one more variable such that may be used.
435 //digi.SetElink(itHitIn->GetElink());
436
437 return digi;
438}
439// -------------------------------------------------------------------------
440
441Bool_t CbmMcbm2018UnpackerAlgoMuch::ProcessMs(const fles::Timeslice& ts, size_t uMsCompIdx, size_t uMsIdx)
442{
443 auto msDescriptor = ts.descriptor(uMsCompIdx, uMsIdx);
444 fuCurrentEquipmentId = msDescriptor.eq_id;
445 const uint8_t* msContent = reinterpret_cast<const uint8_t*>(ts.content(uMsCompIdx, uMsIdx));
446
447 uint32_t uSize = msDescriptor.size;
448 fulCurrentMsIdx = msDescriptor.idx;
449 // Double_t dMsTime = (1e-9) * static_cast<double>(fulCurrentMsIdx);
450 LOG(debug) << "Microslice: " << fulCurrentMsIdx << " from EqId " << std::hex << fuCurrentEquipmentId << std::dec
451 << " has size: " << uSize;
452
453 if (0 == fvbMaskedComponents.size()) fvbMaskedComponents.resize(ts.num_components(), kFALSE);
454
455 fuCurrDpbId = static_cast<uint32_t>(fuCurrentEquipmentId & 0xFFFF);
456 // fuCurrDpbIdx = fDpbIdIndexMap[ fuCurrDpbId ];
457
459 auto it = fDpbIdIndexMap.find(fuCurrDpbId);
460 if (it == fDpbIdIndexMap.end()) {
461 if (kFALSE == fvbMaskedComponents[uMsCompIdx]) {
462 LOG(info) << "---------------------------------------------------------------";
463 /*
464 LOG(info) << "hi hv eqid flag si sv idx/start crc size offset";
465 LOG(info) << Form( "%02x %02x %04x %04x %02x %02x %016llx %08x %08x %016llx",
466 static_cast<unsigned int>(msDescriptor.hdr_id),
467 static_cast<unsigned int>(msDescriptor.hdr_ver), msDescriptor.eq_id, msDescriptor.flags,
468 static_cast<unsigned int>(msDescriptor.sys_id),
469 static_cast<unsigned int>(msDescriptor.sys_ver), msDescriptor.idx, msDescriptor.crc,
470 msDescriptor.size, msDescriptor.offset );
471*/
472 LOG(info) << FormatMsHeaderPrintout(msDescriptor);
473 LOG(warning) << "Could not find the sDPB index for AFCK id 0x" << std::hex << fuCurrDpbId << std::dec
474 << " in timeslice " << fulCurrentTsIdx << " in microslice " << uMsIdx << " component " << uMsCompIdx
475 << "\n"
476 << "If valid this index has to be added in the MUCH "
477 "parameter file in the DbpIdArray field";
478 fvbMaskedComponents[uMsCompIdx] = kTRUE;
479
482 if (1 == fulCurrentTsIdx) return kTRUE;
483 } // if( kFALSE == fvbMaskedComponents[ uMsComp ] )
484 else
485 return kTRUE;
486
487 return kFALSE;
488 } // if( it == fDpbIdIndexMap.end() )
489 else
491
493 UInt_t uTsMsbCycleHeader = std::floor(fulCurrentMsIdx / (stsxyter::kulTsCycleNbBins * stsxyter::kdClockCycleNs));
494
496 if (kTRUE == fbBinningFw)
498
499 if (0 == uMsIdx) {
500 if (uTsMsbCycleHeader != fvuCurrentTsMsbCycle[fuCurrDpbIdx])
501 LOG(info) << " TS " << std::setw(12) << fulCurrentTsIdx << " MS " << std::setw(12) << fulCurrentMsIdx
502 << " MS Idx " << std::setw(4) << uMsIdx << " Msg Idx " << std::setw(5) << 0 << " DPB " << std::setw(2)
503 << fuCurrDpbIdx << " Old TsMsb " << std::setw(5) << fvulCurrentTsMsb[fuCurrDpbIdx] << " Old MsbCy "
504 << std::setw(5) << fvuCurrentTsMsbCycle[fuCurrDpbIdx] << " New MsbCy " << uTsMsbCycleHeader;
505 fvuCurrentTsMsbCycle[fuCurrDpbIdx] = uTsMsbCycleHeader;
507 } // if( 0 == uMsIdx )
508 else if (uTsMsbCycleHeader != fvuCurrentTsMsbCycle[fuCurrDpbIdx]) {
509 if (4194303 == fvulCurrentTsMsb[fuCurrDpbIdx]) {
510 LOG(info) << " TS " << std::setw(12) << fulCurrentTsIdx << " MS " << std::setw(12) << fulCurrentMsIdx
511 << " MS Idx " << std::setw(4) << uMsIdx << " Msg Idx " << std::setw(5) << 0 << " DPB " << std::setw(2)
512 << fuCurrDpbIdx << " Old TsMsb " << std::setw(5) << fvulCurrentTsMsb[fuCurrDpbIdx] << " Old MsbCy "
513 << std::setw(5) << fvuCurrentTsMsbCycle[fuCurrDpbIdx] << " New MsbCy " << uTsMsbCycleHeader;
514 }
515 else {
516 LOG(warning) << "TS MSB cycle from MS header does not match current cycle from data "
517 << "for TS " << std::setw(12) << fulCurrentTsIdx << " MS " << std::setw(12) << fulCurrentMsIdx
518 << " MsInTs " << std::setw(3) << uMsIdx << " ====> " << fvuCurrentTsMsbCycle[fuCurrDpbIdx] << " VS "
519 << uTsMsbCycleHeader;
520 } // else of if( 4194303 == fvulCurrentTsMsb[fuCurrDpbIdx] )
521 fvuCurrentTsMsbCycle[fuCurrDpbIdx] = uTsMsbCycleHeader;
522 } // else if( uTsMsbCycleHeader != fvuCurrentTsMsbCycle[ fuCurrDpbIdx ] )
523
524 // If not integer number of message in input buffer, print warning/error
525 if (0 != (uSize % sizeof(stsxyter::Message)))
526 LOG(error) << "The input microslice buffer does NOT "
527 << "contain only complete sDPB messages!";
528
529 // Compute the number of complete messages in the input microslice buffer
530 uint32_t uNbMessages = (uSize - (uSize % sizeof(stsxyter::Message))) / sizeof(stsxyter::Message);
531
532 // Prepare variables for the loop on contents
533 const stsxyter::Message* pMess = reinterpret_cast<const stsxyter::Message*>(msContent);
534 for (uint32_t uIdx = 0; uIdx < uNbMessages; uIdx++) {
536 stsxyter::MessType typeMess = pMess[uIdx].GetMessType();
537 fmMsgCounter[typeMess]++;
538 // fhStsMessType->Fill( static_cast< uint16_t > (typeMess) );
539 // fhStsMessTypePerDpb->Fill( fuCurrDpbIdx, static_cast< uint16_t > (typeMess) );
540
541 switch (typeMess) {
543 // Extract the eLink and Asic indices => Should GO IN the fill method now that obly hits are link/asic specific!
544 UShort_t usElinkIdx = pMess[uIdx].GetLinkIndex();
546 if (kTRUE == fbBinningFw) usElinkIdx = pMess[uIdx].GetLinkIndexHitBinning();
547 // fhStsMessTypePerElink->Fill( usElinkIdx, static_cast< uint16_t > (typeMess) );
548 // fhStsHitsElinkPerDpb->Fill( fuCurrDpbIdx, usElinkIdx );
549
550 UInt_t uCrobIdx = usElinkIdx / fUnpackPar->GetNbElinkPerCrob();
551 Int_t uFebIdx = fUnpackPar->ElinkIdxToFebIdx(usElinkIdx);
552
553 if (-1 == uFebIdx) {
554 LOG(warning) << "CbmMcbm2018UnpackerAlgoMuch::DoUnpack => "
555 << "Wrong elink Idx! Elink raw " << Form("%d remap %d", usElinkIdx, uFebIdx);
556 continue;
557 } // if( -1 == uFebIdx )
558
559 UInt_t uAsicIdx = (fuCurrDpbIdx * fUnpackPar->GetNbCrobsPerDpb() + uCrobIdx) * fUnpackPar->GetNbAsicsPerCrob()
560 + fUnpackPar->ElinkIdxToAsicIdx(usElinkIdx);
561
562 ProcessHitInfo(pMess[uIdx], usElinkIdx, uAsicIdx, uMsIdx);
563 break;
564 } // case stsxyter::MessType::Hit :
566 // fhStsMessTypePerElink->Fill( fuCurrDpbIdx * fUnpackPar->GetNbElinkPerDpb(), static_cast< uint16_t > (typeMess) );
567
568 ProcessTsMsbInfo(pMess[uIdx], uIdx, uMsIdx);
569 break;
570 } // case stsxyter::MessType::TsMsb :
572 // fhStsMessTypePerElink->Fill( fuCurrDpbIdx * fUnpackPar->GetNbElinkPerDpb(), static_cast< uint16_t > (typeMess) );
573
574 // The first message in the TS is a special ones: EPOCH
575 ProcessEpochInfo(pMess[uIdx]);
576
577 if (0 < uIdx)
578 LOG(info) << "CbmMcbm2018UnpackerAlgoMuch::DoUnpack => "
579 << "EPOCH message at unexpected position in MS: message " << uIdx << " VS message 0 expected!";
580 break;
581 } // case stsxyter::MessType::TsMsb :
584 UShort_t usElinkIdx = pMess[uIdx].GetStatusLink();
585 UInt_t uCrobIdx = usElinkIdx / fUnpackPar->GetNbElinkPerCrob();
586 Int_t uFebIdx = fUnpackPar->ElinkIdxToFebIdx(usElinkIdx);
587
588 if (-1 == uFebIdx) {
589 LOG(warning) << "CbmMcbm2018UnpackerAlgoMuch::DoUnpack => "
590 << "Wrong elink Idx! Elink raw " << Form("%d remap %d", usElinkIdx, uFebIdx);
591 continue;
592 } // if( -1 == uFebIdx )
593
594 UInt_t uAsicIdx = (fuCurrDpbIdx * fUnpackPar->GetNbCrobsPerDpb() + uCrobIdx) * fUnpackPar->GetNbAsicsPerCrob()
595 + fUnpackPar->ElinkIdxToAsicIdx(usElinkIdx);
596
597 ProcessStatusInfo(pMess[uIdx], uAsicIdx);
598 break;
599 } // case stsxyter::MessType::Status
601 // fhStsMessTypePerElink->Fill( fuCurrDpbIdx * fUnpackPar->GetNbElinkPerDpb(), static_cast< uint16_t > (typeMess) );
602 // FillTsMsbInfo( pMess[uIdx] );
603 break;
604 } // case stsxyter::MessType::Empty :
606 if (pMess[uIdx].IsMsErrorFlagOn()) {
607 fErrVect.push_back(
608 CbmErrorMessage(ECbmModuleId::kMuch, fulCurrentMsIdx, fuCurrDpbIdx, 0x20, pMess[uIdx].GetMsErrorType()));
609 } // if( pMess[uIdx].IsMsErrorFlagOn() )
610 break;
611 } // case stsxyter::MessType::EndOfMs :
613 // fhStsMessTypePerElink->Fill( fuCurrDpbIdx * fUnpackPar->GetNbElinkPerDpb(), static_cast< uint16_t > (typeMess) );
614 break;
615 } // case stsxyter::MessType::Dummy / ReadDataAck / Ack :
616 default: {
617 LOG(fatal) << "CbmMcbm2018UnpackerAlgoMuch::DoUnpack => "
618 << "Unknown message type, should never happen, stopping "
619 "here! Type found was: "
620 << static_cast<int>(typeMess);
621 }
622 } // switch( typeMess )
623 } // for (uint32_t uIdx = 0; uIdx < uNbMessages; uIdx ++)
624
625 return kTRUE;
626}
627
628// -------------------------------------------------------------------------
629void CbmMcbm2018UnpackerAlgoMuch::ProcessHitInfo(const stsxyter::Message& mess, const UShort_t& usElinkIdx,
630 const UInt_t& uAsicIdx, const UInt_t& /*uMsIdx*/)
631{
632 UShort_t usChan = mess.GetHitChannel();
633 UShort_t usRawAdc = mess.GetHitAdc();
634 // UShort_t usFullTs = mess.GetHitTimeFull();
635 // UShort_t usTsOver = mess.GetHitTimeOver();
636 UShort_t usRawTs = mess.GetHitTime();
637
639 if (kTRUE == fbBinningFw) usRawTs = mess.GetHitTimeBinning();
640
642 // usChan = 127 - usChan;
643
644 /*
645 fhStsChanCntRaw[ uAsicIdx ]->Fill( usChan );
646 fhStsChanAdcRaw[ uAsicIdx ]->Fill( usChan, usRawAdc );
647 fhStsChanAdcRawProf[ uAsicIdx ]->Fill( usChan, usRawAdc );
648 fhStsChanRawTs[ uAsicIdx ]->Fill( usChan, usRawTs );
649 fhStsChanMissEvt[ uAsicIdx ]->Fill( usChan, mess.IsHitMissedEvts() );
650*/
651 UInt_t uCrobIdx = usElinkIdx / fUnpackPar->GetNbElinkPerCrob();
652 // UInt_t uFebIdx = uAsicIdx / fUnpackPar->GetNbAsicsPerFeb();
653 // UInt_t uAsicInFeb = uAsicIdx % fUnpackPar->GetNbAsicsPerFeb();
654 // UInt_t uChanInFeb = usChan + fUnpackPar->GetNbChanPerAsic() * (uAsicIdx % fUnpackPar->GetNbAsicsPerFeb());
655
657 if (usRawTs == fvvusLastTsChan[uAsicIdx][usChan] &&
658 // usRawAdc == fvvusLastAdcChan[ uAsicIdx ][ usChan ] &&
662 return;
663 } // if same TS, ADC, TS MSB, TS MSB cycle!
664 fvvusLastTsChan[uAsicIdx][usChan] = usRawTs;
665 fvvusLastAdcChan[uAsicIdx][usChan] = usRawAdc;
668
669 /*
670 Double_t dCalAdc = fvdFebAdcOffs[ fuCurrDpbIdx ][ uCrobIdx ][ uFebIdx ]
671 + (usRawAdc - 1)* fvdFebAdcGain[ fuCurrDpbIdx ][ uCrobIdx ][ uFebIdx ];
672*/
673 /*
674 fhStsFebChanCntRaw[ uFebIdx ]->Fill( uChanInFeb );
675 fhStsFebChanAdcRaw[ uFebIdx ]->Fill( uChanInFeb, usRawAdc );
676 fhStsFebChanAdcRawProf[ uFebIdx ]->Fill( uChanInFeb, usRawAdc );
677 fhStsFebChanAdcCal[ uFebIdx ]->Fill( uChanInFeb, dCalAdc );
678 fhStsFebChanAdcCalProf[ uFebIdx ]->Fill( uChanInFeb, dCalAdc );
679 fhStsFebChanRawTs[ uFebIdx ]->Fill( usChan, usRawTs );
680 fhStsFebChanMissEvt[ uFebIdx ]->Fill( usChan, mess.IsHitMissedEvts() );
681*/
682 // Compute the Full time stamp
683 // Use TS w/o overlap bits as they will anyway come from the TS_MSB
684 Long64_t ulHitTime = usRawTs;
685 ulHitTime +=
686 static_cast<ULong64_t>(stsxyter::kuHitNbTsBins) * static_cast<ULong64_t>(fvulCurrentTsMsb[fuCurrDpbIdx])
687 + static_cast<ULong64_t>(stsxyter::kulTsCycleNbBins) * static_cast<ULong64_t>(fvuCurrentTsMsbCycle[fuCurrDpbIdx]);
688
690 if (kTRUE == fbBinningFw)
691 ulHitTime =
692 usRawTs
693 + static_cast<ULong64_t>(stsxyter::kuHitNbTsBinsBinning) * static_cast<ULong64_t>(fvulCurrentTsMsb[fuCurrDpbIdx])
694 + static_cast<ULong64_t>(stsxyter::kulTsCycleNbBinsBinning)
695 * static_cast<ULong64_t>(fvuCurrentTsMsbCycle[fuCurrDpbIdx]);
696
697 // Convert the Hit time in bins to Hit time in ns
698 Double_t dHitTimeNs = ulHitTime * stsxyter::kdClockCycleNs;
699
701 //if( 0 != fviFebAddress[ uFebIdx ] && fdAdcCut < usRawAdc )
702 if (fdAdcCut < usRawAdc) {
705 if (kFALSE == fbUseChannelMask || false == fvvbMaskedChannels[uAsicIdx][usChan])
706 fvmHitsInMs.push_back(stsxyter::FinalHit(ulHitTime, usRawAdc, uAsicIdx, usChan, fuCurrDpbIdx, uCrobIdx));
707 } // if( 0 != fviFebAddress[ uFebIdx ] )
708
709 // fvmHitsInMs.push_back( stsxyter::FinalHit( ulHitTime, usRawAdc, uAsicIdx, usChan, fuCurrDpbIdx, uCrobIdx ) );
710
712 if (mess.IsHitMissedEvts())
713 fErrVect.push_back(
714 CbmErrorMessage(ECbmModuleId::kMuch, dHitTimeNs, uAsicIdx, 1 << stsxyter::kusLenStatStatus, usChan));
715
716 // Check Starting point of histos with time as X axis
717 if (-1 == fdStartTime) fdStartTime = dHitTimeNs;
718 /*
719 // Fill histos with time as X axis
720 Double_t dTimeSinceStartSec = (dHitTimeNs - fdStartTime)* 1e-9;
721 Double_t dTimeSinceStartMin = dTimeSinceStartSec / 60.0;
722
723 fviFebCountsSinceLastRateUpdate[uFebIdx]++;
724 fvdFebChanCountsSinceLastRateUpdate[uFebIdx][uChanInFeb] += 1;
725
726 fhStsFebChanHitRateEvo[ uFebIdx ]->Fill( dTimeSinceStartSec, uhanInFeb );
727 fhStsFebAsicHitRateEvo[ uFebIdx ]->Fill( dTimeSinceStartSec, uAsicInFeb );
728 fhStsFebHitRateEvo[ uFebIdx ]->Fill( dTimeSinceStartSec );
729 fhStsFebChanHitRateEvoLong[ uFebIdx ]->Fill( dTimeSinceStartMin, uChanInFeb, 1.0/60.0 );
730 fhStsFebAsicHitRateEvoLong[ uFebIdx ]->Fill( dTimeSinceStartMin, uAsicInFeb, 1.0/60.0 );
731 fhStsFebHitRateEvoLong[ uFebIdx ]->Fill( dTimeSinceStartMin, 1.0/60.0 );
732 if( mess.IsHitMissedEvts() )
733 {
734 fhStsFebChanMissEvtEvo[ uFebIdx ]->Fill( dTimeSinceStartSec, uChanInFeb );
735 fhStsFebAsicMissEvtEvo[ uFebIdx ]->Fill( dTimeSinceStartSec, uAsicInFeb );
736 fhStsFebMissEvtEvo[ uFebIdx ]->Fill( dTimeSinceStartSec );
737 } // if( mess.IsHitMissedEvts() )
738*/
739 /*
740 if( kTRUE == fbLongHistoEnable )
741 {
742 std::chrono::steady_clock::time_point tNow = std::chrono::steady_clock::now();
743 Double_t dUnixTimeInRun = std::chrono::duration_cast< std::chrono::seconds >(tNow - ftStartTimeUnix).count();
744 fhFebRateEvoLong[ uAsicIdx ]->Fill( dUnixTimeInRun , 1.0 / fuLongHistoBinSizeSec );
745 fhFebChRateEvoLong[ uAsicIdx ]->Fill( dUnixTimeInRun , usChan, 1.0 / fuLongHistoBinSizeSec );
746 } // if( kTRUE == fbLongHistoEnable )
747*/
748}
749
750void CbmMcbm2018UnpackerAlgoMuch::ProcessTsMsbInfo(const stsxyter::Message& mess, UInt_t uMessIdx, UInt_t uMsIdx)
751{
752 UInt_t uVal = mess.GetTsMsbVal();
753
755 if (kTRUE == fbBinningFw) uVal = mess.GetTsMsbValBinning();
756
757 /*
758 if( (uVal != fvulCurrentTsMsb[fuCurrDpbIdx] + 1) && 0 < uVal &&
759 !( 1 == uMessIdx && usVal == fvulCurrentTsMsb[fuCurrDpbIdx] ) ) // 1st TS_MSB in MS is always a repeat of the last one in previous MS!
760 {
761 LOG(info) << "TS MSB not increasing by 1! TS " << std::setw( 12 ) << fulCurrentTsIdx
762 << " MS " << std::setw( 12 ) << fulCurrentMsIdx
763 << " MsInTs " << std::setw( 3 ) << uMsIdx
764 << " DPB " << std::setw( 2 ) << fuCurrDpbIdx
765 << " Mess " << std::setw( 5 ) << uMessIdx
766 << " Old TsMsb " << std::setw( 5 ) << fvulCurrentTsMsb[fuCurrDpbIdx]
767 << " new TsMsb " << std::setw( 5 ) << uVal
768 << " Diff " << std::setw( 5 ) << uVal - fvulCurrentTsMsb[fuCurrDpbIdx]
769 << " Old MsbCy " << std::setw( 5 ) << fvuCurrentTsMsbCycle[fuCurrDpbIdx];
770 } // if( (uVal != fvulCurrentTsMsb[fuCurrDpbIdx] + 1) && 0 < uVal )
771*/
772
773 // Update Status counters
774 if (uVal < fvulCurrentTsMsb[fuCurrDpbIdx]) {
775
776 LOG(info) << " TS " << std::setw(12) << fulCurrentTsIdx << " MS " << std::setw(12) << fulCurrentMsIdx << " MS Idx "
777 << std::setw(4) << uMsIdx << " Msg Idx " << std::setw(5) << uMessIdx << " DPB " << std::setw(2)
778 << fuCurrDpbIdx << " Old TsMsb " << std::setw(5) << fvulCurrentTsMsb[fuCurrDpbIdx] << " Old MsbCy "
779 << std::setw(5) << fvuCurrentTsMsbCycle[fuCurrDpbIdx] << " new TsMsb " << std::setw(5) << uVal;
780
782 } // if( uVal < fvulCurrentTsMsb[fuCurrDpbIdx] )
783 if (
784 uVal != fvulCurrentTsMsb[fuCurrDpbIdx] + 1 && !(0 == uVal && 4194303 == fvulCurrentTsMsb[fuCurrDpbIdx])
785 &&
786 1 != uMessIdx &&
787 !(0 == uVal && 0 == fvulCurrentTsMsb[fuCurrDpbIdx] && 2 == uMessIdx) &&
788 uVal < fvulCurrentTsMsb
789 [fuCurrDpbIdx]
790 ) {
791 LOG(info) << "TS MSb Jump in "
792 << " TS " << std::setw(12) << fulCurrentTsIdx << " MS " << std::setw(12) << fulCurrentMsIdx << " MS Idx "
793 << std::setw(4) << uMsIdx << " Msg Idx " << std::setw(5) << uMessIdx << " DPB " << std::setw(2)
794 << fuCurrDpbIdx << " => Old TsMsb " << std::setw(5) << fvulCurrentTsMsb[fuCurrDpbIdx] << " new TsMsb "
795 << std::setw(5) << uVal;
796 } // if( uVal + 1 != fvulCurrentTsMsb[fuCurrDpbIdx] && 4194303 != uVal && 0 != fvulCurrentTsMsb[fuCurrDpbIdx] )
797
800 if (4194303 == uVal && 1 == uMessIdx) fvulCurrentTsMsb[fuCurrDpbIdx] = 0;
801 else
803 /*
804 if( 1 < uMessIdx )
805 {
806 fhStsDpbRawTsMsb->Fill( fuCurrDpbIdx, fvulCurrentTsMsb[fuCurrDpbIdx] );
807 fhStsDpbRawTsMsbSx->Fill( fuCurrDpbIdx, ( fvulCurrentTsMsb[fuCurrDpbIdx] & 0x1F ) );
808 fhStsDpbRawTsMsbDpb->Fill( fuCurrDpbIdx, ( fvulCurrentTsMsb[fuCurrDpbIdx] >> 5 ) );
809 } // if( 0 < uMessIdx )
810*/
811 // fhStsAsicTsMsb->Fill( fvulCurrentTsMsb[fuCurrDpbIdx], uAsicIdx );
812 /*
813 ULong64_t ulNewTsMsbTime = static_cast< ULong64_t >( stsxyter::kuHitNbTsBins )
814 * static_cast< ULong64_t >( fvulCurrentTsMsb[fuCurrDpbIdx])
815 + static_cast< ULong64_t >( stsxyter::kulTsCycleNbBins )
816 * static_cast< ULong64_t >( fvuCurrentTsMsbCycle[fuCurrDpbIdx] );
817*/
818}
819
821{
822 // UInt_t uVal = mess.GetEpochVal();
823 // UInt_t uCurrentCycle = uVal % stsxyter::kulTsCycleNbBins;
824
825 /*
826 // Update Status counters
827 if( usVal < fvulCurrentTsMsb[fuCurrDpbIdx] )
828 fvuCurrentTsMsbCycle[fuCurrDpbIdx] ++;
829 fvulCurrentTsMsb[fuCurrDpbIdx] = usVal;
830
831// fhStsAsicTsMsb->Fill( fvulCurrentTsMsb[fuCurrDpbIdx], uAsicIdx );
832*/
833}
834
836{
837 /*
838 UShort_t usStatusField = mess.GetStatusStatus();
839
840 fhPulserStatusMessType->Fill( uAsicIdx, usStatusField );
842 if( fbPrintMessages )
843 {
844 std::cout << Form("DPB %2u TS %12u MS %12u mess %5u ", fuCurrDpbIdx, fulCurrentTsIdx, fulCurrentMsIdx, uIdx );
845 mess.PrintMess( std::cout, fPrintMessCtrl );
846 } // if( fbPrintMessages )
847*/
850 Long64_t ulTime =
851 static_cast<ULong64_t>(stsxyter::kuHitNbTsBins) * static_cast<ULong64_t>(fvulCurrentTsMsb[fuCurrDpbIdx])
852 + static_cast<ULong64_t>(stsxyter::kulTsCycleNbBins) * static_cast<ULong64_t>(fvuCurrentTsMsbCycle[fuCurrDpbIdx]);
853
855 if (kTRUE == fbBinningFw)
856 ulTime =
857 static_cast<ULong64_t>(stsxyter::kuHitNbTsBinsBinning) * static_cast<ULong64_t>(fvulCurrentTsMsb[fuCurrDpbIdx])
858 + static_cast<ULong64_t>(stsxyter::kulTsCycleNbBinsBinning)
859 * static_cast<ULong64_t>(fvuCurrentTsMsbCycle[fuCurrDpbIdx]);
860
862 Double_t dTimeNs = ulTime * stsxyter::kdClockCycleNs;
863
864 fErrVect.push_back(CbmErrorMessage(ECbmModuleId::kMuch, dTimeNs, uAsicIdx, mess.GetStatusStatus(), mess.GetData()));
865}
866
867// -------------------------------------------------------------------------
868
869// -------------------------------------------------------------------------
870
872{
875 new TH1I("hMuchDigisTimeInRun", "Digis Nb vs Time in Run; Time in run [s]; Digis Nb []", 36000, 0, 3600);
877
878 fhVectorSize = new TH1I("fhVectorSize", "Size of the vector VS TS index; TS index; Size [bytes]", 10000, 0., 10000.);
880 new TH1I("fhVectorCapacity", "Size of the vector VS TS index; TS index; Size [bytes]", 10000, 0., 10000.);
883
884 /*
886 for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
887 {
888 UInt_t uSector = fUnpackPar->GetGdpbToSectorOffset() + uGdpb;
889 std::string sFolder = Form( "sector%2u", uSector);
890
891 LOG(info) << "gDPB " << uGdpb << " is " << sFolder;
892
893 fvhHitsTimeToTriggerRaw.push_back( new TH1D(
894 Form( "hHitsTimeToTriggerRawSect%2u", uSector ),
895 Form( "Time to trigger for all neighboring hits in sector %2u; t - Ttrigg [ns]; Hits []", uSector ),
896 2000, -5000, 5000 ) );
897
898 UInt_t uNbBinsDtSel = fdStarTriggerWinSize[ uGdpb ];
899 Double_t dMaxDtSel = fdStarTriggerDelay[ uGdpb ] + fdStarTriggerWinSize[ uGdpb ];
900 fvhHitsTimeToTriggerSel.push_back( new TH1D(
901 Form( "hHitsTimeToTriggerSelSect%2u", uSector ),
902 Form( "Time to trigger for all selected hits in sector %2u; t - Ttrigg [ns]; Hits []", uSector ),
903 uNbBinsDtSel, fdStarTriggerDelay[ uGdpb ], dMaxDtSel ) );
904
906 AddHistoToVector( fvhHitsTimeToTriggerRaw[ uGdpb ], sFolder );
907 AddHistoToVector( fvhHitsTimeToTriggerSel[ uGdpb ], sFolder );
908
909 if( kTRUE == fbDebugMonitorMode )
910 {
911 fvhHitsTimeToTriggerSelVsDaq.push_back( new TH2D(
912 Form( "hHitsTimeToTriggerSelVsDaqSect%2u", uSector ),
913 Form( "Time to trigger for all selected hits vs DAQ CMD in sector %2u; t - Ttrigg [ns]; DAQ CMD []; Hits []", uSector ),
914 uNbBinsDtSel, fdStarTriggerDelay[ uGdpb ], dMaxDtSel,
915 16, 0., 16. ) );
916
917 fvhHitsTimeToTriggerSelVsTrig.push_back( new TH2D(
918 Form( "hHitsTimeToTriggerSelVsTrigSect%2u", uSector ),
919 Form( "Time to trigger for all selected hits vs TRIG CMD in sector %2u; t - Ttrigg [ns]; TRIG CMD []; Hits []", uSector ),
920 uNbBinsDtSel, fdStarTriggerDelay[ uGdpb ], dMaxDtSel,
921 16, 0., 16. ) );
922
923 fvhTriggerDt.push_back( new TH1I(
924 Form( "hTriggerDtSect%2u", uSector ),
925 Form( "Trigger time difference between sector %2u and the first sector, full events only; Ttrigg%2u - TtriggRef [Clk]; events []",
926 uSector, uSector ),
927 200, -100, 100 ) );
928
931 UInt_t uNbBinsInTs = fdMsSizeInNs * 111 / 1000. / 10.;
932 UInt_t uNbBinsInMs = fdMsSizeInNs * 20 / 1000. / 10.;
933
934 fvhTriggerDistributionInTs.push_back( new TH1I( Form( "hTriggerDistInTsSect%2u", uSector ),
935 Form( "Trigger distribution inside TS in sector %2u; Time in TS [us]; Trigger [];", uSector ),
936 uNbBinsInTs, -0.5 - fdMsSizeInNs * 10 / 1000., fdMsSizeInNs * 101 / 1000. - 0.5 ) );
937
938 fvhTriggerDistributionInMs.push_back( new TH1I( Form( "hTriggerDistInMsSect%2u", uSector ),
939 Form( "Trigger distribution inside MS in sector %2u; Time in MS [us]; Trigger [];", uSector ),
940 uNbBinsInMs, -0.5 - fdMsSizeInNs * 10 / 1000., fdMsSizeInNs * 10 / 1000. - 0.5 ) );
941
942 fvhMessDistributionInMs.push_back( new TH1I( Form( "hMessDistInMsSect%2u", uSector ),
943 Form( "Messages distribution inside MS in sector %2u; Time in MS [us]; Trigger [];", uSector ),
944 uNbBinsInMs, -0.5 - fdMsSizeInNs * 10 / 1000., fdMsSizeInNs * 10 / 1000. - 0.5 ) );
945
947 AddHistoToVector( fvhHitsTimeToTriggerSelVsDaq[ uGdpb ], sFolder );
948 AddHistoToVector( fvhHitsTimeToTriggerSelVsTrig[ uGdpb ], sFolder );
949 AddHistoToVector( fvhTriggerDt[ uGdpb ], sFolder );
950 AddHistoToVector( fvhTriggerDistributionInTs[ uGdpb ], sFolder );
951 AddHistoToVector( fvhTriggerDistributionInMs[ uGdpb ], sFolder );
952 AddHistoToVector( fvhMessDistributionInMs[ uGdpb ], sFolder );
953 } // if( kTRUE == fbDebugMonitorMode )
954 } // for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
955
957 fhEventNbPerTs = new TH1I( "hEventNbPerTs",
958 "Number of Events per TS; Events []; TS []",
959 1000 , 0, 1000 );
960
961 fhEventSizeDistribution = new TH1I( "hEventSizeDistribution",
962 "Event size distribution; Event size [byte]; Events []",
963 CbmTofStarSubevent2019::GetMaxOutputSize()/8 , 0, CbmTofStarSubevent2019::GetMaxOutputSize() );
964
965 fhEventSizeEvolution = new TProfile( "hEventSizeEvolution",
966 "Event size evolution; Time in run [min]; mean Event size [byte];",
967 14400, 0, 14400 );
968
969 fhEventNbEvolution = new TH1I( "hEventNbEvolution",
970 "Event number evolution; Time in run [min]; Events [];",
971 14400, 0, 14400 );
972
974 AddHistoToVector( fhEventNbPerTs, "eventbuilder" );
975 AddHistoToVector( fhEventSizeDistribution, "eventbuilder" );
976 AddHistoToVector( fhEventSizeEvolution, "eventbuilder" );
977 AddHistoToVector( fhEventNbEvolution, "eventbuilder" );
978
979 if( kTRUE == fbDebugMonitorMode )
980 {
983 UInt_t uNbBinsInTs = fdMsSizeInNs * 101 / 1000. / 10.;
984
985 fhEventNbDistributionInTs = new TH1I( "hEventNbDistributionInTs",
986 "Event number distribution inside TS; Time in TS [us]; Events [];",
987 uNbBinsInTs, -0.5, fdMsSizeInNs * 101 / 1000. - 0.5 );
988
989 fhEventSizeDistributionInTs = new TProfile( "hEventSizeDistributionInTs",
990 "Event size distribution inside TS; Time in TS [us]; mean Event size [Byte];",
991 uNbBinsInTs, -0.5, fdMsSizeInNs * 101 / 1000. - 0.5 );
992
993 fhRawTriggersStats = new TH2I(
994 "hRawTriggersStats",
995 "Raw triggers statistics per sector; ; Sector []; Messages []",
996 5, 0, 5,
997 12, 13, 25 );
998 fhRawTriggersStats->GetXaxis()->SetBinLabel( 1, "A");
999 fhRawTriggersStats->GetXaxis()->SetBinLabel( 2, "B");
1000 fhRawTriggersStats->GetXaxis()->SetBinLabel( 3, "C");
1001 fhRawTriggersStats->GetXaxis()->SetBinLabel( 4, "D");
1002 fhRawTriggersStats->GetXaxis()->SetBinLabel( 5, "F");
1003
1004 fhMissingTriggersEvolution = new TH2I(
1005 "hMissingTriggersEvolution",
1006 "Missing trigger counts per sector vs time in run; Time in run [min]; Sector []; Missing triggers []",
1007 14400, 0, 14400,
1008 12, 13, 25 );
1009
1011 AddHistoToVector( fhEventNbDistributionInTs, "eventbuilder" );
1012 AddHistoToVector( fhEventSizeDistributionInTs, "eventbuilder" );
1013 AddHistoToVector( fhRawTriggersStats, "eventbuilder" );
1014 AddHistoToVector( fhMissingTriggersEvolution, "eventbuilder" );
1015 } // if( kTRUE == fbDebugMonitorMode )
1016
1018 Double_t w = 10;
1019 Double_t h = 10;
1020
1022 fcTimeToTrigRaw = new TCanvas( "cTimeToTrigRaw", "Raw Time to trig for all sectors", w, h);
1023 fcTimeToTrigRaw->Divide( 2, fuNrOfGdpbs / 2 );
1024 for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
1025 {
1026 fcTimeToTrigRaw->cd( 1 + uGdpb );
1027 gPad->SetGridx();
1028 gPad->SetGridy();
1029 gPad->SetLogy();
1030 fvhHitsTimeToTriggerRaw[ uGdpb ]->Draw();
1031 } // for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
1032
1034 fcTimeToTrigSel = new TCanvas( "cTimeToTrigSel", "Selected Time to trig for all sectors", w, h);
1035 fcTimeToTrigSel->Divide( 2, fuNrOfGdpbs / 2 );
1036 for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
1037 {
1038 fcTimeToTrigSel->cd( 1 + uGdpb );
1039 gPad->SetGridx();
1040 gPad->SetGridy();
1041 gPad->SetLogy();
1042 fvhHitsTimeToTriggerSel[ uGdpb ]->Draw();
1043 } // for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
1044
1045 if( kTRUE == fbDebugMonitorMode )
1046 {
1048 fcTrigDistMs = new TCanvas( "cTrigDistMs", "Trigger time to MS start for all sectors", w, h);
1049 fcTrigDistMs->Divide( 2, fuNrOfGdpbs / 2 );
1050 for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
1051 {
1052 fcTrigDistMs->cd( 1 + uGdpb );
1053 gPad->SetGridx();
1054 gPad->SetGridy();
1055 gPad->SetLogy();
1056 fvhTriggerDistributionInMs[ uGdpb ]->Draw();
1057 } // for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
1058
1060 fcMessDistMs = new TCanvas( "cMessDistMs", "Message time to MS start for all sectors", w, h);
1061 fcMessDistMs->Divide( 2, fuNrOfGdpbs / 2 );
1062 for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
1063 {
1064 fcMessDistMs->cd( 1 + uGdpb );
1065 gPad->SetGridx();
1066 gPad->SetGridy();
1067 gPad->SetLogy();
1068 fvhMessDistributionInMs[ uGdpb ]->Draw();
1069 } // for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
1070 } // if( kTRUE == fbDebugMonitorMode )
1071
1073 fcEventBuildStats = new TCanvas( "cEvtBuildStats", "Event building statistics", w, h);
1074 if( kTRUE == fbDebugMonitorMode )
1075 fcEventBuildStats->Divide( 2, 3 );
1076 else fcEventBuildStats->Divide( 2, 2 );
1077
1078 fcEventBuildStats->cd( 1 );
1079 gPad->SetGridx();
1080 gPad->SetGridy();
1081 gPad->SetLogy();
1082 fhEventNbPerTs->Draw();
1083
1084 fcEventBuildStats->cd( 2 );
1085 gPad->SetGridx();
1086 gPad->SetGridy();
1087 gPad->SetLogy();
1088 fhEventSizeDistribution->Draw();
1089
1090 fcEventBuildStats->cd( 3 );
1091 gPad->SetGridx();
1092 gPad->SetGridy();
1093 gPad->SetLogy();
1094 fhEventSizeEvolution->Draw();
1095
1096 fcEventBuildStats->cd( 4 );
1097 gPad->SetGridx();
1098 gPad->SetGridy();
1099 gPad->SetLogy();
1100 fhEventNbEvolution->Draw();
1101
1102 if( kTRUE == fbDebugMonitorMode )
1103 {
1104 fcEventBuildStats->cd( 5 );
1105 gPad->SetGridx();
1106 gPad->SetGridy();
1107 gPad->SetLogy();
1108 fhEventNbDistributionInTs->Draw();
1109
1110 fcEventBuildStats->cd( 6 );
1111 gPad->SetGridx();
1112 gPad->SetGridy();
1113 gPad->SetLogy();
1114 fhEventSizeDistributionInTs->Draw();
1115 } // if( kTRUE == fbDebugMonitorMode )
1116
1117 AddCanvasToVector( fcEventBuildStats, "canvases" );
1118*/
1119 return kTRUE;
1120}
1122{
1123 for (auto itHit = fDigiVect.begin(); itHit != fDigiVect.end(); ++itHit) {
1124 fhDigisTimeInRun->Fill(itHit->GetTime() * 1e-9);
1125 } // for( auto itHit = fDigiVect.begin(); itHit != fDigiVect.end(); ++itHit)
1126 return kTRUE;
1127}
1129{
1130 fhDigisTimeInRun->Reset();
1131 /*
1132 for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
1133 {
1134 fvhHitsTimeToTriggerRaw[ uGdpb ]->Reset();
1135 fvhHitsTimeToTriggerSel[ uGdpb ]->Reset();
1136
1137 if( kTRUE == fbDebugMonitorMode )
1138 {
1139 fvhHitsTimeToTriggerSelVsDaq[ uGdpb ]->Reset();
1140 fvhHitsTimeToTriggerSelVsTrig[ uGdpb ]->Reset();
1141 fvhTriggerDt[ uGdpb ]->Reset();
1142 fvhTriggerDistributionInTs[ uGdpb ]->Reset();
1143 fvhTriggerDistributionInMs[ uGdpb ]->Reset();
1144 fvhMessDistributionInMs[ uGdpb ]->Reset();
1145 } // if( kTRUE == fbDebugMonitorMode )
1146 } // for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
1147
1149 fhEventNbPerTs->Reset();
1150 fhEventSizeDistribution->Reset();
1151 fhEventSizeEvolution->Reset();
1152 fhEventNbEvolution->Reset();
1153
1154 if( kTRUE == fbDebugMonitorMode )
1155 {
1156 fhEventNbDistributionInTs->Reset();
1157 fhEventSizeDistributionInTs->Reset();
1158 fhRawTriggersStats->Reset();
1159 fhMissingTriggersEvolution->Reset();
1160 } // if( kTRUE == fbDebugMonitorMode )
1161*/
1162 return kTRUE;
1163}
1164// -------------------------------------------------------------------------
1165
1166void CbmMcbm2018UnpackerAlgoMuch::SetTimeOffsetNsAsic(UInt_t uAsicIdx, Double_t dOffsetIn)
1167{
1168 if (uAsicIdx >= fvdTimeOffsetNsAsics.size()) {
1169 fvdTimeOffsetNsAsics.resize(uAsicIdx + 1, 0.0);
1170 } // if( uAsicIdx < fvdTimeOffsetNsAsics.size() )
1171
1172 fvdTimeOffsetNsAsics[uAsicIdx] = dOffsetIn;
1173}
1174// -------------------------------------------------------------------------
1175void CbmMcbm2018UnpackerAlgoMuch::MaskNoisyChannel(UInt_t uFeb, UInt_t uChan, Bool_t bMasked)
1176{
1177 if (kFALSE == fbUseChannelMask) {
1178 fbUseChannelMask = kTRUE;
1180 for (UInt_t uFebIdx = 0; uFebIdx < fuNbFebs; ++uFebIdx) {
1181 fvvbMaskedChannels[uFebIdx].resize(fUnpackPar->GetNbChanPerFeb(), false);
1182 } // for( UInt_t uFeb = 0; uFeb < fuNbFebs; ++uFeb )
1183 } // if( kFALSE == fbUseChannelMask )
1184 if (uFeb < fuNbFebs && uChan < fUnpackPar->GetNbChanPerFeb()) fvvbMaskedChannels[uFeb][uChan] = bMasked;
1185 else
1186 LOG(fatal) << "CbmMcbm2018UnpackerAlgoMuch::MaskNoisyChannel => Invalid "
1187 "FEB and/or CHAN index:"
1188 << Form(" %u vs %u and %u vs %u", uFeb, fuNbFebs, uChan, fUnpackPar->GetNbChanPerFeb());
1189}
1190// -------------------------------------------------------------------------
@ kMuch
Muon detection system.
Definition CbmDefs.h:50
std::string FormatMsHeaderPrintout(const fles::MicrosliceDescriptor &msDescriptor)
int Int_t
bool Bool_t
CbmRoot (+externals) headers.
std::map< UInt_t, UInt_t > fDpbIdIndexMap
Total number of STS DPBs in system.
std::vector< std::vector< UShort_t > > fvvusLastTsChan
Limit how many different TS_MSB are checked for same duplicate hit => set to 1 uS.
std::chrono::steady_clock::time_point ftStartTimeUnix
CbmMuchBeamTimeDigi * CreateMuchDigi(stsxyter::FinalHit *)
void ProcessTsMsbInfo(const stsxyter::Message &mess, UInt_t uMessIdx=0, UInt_t uMsIdx=0)
CbmMcbm2018MuchPar * fUnpackPar
Switch to smx2.0/smx2.1 data-> fiFlag = 0 for 2.0 and fiFlag = 1 for 2.1.
Int_t fiBinSizeDatePlots
Start of run time since "epoch" in s, for the plots with date as X axis.
Bool_t fbBinningFw
=> Quick and dirty hack for binning FW!!!
void MaskNoisyChannel(UInt_t uFeb, UInt_t uChan, Bool_t bMasked=kTRUE)
Bool_t ProcessMs(const fles::Timeslice &ts, size_t uMsCompIdx, size_t uMsIdx)
std::vector< Bool_t > fvbMaskedComponents
Switch ON the filling of a additional set of histograms.
UInt_t fuNbFebs
Array to hold the active flag for all CROBs, [ NbDpb ][ NbCrobPerDpb ].
std::vector< std::vector< UShort_t > > fvvusLastTsMsbChan
ADC of last hit message for each channel, [ AsicIdx ][ Chan ].
Double_t fdTimeOffsetNs
Number of StsXyter ASICs.
std::vector< std::vector< UShort_t > > fvvusLastTsMsbCycleChan
TS MSB of last hit message for each channel, [ AsicIdx ][ Chan ].
UInt_t fuCurrDpbId
Current equipment ID, tells from which DPB the current MS is originating.
void SetTimeOffsetNsAsic(UInt_t uAsicIdx, Double_t dOffsetIn=0.0)
std::vector< UInt_t > fvuCurrentTsMsbCycle
Current TS MSB for each DPB.
UInt_t fuNbStsXyters
Number of FEBs with StsXyter ASICs.
std::vector< std::vector< Bool_t > > fvbCrobActiveFlag
Map of DPB Identifier to DPB index.
std::vector< std::vector< bool > > fvvbMaskedChannels
Bool_t fbDebugMonitorMode
Switch ON the filling of a minimal set of histograms.
Bool_t ProcessTs(const fles::Timeslice &ts)
std::vector< Double_t > fvdTimeOffsetNsAsics
std::vector< stsxyter::FinalHit > fvmHitsInMs
Hits time-sorting.
std::map< stsxyter::MessType, UInt_t > fmMsgCounter
Int_t fiRunStartDateTimeSec
Index of the DPB from which the MS currently unpacked is coming.
void AddMsComponentToList(size_t component, UShort_t usDetectorId)
void ProcessStatusInfo(const stsxyter::Message &mess, const UInt_t &uAsicIdx)
UInt_t fdAdcCut
Vector of channel masks, [ NbFeb ][ NbCHanInFeb ], used only if fbUseChannelMask is true.
TH1 * fhDigisTimeInRun
TS MSB cycle of last hit message for each channel, [ AsicIdx ][ Chan ].
std::vector< std::vector< UShort_t > > fvvusLastAdcChan
TS 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)
void ProcessEpochInfo(const stsxyter::Message &mess)
static const UInt_t kuMaxTsMsbDiffDuplicates
All hits (time in bins, ADC in bins, asic, channel) in last MS, sorted with "<" operator.
UInt_t fuMsIndex
Start Time in ns of current MS from its index field in header.
UInt_t fuCurrDpbIdx
Temp holder until Current equipment ID is properly filled in MS.
Double_t fdMsTime
End Time in ns of current TS Core from the index of the first MS first component.
Double_t fdTsStopTimeCore
Time in ns of current TS from the index of the first MS first component.
std::vector< ULong64_t > fvulCurrentTsMsb
static uint32_t GetAddress(int32_t station=0, int32_t layer=0, int32_t side=0, int32_t module=0, int32_t sector=0, int32_t channel=0)
void SetNxId(int32_t nxId)
void SetPadX(int32_t padX)
void SetPadY(int32_t padY)
void SetRocId(int32_t rocId)
void SetNxCh(int32_t nxCh)
std::vector< CbmErrorMessage > fErrVect
void AddHistoToVector(TNamed *pointer, std::string sFolder="")
std::vector< CbmMuchBeamTimeDigi > fDigiVect
uint16_t GetDpb() const
uint16_t GetAsic() const
uint16_t GetChan() const
uint64_t GetTs() const
uint16_t GetAdc() const
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