CbmRoot
Loading...
Searching...
No Matches
CbmMcbm2018MonitorAlgoPsd.cxx
Go to the documentation of this file.
1/* Copyright (C) 2019-2021 GSI Helmholtzzentrum fuer Schwerionenforschung, Darmstadt
2 SPDX-License-Identifier: GPL-3.0-only
3 Authors: David Emschermann [committer], Nikolay Karpushkin */
4
5// -----------------------------------------------------------------------------
6// ----- -----
7// ----- CbmMcbm2018MonitorAlgoPsd -----
8// ----- Created 26.09.2019 by N.Karpushkin -----
9// ----- based on CbmMcbm2018MonitorAlgoT0 by P.-A. Loizeau -----
10// ----- -----
11// -----------------------------------------------------------------------------
12
14
15#include "CbmFlesHistosTools.h"
17#include "CbmMcbm2018PsdPar.h"
18
19#include "FairRootManager.h"
20#include "FairRun.h"
21#include "FairRunOnline.h"
22#include "FairRuntimeDb.h"
23#include <Logger.h>
24
25#include "TCanvas.h"
26#include "TGraph.h"
27#include "TH1.h"
28#include "TH2.h"
29#include "TList.h"
30#include "TPaveStats.h"
31#include "TProfile.h"
32#include "TROOT.h"
33#include "TString.h"
34
35#include <fstream>
36#include <iomanip>
37#include <iostream>
38
39#include <stdint.h>
40
41// -------------------------------------------------------------------------
46
51
52// -------------------------------------------------------------------------
54{
55 LOG(info) << "Initializing mCBM Psd 2019 monitor algo";
56
57 return kTRUE;
58}
66
67// -------------------------------------------------------------------------
69{
70 LOG(info) << "Init parameter containers for CbmMcbm2018MonitorAlgoPsd";
71 Bool_t initOK = ReInitContainers();
72
73 return initOK;
74}
76{
77 LOG(info) << "**********************************************";
78 LOG(info) << "ReInit parameter containers for CbmMcbm2018MonitorAlgoPsd";
79
80 fUnpackPar = (CbmMcbm2018PsdPar*) fParCList->FindObject("CbmMcbm2018PsdPar");
81 if (nullptr == fUnpackPar) return kFALSE;
82
83 Bool_t initOK = InitParameters();
84
85 return initOK;
86}
88{
89 if (nullptr == fParCList) fParCList = new TList();
90 fUnpackPar = new CbmMcbm2018PsdPar("CbmMcbm2018PsdPar");
92
93 return fParCList;
94}
96{
98 LOG(info) << "Data Version: " << fuRawDataVersion;
99
101 LOG(info) << "Nr. of Tof GDPBs: " << fuNrOfGdpbs;
102
104 LOG(info) << "Nr. of FEEs per Psd GDPB: " << fuNrOfFeePerGdpb;
105
107 LOG(info) << "Nr. of channels per FEE: " << fuNrOfChannelsPerFee;
108
110 LOG(info) << "Nr. of channels per GDPB: " << fuNrOfChannelsPerGdpb;
111
112 fGdpbIdIndexMap.clear();
113 for (UInt_t i = 0; i < fuNrOfGdpbs; ++i) {
115 LOG(info) << "GDPB Id of PSD " << i << " : " << std::hex << fUnpackPar->GetGdpbId(i) << std::dec;
116 } // for( UInt_t i = 0; i < fuNrOfGdpbs; ++i )
117
118
119 return kTRUE;
120}
121// -------------------------------------------------------------------------
122
123void CbmMcbm2018MonitorAlgoPsd::AddMsComponentToList(size_t component, UShort_t usDetectorId)
124{
126 for (UInt_t uCompIdx = 0; uCompIdx < fvMsComponentsList.size(); ++uCompIdx)
127 if (component == fvMsComponentsList[uCompIdx]) return;
128
130 fvMsComponentsList.push_back(component);
131
132 LOG(info) << "CbmMcbm2018MonitorAlgoPsd::AddMsComponentToList => Component " << component << " with detector ID 0x"
133 << std::hex << usDetectorId << std::dec << " added to list";
134}
135// -------------------------------------------------------------------------
136
137Bool_t CbmMcbm2018MonitorAlgoPsd::ProcessTs(const fles::Timeslice& ts)
138{
139 fulCurrentTsIdx = ts.index();
140 fdTsStartTime = static_cast<Double_t>(ts.descriptor(0, 0).idx);
141
143 if (0 == fulCurrentTsIdx) {
144 LOG(info) << "Reseting Histos for a new run";
146 return kTRUE;
147 }
148
150 if (-1.0 == fdTsCoreSizeInNs) {
151 fuNbCoreMsPerTs = ts.num_core_microslices();
152 fuNbOverMsPerTs = ts.num_microslices(0) - ts.num_core_microslices();
155 LOG(info) << "Timeslice parameters: each TS has " << fuNbCoreMsPerTs << " Core MS and " << fuNbOverMsPerTs
156 << " Overlap MS, for a core duration of " << fdTsCoreSizeInNs << " ns and a full duration of "
157 << fdTsFullSizeInNs << " ns";
158
162 LOG(info) << "In each TS " << fuNbMsLoop << " MS will be looped over";
163 } // if( -1.0 == fdTsCoreSizeInNs )
164
166 for (fuMsIndex = 0; fuMsIndex < fuNbMsLoop; fuMsIndex++) {
168 for (UInt_t uMsCompIdx = 0; uMsCompIdx < fvMsComponentsList.size(); ++uMsCompIdx) {
169 UInt_t uMsComp = fvMsComponentsList[uMsCompIdx];
170
171 if (kFALSE == ProcessMs(ts, uMsComp, fuMsIndex)) {
172 LOG(error) << "Failed to process ts " << fulCurrentTsIdx << " MS " << fuMsIndex << " for component " << uMsComp;
173 return kFALSE;
174 } // if( kFALSE == ProcessMs( ts, uMsCompIdx, fuMsIndex ) )
175 } // for( UInt_t uMsCompIdx = 0; uMsCompIdx < fvMsComponentsList.size(); ++uMsCompIdx )
176
178
179 } // for( fuMsIndex = 0; fuMsIndex < uNbMsLoop; fuMsIndex ++ )
180
182
183
185 if (fbMonitorMode) {
186 if (kFALSE == FillHistograms()) {
187 LOG(error) << "Failed to fill histos in ts " << fulCurrentTsIdx;
188 return kFALSE;
189 } // if( kFALSE == FillHistograms() )
190 } // if( fbMonitorMode )
191
192 return kTRUE;
193}
194
195Bool_t CbmMcbm2018MonitorAlgoPsd::ProcessMs(const fles::Timeslice& ts, size_t uMsCompIdx, size_t uMsIdx)
196{
197 auto msDescriptor = ts.descriptor(uMsCompIdx, uMsIdx);
198 fuCurrentEquipmentId = msDescriptor.eq_id;
199 const uint8_t* msContent = reinterpret_cast<const uint8_t*>(ts.content(uMsCompIdx, uMsIdx));
200
201 uint32_t uSize = msDescriptor.size;
202 fulCurrentMsIdx = msDescriptor.idx;
203
204 fdMsTime = (1e-9) * static_cast<double>(fulCurrentMsIdx);
205
206 LOG(debug) << "Microslice: " << fulCurrentMsIdx << " from EqId " << std::hex << fuCurrentEquipmentId << std::dec
207 << " has size: " << uSize;
208
209 if (0 == fvbMaskedComponents.size()) fvbMaskedComponents.resize(ts.num_components(), kFALSE);
210
211 fuCurrDpbId = static_cast<uint32_t>(fuCurrentEquipmentId & 0xFFFF);
212
214 auto it = fGdpbIdIndexMap.find(fuCurrDpbId);
215 if (it == fGdpbIdIndexMap.end()) {
216 if (kFALSE == fvbMaskedComponents[uMsCompIdx]) {
217 LOG(info) << "---------------------------------------------------------------";
218
219 LOG(info) << FormatMsHeaderPrintout(msDescriptor);
220 LOG(warning) << "Could not find the gDPB index for AFCK id 0x" << std::hex << fuCurrDpbId << std::dec
221 << " in timeslice " << fulCurrentTsIdx << " in microslice " << uMsIdx << " component " << uMsCompIdx
222 << "\n"
223 << "If valid this index has to be added in the PSD "
224 "parameter file in the DbpIdArray field";
225 fvbMaskedComponents[uMsCompIdx] = kTRUE;
226 } // if( kFALSE == fvbMaskedComponents[ uMsComp ] )
227 else
228 return kTRUE;
229
232
233
234 return kFALSE;
235
236 } // if( it == fGdpbIdIndexMap.end() )
237 else
239
241 if (0 == fuCurrDpbIdx) {
243 if (1.0 < fdMsTime - fdLastSecondTime) {
246 fbSpillOn = kFALSE;
249 } // if( fbSpillOn && fuCountsLastSecond < kuOffSpillCountLimit )
251 fbSpillOn = kTRUE;
252
255 } // if( 1 < fdMsTime - fdLastSecondTime )
256 } // if( 0 == fuCurrDpbIdx )
257
260
265 } // else if( fuHistoryHistoSize < fdMsTime - fdStartTime )
266
267 // If MS time is less than start time print error and return
268 if (fdMsTime - fdStartTime < 0) {
269 LOG(error) << "negative time! ";
271 return kFALSE;
272 }
273
274 // If not integer number of message in input buffer, print warning/error
275 if (0 != (uSize % kuBytesPerMessage))
276 LOG(error) << "The input microslice buffer does NOT "
277 << "contain only complete nDPB messages!";
278
279 // Compute the number of complete messages in the input microslice buffer
280 uint32_t uNbMessages = (uSize - (uSize % kuBytesPerMessage)) / kuBytesPerMessage;
281
282 // Prepare variables for the loop on contents
283 const uint64_t* pInBuff = reinterpret_cast<const uint64_t*>(msContent);
284
285 if (fair::Logger::Logging(fair::Severity::debug)) {
286 if (uNbMessages != 0) printf("\n\n");
287 for (uint32_t uIdx = 0; uIdx < uNbMessages; uIdx++) {
288 // Fill message
289 uint64_t ulData = static_cast<uint64_t>(pInBuff[uIdx]);
290 printf("%016llx\n", (long long int) ulData);
291 }
292 }
293
294 switch (fuRawDataVersion) {
295 // --------------------------------------------------------------------------------------------------
296 case 100: {
297
298 PsdDataV100::PsdGbtReader PsdReader(pInBuff);
299 if (fair::Logger::Logging(fair::Severity::debug)) PsdReader.SetPrintOutMode(true);
300 // every 80bit gbt word is decomposed into two 64bit words
301 if (uNbMessages > 1) {
302 while (PsdReader.GetTotalGbtWordsRead() < uNbMessages) {
303 int ReadResult = PsdReader.ReadMs();
304 if (ReadResult == 0) {
307
308 //hit loop
309 for (uint64_t hit_iter = 0; hit_iter < PsdReader.VectHitHdr.size(); hit_iter++) {
310 if (PsdReader.VectPackHdr.size() != PsdReader.VectHitHdr.size()) {
311 LOG(error) << "Different vector headers sizes!"
312 << " in VectPackHdr " << PsdReader.VectPackHdr.size() << " in VectHitHdr "
313 << PsdReader.VectHitHdr.size() << "\n";
314 break;
315 }
316
317 uint16_t uHitChannel = PsdReader.VectHitHdr.at(hit_iter).uHitChannel;
318 //uint8_t uLinkIndex = PsdReader.VectPackHdr.at(hit_iter).uLinkIndex;
319 uint32_t uSignalCharge = PsdReader.VectHitHdr.at(hit_iter).uSignalCharge;
320 uint16_t uZeroLevel = PsdReader.VectHitHdr.at(hit_iter).uZeroLevel;
321 uint32_t uAccum = PsdReader.VectHitHdr.at(hit_iter).uFeeAccum;
322 //double dHitTime = (double) fulCurrentMsIdx + PsdReader.VectPackHdr.at(hit_iter).uAdcTime * 12.5; //in ns
323 //double dHitTime = PsdReader.MsHdr.ulMicroSlice*1000. + PsdReader.VectPackHdr.at(hit_iter).uAdcTime*12.5; //in ns
324 std::vector<uint16_t> uWfm = PsdReader.VectHitData.at(hit_iter).uWfm;
325 uSignalCharge /= (int) kfAdc_to_mV; // ->now in mV
326
327 fhAdcTime->Fill(PsdReader.VectPackHdr.at(hit_iter).uAdcTime);
328
329 //uHitChannel numerated from 0
330 if (uHitChannel >= kuNbChanPsd) {
331 LOG(error) << "hit channel number out of range! channel index: " << uHitChannel
332 << " max: " << GetNbChanPsd();
333 break;
334 }
335
336 //Pack header
337 fhHitChargeMap->Fill(uHitChannel, uSignalCharge);
338 fhChanHitMapEvo->Fill(uHitChannel, fdMsTime - fdStartTime); //should be placed map(channel)
339 if (fbMonitorChanMode) {
340 fvhHitChargeChan[uHitChannel]->Fill(uSignalCharge);
341 fvhHitZeroLevelChan[uHitChannel]->Fill(uZeroLevel);
342 fvhHitZLChanEvo[uHitChannel]->Fill(fdMsTime - fdStartTime, uZeroLevel);
343 fvhHitFAChanEvo[uHitChannel]->Fill(fdMsTime - fdStartTime, uAccum);
344
345 //Hit data
346 double dHitAmlpitude = 0;
347 double dHitChargeWfm = 0;
348 if (fbMonitorWfmMode) fvhHitWfmChan[uHitChannel]->Reset();
349 if (fbMonitorFitMode) fvhHitFitWfmChan[uHitChannel]->Reset();
350
351 if (!uWfm.empty()) {
352 dHitChargeWfm = std::accumulate(uWfm.begin(), uWfm.end(), 0);
353 dHitChargeWfm -= uZeroLevel * uWfm.size();
354 auto const max_iter = std::max_element(uWfm.begin(), uWfm.end());
355 assert(max_iter != uWfm.end());
356 if (max_iter == uWfm.end()) break;
357 //uint8_t hit_time_max = std::distance(uWfm.begin(), max_iter);
358 dHitAmlpitude = *max_iter - uZeroLevel;
359 dHitAmlpitude /= kfAdc_to_mV;
360 dHitChargeWfm /= kfAdc_to_mV;
361 fvhHitAmplChan[uHitChannel]->Fill(dHitAmlpitude);
362 fvhHitChargeByWfmChan[uHitChannel]->Fill(dHitChargeWfm);
363
364 if (fbMonitorWfmMode) {
365 fvhHitLPChanEvo[uHitChannel]->Fill(fdMsTime - fdStartTime, uWfm.back());
366 for (UInt_t wfm_iter = 0; wfm_iter < uWfm.size(); wfm_iter++)
367 fvhHitWfmChan[uHitChannel]->Fill(wfm_iter, uWfm.at(wfm_iter));
368 fvhHitWfmChan[uHitChannel]->SetTitle(
369 Form("Waveform channel %03u charge %0u zero level %0u; Time [adc "
370 "counts]; Amplitude [adc counts]",
371 uHitChannel, uSignalCharge, uZeroLevel));
372 for (uint8_t i = 0; i < kuNbWfmRanges; ++i) {
373 if (uSignalCharge > kvuWfmRanges.at(i) && uSignalCharge < kvuWfmRanges.at(i + 1)) {
374 UInt_t uFlatIndexOfChange = i * kuNbChanPsd + uHitChannel;
375
376 UInt_t uWfmExampleIter = kvuWfmInRangeToChangeChan.at(uFlatIndexOfChange);
377 UInt_t uFlatIndexHisto =
378 uWfmExampleIter * kuNbWfmRanges * kuNbChanPsd + i * kuNbChanPsd + uHitChannel;
379 fv3hHitWfmFlattenedChan[uFlatIndexHisto]->Reset();
380
381 for (UInt_t wfm_iter = 0; wfm_iter < uWfm.size(); wfm_iter++)
382 fv3hHitWfmFlattenedChan[uFlatIndexHisto]->Fill(wfm_iter, uWfm.at(wfm_iter));
383 fv3hHitWfmFlattenedChan[uFlatIndexHisto]->SetTitle(
384 Form("Waveform channel %03u charge %0u zero level %0u; Time "
385 "[adc counts]; Amplitude [adc counts]",
386 uHitChannel, uSignalCharge, uZeroLevel));
387
388 kvuWfmInRangeToChangeChan.at(uFlatIndexOfChange)++;
389 if (kvuWfmInRangeToChangeChan.at(uFlatIndexOfChange) == kuNbWfmExamples)
390 kvuWfmInRangeToChangeChan.at(uFlatIndexOfChange) = 0;
391
392 } // if( uSignalCharge > kvuWfmRanges.at(i) && uSignalCharge < kvuWfmRanges.at(i+1) )
393 } // for (uint8_t i=0; i<kuNbWfmRanges; ++i)
394 } //if (fbMonitorWfmMode)
395 } //if (!uWfm.empty())
396
397
398 if (fbMonitorFitMode && !uWfm.empty()) {
399 int gate_beg = 0;
400 int gate_end = 10; //uWfm.size() - 1;
401 PsdSignalFitting::PronyFitter Pfitter(2, 2, gate_beg, gate_end);
402
403 Pfitter.SetDebugMode(0);
404 Pfitter.SetWaveform(uWfm, uZeroLevel);
405 int SignalBeg = 4;
406 std::complex<float> first_fit_harmonic = {0.72, 0.0};
407 std::complex<float> second_fit_harmonic = {0.38, -0.0};
408 Pfitter.SetExternalHarmonics(first_fit_harmonic, second_fit_harmonic);
409 Int_t best_signal_begin = Pfitter.ChooseBestSignalBegin(SignalBeg - 1, SignalBeg + 1);
410 Pfitter.SetSignalBegin(best_signal_begin);
411 Pfitter.CalculateFitAmplitudes();
412
413 Pfitter.SetSignalBegin(best_signal_begin);
414 Pfitter.CalculateFitHarmonics();
415 Pfitter.CalculateFitAmplitudes();
416
417 Float_t fit_integral = Pfitter.GetIntegral(gate_beg, uWfm.size() - 1) / kfAdc_to_mV;
418 Float_t fit_R2 = Pfitter.GetRSquare(gate_beg, uWfm.size() - 1);
419
420 std::complex<float>* harmonics = Pfitter.GetHarmonics();
421 std::vector<uint16_t> uFitWfm = Pfitter.GetFitWfm();
422 for (UInt_t wfm_iter = 0; wfm_iter < uFitWfm.size(); wfm_iter++)
423 fvhHitFitWfmChan[uHitChannel]->Fill(wfm_iter, uFitWfm.at(wfm_iter));
424 fvhHitWfmChan[uHitChannel]->SetTitle(Form("Waveform channel %03u charge %0u zero level %0u R2 %.5f; "
425 "Time [adc counts]; Amplitude [adc counts]",
426 uHitChannel, uSignalCharge, uZeroLevel, fit_R2));
427
428 fvhFitQaChan[uHitChannel]->Fill(fit_integral, fit_R2);
429
430 if (fit_R2 > 0.02) continue;
431 fvhFitHarmonic1Chan[uHitChannel]->Fill(std::real(harmonics[1]), std::imag(harmonics[1]));
432 fvhFitHarmonic2Chan[uHitChannel]->Fill(std::real(harmonics[2]), std::imag(harmonics[2]));
433 } //if (fbMonitorFitMode && !uWfm.empty())
434 } //if (fbMonitorChanMode)
435 } // for(int hit_iter = 0; hit_iter < PsdReader.EvHdrAb.uHitsNumber; hit_iter++)
436 }
437 else if (ReadResult == 1) {
438 LOG(error) << "no pack headers in message!";
439 break;
440 }
441 else if (ReadResult == 2) {
442 LOG(error) << "wrong channel! In header: " << PsdReader.HitHdr.uHitChannel;
443 break;
444 }
445 else if (ReadResult == 3) {
446 LOG(error) << "check number of waveform points! In header: " << PsdReader.HitHdr.uWfmWords - 1;
447 break;
448 }
449 else {
450 LOG(error) << "PsdGbtReader.ReadEventFles() didn't return expected values";
451 break;
452 }
453
454 } // while(PsdReader.GetTotalGbtWordsRead()<uNbMessages)
455
456
457 if (uNbMessages != PsdReader.GetTotalGbtWordsRead()) {
458 fbPsdMissedData = kTRUE;
459 LOG(error) << "Wrong amount of messages read!"
460 << " in microslice " << uNbMessages << " by PsdReader " << PsdReader.GetTotalGbtWordsRead();
461
463 std::ofstream error_log(Form("%llu_errorlog.txt", fulCurrentMsIdx), std::ofstream::out);
464 for (uint32_t uIdx = 0; uIdx < uNbMessages; uIdx++) {
465 uint64_t ulData = static_cast<uint64_t>(pInBuff[uIdx]);
466 error_log << Form("%016llx\n", (long long int) ulData);
467 }
468 error_log.close();
469 fbFirstPackageError = kFALSE;
470 printf("Written file %llu_errorlog.txt\n", fulCurrentMsIdx);
471 }
472 }
473
474 fuMsgsCntInMs += uNbMessages;
476 fuLostMsgsCntInMs += uNbMessages - PsdReader.GetTotalGbtWordsRead();
477 } //if (uNbMessages > 1)
478
481
482 break;
483
484 } // case 1
485 // --------------------------------------------------------------------------------------------------
486 // --------------------------------------------------------------------------------------------------
487 case 000: {
488
489 PsdDataV000::PsdGbtReader PsdReader(pInBuff);
490 if (fair::Logger::Logging(fair::Severity::debug)) PsdReader.SetPrintOutMode(true);
491 if (uNbMessages > 1) {
492 while (PsdReader.GetTotalGbtWordsRead() < uNbMessages) {
493 int ReadResult = PsdReader.ReadEventFles();
494 if (PsdReader.EvHdrAb.uHitsNumber > kuNbChanPsd) {
495 LOG(error) << "too many triggered channels! In header: " << PsdReader.EvHdrAb.uHitsNumber
496 << " in PSD: " << GetNbChanPsd();
497 break;
498 }
499
500 if (ReadResult == 0) {
502 fhAdcTime->Fill(PsdReader.EvHdrAc.uAdcTime);
504
505 //hit loop
506 for (int hit_iter = 0; hit_iter < PsdReader.EvHdrAb.uHitsNumber; hit_iter++) {
507 UInt_t uHitChannel = PsdReader.VectHitHdr.at(hit_iter).uHitChannel;
508 UInt_t uSignalCharge = PsdReader.VectHitHdr.at(hit_iter).uSignalCharge;
509 UInt_t uZeroLevel = PsdReader.VectHitHdr.at(hit_iter).uZeroLevel;
510 std::vector<uint16_t> uWfm = PsdReader.VectHitData.at(hit_iter).uWfm;
511 uSignalCharge /= (int) kfAdc_to_mV; // ->now in mV
512
513 if (uHitChannel >= kuNbChanPsd) //uHitChannel numerated from 0
514 {
515 LOG(error) << "hit channel number out of range! channel index: " << uHitChannel
516 << " max: " << GetNbChanPsd();
517 break;
518 }
519 //Hit header
520 fhHitChargeMap->Fill(uHitChannel, uSignalCharge);
521 fhHitMapEvo->Fill(uHitChannel, fdMsTime - fdStartTime);
522 fhChanHitMapEvo->Fill(uHitChannel, fdMsTime - fdStartTime); //should be placed map(ch)
523
524 if (fbMonitorChanMode) {
525
526 fvhHitChargeChan[uHitChannel]->Fill(uSignalCharge);
527 fvhHitZeroLevelChan[uHitChannel]->Fill(uZeroLevel);
528
529 //Hit data
530 double dHitAmlpitude = 0;
531 double dHitChargeWfm = 0;
532 if (fbMonitorWfmMode) fvhHitWfmChan[uHitChannel]->Reset();
533 if (fbMonitorFitMode) fvhHitFitWfmChan[uHitChannel]->Reset();
534
535 dHitChargeWfm = std::accumulate(uWfm.begin(), uWfm.end(), 0);
536 dHitChargeWfm -= uZeroLevel * uWfm.size();
537 auto const max_iter = std::max_element(uWfm.begin(), uWfm.end());
538 assert(max_iter != uWfm.end());
539 if (max_iter == uWfm.end()) break;
540 //uint8_t hit_time_max = std::distance(uWfm.begin(), max_iter);
541 dHitAmlpitude = *max_iter - uZeroLevel;
542 dHitAmlpitude /= kfAdc_to_mV;
543 dHitChargeWfm /= kfAdc_to_mV;
544 fvhHitAmplChan[uHitChannel]->Fill(dHitAmlpitude);
545 fvhHitChargeByWfmChan[uHitChannel]->Fill(dHitChargeWfm);
546
547 if (fbMonitorWfmMode) {
548 for (UInt_t wfm_iter = 0; wfm_iter < uWfm.size(); wfm_iter++)
549 fvhHitWfmChan[uHitChannel]->Fill(wfm_iter, uWfm.at(wfm_iter));
550 fvhHitWfmChan[uHitChannel]->SetTitle(
551 Form("Waveform channel %03u charge %0u zero level %0u; Time [adc "
552 "counts]; Amplitude [adc counts]",
553 uHitChannel, uSignalCharge, uZeroLevel));
554 for (uint8_t i = 0; i < kuNbWfmRanges; ++i) {
555 if (uSignalCharge > kvuWfmRanges.at(i) && uSignalCharge < kvuWfmRanges.at(i + 1)) {
556 UInt_t uFlatIndexOfChange = i * kuNbChanPsd + uHitChannel;
557
558 UInt_t uWfmExampleIter = kvuWfmInRangeToChangeChan.at(uFlatIndexOfChange);
559 UInt_t uFlatIndexHisto =
560 uWfmExampleIter * kuNbWfmRanges * kuNbChanPsd + i * kuNbChanPsd + uHitChannel;
561 fv3hHitWfmFlattenedChan[uFlatIndexHisto]->Reset();
562
563 for (UInt_t wfm_iter = 0; wfm_iter < uWfm.size(); wfm_iter++)
564 fv3hHitWfmFlattenedChan[uFlatIndexHisto]->Fill(wfm_iter, uWfm.at(wfm_iter));
565 fv3hHitWfmFlattenedChan[uFlatIndexHisto]->SetTitle(
566 Form("Waveform channel %03u charge %0u zero level %0u; Time "
567 "[adc counts]; Amplitude [adc counts]",
568 uHitChannel, uSignalCharge, uZeroLevel));
569
570 kvuWfmInRangeToChangeChan.at(uFlatIndexOfChange)++;
571 if (kvuWfmInRangeToChangeChan.at(uFlatIndexOfChange) == kuNbWfmExamples)
572 kvuWfmInRangeToChangeChan.at(uFlatIndexOfChange) = 0;
573
574 } // if( uSignalCharge > kvuWfmRanges.at(i) && uSignalCharge < kvuWfmRanges.at(i+1) )
575 } // for (uint8_t i=0; i<kuNbWfmRanges; ++i)
576 } //if (fbMonitorWfmMode)
577
578 if (fbMonitorFitMode) {
579 int gate_beg = 0;
580 int gate_end = uWfm.size() - 1;
581 PsdSignalFitting::PronyFitter Pfitter(2, 2, gate_beg, gate_end);
582
583 Pfitter.SetDebugMode(0);
584 Pfitter.SetWaveform(uWfm, uZeroLevel);
585 int SignalBeg = 2;
586 Int_t best_signal_begin = Pfitter.ChooseBestSignalBeginHarmonics(SignalBeg - 1, SignalBeg + 1);
587
588 Pfitter.SetSignalBegin(best_signal_begin);
589 Pfitter.CalculateFitHarmonics();
590 Pfitter.CalculateFitAmplitudes();
591
592 Float_t fit_integral = Pfitter.GetIntegral(gate_beg, gate_end) / kfAdc_to_mV;
593 Float_t fit_R2 = Pfitter.GetRSquare(gate_beg, gate_end);
594
595 std::complex<float>* harmonics = Pfitter.GetHarmonics();
596 std::vector<uint16_t> uFitWfm = Pfitter.GetFitWfm();
597 for (UInt_t wfm_iter = 0; wfm_iter < uFitWfm.size(); wfm_iter++)
598 fvhHitFitWfmChan[uHitChannel]->Fill(wfm_iter, uFitWfm.at(wfm_iter));
599 fvhHitWfmChan[uHitChannel]->SetTitle(Form("Waveform channel %03u charge %0u zero level %0u R2 %.5f; "
600 "Time [adc counts]; Amplitude [adc counts]",
601 uHitChannel, uSignalCharge, uZeroLevel, fit_R2));
602
603 fvhFitQaChan[uHitChannel]->Fill(fit_integral, fit_R2);
604
605 if (fit_R2 > 0.02) continue;
606 fvhFitHarmonic1Chan[uHitChannel]->Fill(std::real(harmonics[1]), std::imag(harmonics[1]));
607 fvhFitHarmonic2Chan[uHitChannel]->Fill(std::real(harmonics[2]), std::imag(harmonics[2]));
608 } //if (fbMonitorFitMode)
609 } //if (fbMonitorChanMode)
610
611 } // for(int hit_iter = 0; hit_iter < PsdReader.EvHdrAb.uHitsNumber; hit_iter++)
612 }
613 else if (ReadResult == 1) {
614 LOG(error) << "no event headers in message!";
615 break;
616 }
617 else if (ReadResult == 2) {
618 LOG(error) << "check number of waveform points! In header: " << PsdReader.HitHdr.uWfmPoints
619 << " should be: " << 8;
620 break;
621 }
622 else if (ReadResult == 3) {
623 LOG(error) << "wrong amount of hits read! In header: " << PsdReader.EvHdrAb.uHitsNumber
624 << " in hit vector: " << PsdReader.VectHitHdr.size();
625 break;
626 }
627 else {
628 LOG(error) << "PsdGbtReader.ReadEventFles() didn't return expected values";
629 break;
630 }
631
632 } // while(PsdReader.GetTotalGbtWordsRead()<uNbMessages)
633
634 if (uNbMessages != PsdReader.GetTotalGbtWordsRead()) {
635 fbPsdMissedData = kTRUE;
636 LOG(error) << "Wrong amount of messages read!"
637 << " in microslice " << uNbMessages << " by PsdReader " << PsdReader.GetTotalGbtWordsRead();
638
640 std::ofstream error_log(Form("%llu_errorlog.txt", fulCurrentMsIdx), std::ofstream::out);
641 for (uint32_t uIdx = 0; uIdx < uNbMessages; uIdx++) {
642 uint64_t ulData = static_cast<uint64_t>(pInBuff[uIdx]);
643 error_log << Form("%016llx\n", (long long int) ulData);
644 }
645 error_log.close();
646 fbFirstPackageError = kFALSE;
647 printf("Written file %llu_errorlog.txt\n", fulCurrentMsIdx);
648 }
649 }
650
651 if (fulCurrentMsIdx != PsdReader.EvHdrAb.ulMicroSlice)
652 LOG(error) << "Wrong MS index!"
653 << " in microslice " << fulCurrentMsIdx << " by PsdReader " << PsdReader.EvHdrAb.ulMicroSlice
654 << "\n";
655
656 fuMsgsCntInMs += uNbMessages;
658 fuLostMsgsCntInMs += uNbMessages - PsdReader.GetTotalGbtWordsRead();
659
660 } //if(uNbMessages > 1)
661
662
665
666 break;
667 } // case 0
668 // --------------------------------------------------------------------------------------------------
669
670 } // switch
671
673 else {
676 }
677
678 return kTRUE;
679}
680
682{
683 std::string sFolder = "MoniPsd";
684 std::string sFitFolder = "PronyFit";
685 LOG(info) << "create Histos for PSD monitoring ";
686
688 uint32_t iNbBinsLog = 0;
690 std::vector<double> dBinsLogVector = GenerateLogBinArray(2, 3, 1, iNbBinsLog);
691 double* dBinsLog = dBinsLogVector.data();
692
693 fhHitChargeMap = new TH2I("hHitChargeMap", "Map of hits charges in PSD detector; Chan; Charge [mV]", kuNbChanPsd, 0.,
695 fhHitMapEvo = new TH2I("hHitMapEvo",
696 "Map of hits in PSD detector electronics vs time in "
697 "run; Chan; Time in run [s]; Hits Count []",
699 fhChanHitMapEvo = new TH2I("hChanHitMapEvo",
700 "Map of hits in PSD detector vs time in run; "
701 "Chan; Time in run [s]; Hits Count []",
703
704
705 fhMissedData = new TH1I("hMissedData", "PSD Missed data", 2, 0, 2);
706
707 fhAdcTime = new TH1I("hAdcTime", "ADC time; Adc time []", 100, 0, 160000);
708
709 fhMsLengthEvo = new TH2I("hMsLengthEvo", "Evolution of MS length vs time in run; Time in run [s]; MS length [ns]",
710 fuHistoryHistoSize, 0, fuHistoryHistoSize, (Int_t) 1e2, 0, 1e8);
711
712 fhMsgsCntPerMsEvo = new TH2I("hMsgsCntPerMsEvo",
713 "Evolution of TotalMsgs counts, per MS vs time in run; Time in "
714 "run [s]; TotalMsgs Count/MS []; MS",
715 fuHistoryHistoSize, 0, fuHistoryHistoSize, iNbBinsLog, dBinsLog);
716 fhReadMsgsCntPerMsEvo = new TH2I("ReadMsgsCntPerMsEvo",
717 "Evolution of ReadMsgs counts, per MS vs time in run; Time in run "
718 "[s]; ReadMsgs Count/MS []; MS",
719 fuHistoryHistoSize, 0, fuHistoryHistoSize, iNbBinsLog, dBinsLog);
720 fhLostMsgsCntPerMsEvo = new TH2I("hLostMsgsCntPerMsEvo",
721 "Evolution of LostMsgs counts, per MS vs time in run; Time in run "
722 "[s]; LostMsgs Count/MS []; MS",
723 fuHistoryHistoSize, 0, fuHistoryHistoSize, iNbBinsLog, dBinsLog);
724 fhReadEvtsCntPerMsEvo = new TH2I("hReadEvtCntPerMsEvo",
725 "Evolution of ReadEvents counts, per MS vs time in run; Time in "
726 "run [s]; ReadEvents Count/MS []; MS",
727 fuHistoryHistoSize, 0, fuHistoryHistoSize, iNbBinsLog, dBinsLog);
728
733
735 AddHistoToVector(fhAdcTime, sFolder);
737
742
743 /*******************************************************************/
744 if (fbMonitorChanMode) {
745
746 for (UInt_t uChan = 0; uChan < kuNbChanPsd; ++uChan) {
747
748 fvhHitZLChanEvo[uChan] = new TH2I(
749 Form("hHitZLChanEvo%03u", uChan),
750 Form("Hits ZeroLevel evolution for channel %03u; Time in run [s]; ZeroLevel [adc counts]", uChan),
752 fvhHitZLChanEvo[uChan]->SetMarkerColor(kRed);
753 fvhHitLPChanEvo[uChan] = new TH2I(
754 Form("hHitLPChanEvo%03u", uChan),
755 Form("Hits LastPoint evolution for channel %03u; Time in run [s]; ZeroLevel [adc counts]", uChan),
757 fvhHitLPChanEvo[uChan]->SetMarkerColor(kBlue);
758 fvhHitFAChanEvo[uChan] = new TH2I(
759 Form("hHitFAChanEvo%03u", uChan),
760 Form("Hits FeeAccumulator evolution for channel %03u; Time in run [s]; ZeroLevel [adc counts]", uChan),
762 fvhHitFAChanEvo[uChan]->SetMarkerColor(kOrange);
763
764
765 fvhHitChargeChan[uChan] = new TH1I(Form("hHitChargeChan%03u", uChan),
766 Form("Hits charge distribution for channel %03u; Charge [mV]", uChan),
768
769 fvhHitZeroLevelChan[uChan] =
770 new TH1I(Form("hHitZeroLevelChan%03u", uChan),
771 Form("Hits zero level distribution for channel %03u; ZeroLevel [adc counts]", uChan),
772 fviHistoZLArgs.at(0), fviHistoZLArgs.at(1), fviHistoZLArgs.at(2));
773
774 fvhHitAmplChan[uChan] = new TH1F(Form("hHitAmplChan%03u", uChan),
775 Form("Hits amplitude distribution for channel %03u; Amplitude [mV]", uChan),
777
778 fvhHitChargeByWfmChan[uChan] =
779 new TH1I(Form("hHitChargeByWfmChan%03u", uChan),
780 Form("Hits charge by waveform distribution for channel %03u; "
781 "Charge [mV]",
782 uChan),
784
785
786 AddHistoToVector(fvhHitZLChanEvo[uChan], sFolder);
787 AddHistoToVector(fvhHitLPChanEvo[uChan], sFolder);
788 AddHistoToVector(fvhHitFAChanEvo[uChan], sFolder);
789 AddHistoToVector(fvhHitChargeChan[uChan], sFolder);
790 AddHistoToVector(fvhHitZeroLevelChan[uChan], sFolder);
791 AddHistoToVector(fvhHitAmplChan[uChan], sFolder);
793
794 if (fbMonitorWfmMode) {
795 fvhHitWfmChan[uChan] = new TH1I(Form("hHitWfmChan%03u", uChan), Form("HitWfmChan%03u", uChan), 32, 0, 32);
796 fvhHitWfmChan[uChan]->SetMarkerStyle(31);
797 fvhHitWfmChan[uChan]->SetMarkerSize(0.5);
798
799 for (UInt_t uWfmRangeIter = 0; uWfmRangeIter < kuNbWfmRanges; uWfmRangeIter++) {
800 for (UInt_t uWfmExampleIter = 0; uWfmExampleIter < kuNbWfmExamples; uWfmExampleIter++) {
801 UInt_t uFlatIndex = uWfmExampleIter * kuNbWfmRanges * kuNbChanPsd + uWfmRangeIter * kuNbChanPsd + uChan;
802 fv3hHitWfmFlattenedChan[uFlatIndex] =
803 new TH1I(Form("hHitWfmChan%03uRange%02uExample%02u", uChan, uWfmRangeIter, uWfmExampleIter),
804 Form("HitWfmChan%03uRange%02uExample%02u", uChan, uWfmRangeIter, uWfmExampleIter), 32, 0, 32);
805
806 } // for( UInt_t uWfmRangeIter = 0; uWfmRangeIter < kuNbWfmRanges; uWfmRangeIter ++)
807 } // for( UInt_t uWfmExampleIter = 0; uWfmExampleIter < kuNbWfmExamples; uWfmExampleIter ++)
808 } // if(fbMonitorWfmMode)
809
810 if (fbMonitorFitMode) {
811
812 fvhHitFitWfmChan[uChan] =
813 new TH1I(Form("hHitFitWfmChan%03u", uChan), Form("HitFitWfmChan%03u", uChan), 32, 0, 32);
814 fvhHitFitWfmChan[uChan]->SetLineColor(kRed);
815 fvhHitFitWfmChan[uChan]->SetLineWidth(2);
816
817 fvhFitHarmonic1Chan[uChan] = new TH2I(
818 Form("hFitHarmonic1Chan%03u", uChan),
819 Form("Waveform fit harmonic 1 for channel %03u; Real part []; Imag part []", uChan), 400, -2, 2, 200, -1, 1);
820 fvhFitHarmonic1Chan[uChan]->SetMarkerColor(kRed);
821
822 fvhFitHarmonic2Chan[uChan] = new TH2I(
823 Form("hFitHarmonic2Chan%03u", uChan),
824 Form("Waveform fit harmonic 2 for channel %03u; Real part []; Imag part []", uChan), 400, -2, 2, 200, -1, 1);
825 fvhFitHarmonic2Chan[uChan]->SetMarkerColor(kBlue);
826
827 fvhFitQaChan[uChan] = new TH2I(
828 Form("hFitQaChan%03u", uChan), Form("Waveform fit QA for channel %03u; Integral [adc counts]; R2 []", uChan),
829 fviHistoChargeArgs.at(0), fviHistoChargeArgs.at(1), fviHistoChargeArgs.at(2), 500, 0, 1);
830
831 AddHistoToVector(fvhFitHarmonic1Chan[uChan], sFitFolder);
832 AddHistoToVector(fvhFitHarmonic2Chan[uChan], sFitFolder);
833 AddHistoToVector(fvhFitQaChan[uChan], sFitFolder);
834
835 } // if(fbMonitorFitMode)
836 } // for( UInt_t uChan = 0; uChan < kuNbChanPsd; ++uChan )
837
838 } // if (fbMonitorChanMode)
839 /*******************************************************************/
840
842 Double_t w = 10;
843 Double_t h = 10;
844
845 /*******************************************************************/
847 fcHitMaps = new TCanvas("cHitMaps", "Hit maps", w, h);
848 fcHitMaps->Divide(2);
849
850 fcHitMaps->cd(1);
851 gPad->SetGridx();
852 gPad->SetGridy();
853 gPad->SetLogz();
854 fhChanHitMapEvo->Draw("colz");
855
856 fcHitMaps->cd(2);
857 gPad->SetGridx();
858 gPad->SetGridy();
859 gPad->SetLogz();
860 fhHitChargeMap->Draw("colz");
861
862 AddCanvasToVector(fcHitMaps, "canvases");
863 /*******************************************************************/
864
865 /*******************************************************************/
867 fcSummary = new TCanvas("cSummary", "Hit maps, Hit rate, Error fraction", w, h);
868 fcSummary->Divide(2, 2);
869
870 fcSummary->cd(1);
871 gPad->SetGridx();
872 gPad->SetGridy();
873 gPad->SetLogz();
874 fhHitChargeMap->Draw("colz");
875
876 fcSummary->cd(2);
877 gPad->SetGridx();
878 gPad->SetGridy();
879 gPad->SetLogz();
880 fhChanHitMapEvo->Draw("colz");
881
882 fcSummary->cd(3);
883 gPad->SetGridx();
884 gPad->SetGridy();
885 gPad->SetLogy();
886 fhMissedData->Draw();
887
888 fcSummary->cd(4);
889 gPad->SetGridx();
890 gPad->SetGridy();
891 gPad->SetLogy();
892 gPad->SetLogz();
893 fhMsLengthEvo->Draw("colz");
894
895 AddCanvasToVector(fcSummary, "canvases");
896 /*******************************************************************/
897
898 /*******************************************************************/
900 fcGenCntsPerMs = new TCanvas("cGenCntsPerMs", "Messages and hit cnt per MS, Error and Evt Loss Fract. per MS ", w, h);
901 fcGenCntsPerMs->Divide(2, 2);
902
903 fcGenCntsPerMs->cd(1);
904 gPad->SetGridx();
905 gPad->SetGridy();
906 gPad->SetLogy();
907 gPad->SetLogz();
908 fhMsgsCntPerMsEvo->Draw("colz");
909
910 fcGenCntsPerMs->cd(2);
911 gPad->SetGridx();
912 gPad->SetGridy();
913 gPad->SetLogy();
914 gPad->SetLogz();
915 fhReadMsgsCntPerMsEvo->Draw("colz");
916
917 fcGenCntsPerMs->cd(3);
918 gPad->SetGridx();
919 gPad->SetGridy();
920 gPad->SetLogy();
921 gPad->SetLogz();
922 fhReadEvtsCntPerMsEvo->Draw("colz");
923
924 fcGenCntsPerMs->cd(4);
925 gPad->SetGridx();
926 gPad->SetGridy();
927 gPad->SetLogy();
928 gPad->SetLogz();
929 fhLostMsgsCntPerMsEvo->Draw("colz");
930
932 /*******************************************************************/
933
934 if (fbMonitorChanMode) {
935
936 /*******************************************************************/
938 fcZLevo = new TCanvas("cZLevo", "ZeroLevel evolaution in all channels", w, h);
939 fcZLevo->DivideSquare(kuNbChanPsd);
940
941 for (UInt_t uChan = 0; uChan < kuNbChanPsd; uChan++) {
942 fcZLevo->cd(1 + uChan);
943 fvhHitZLChanEvo[uChan]->Draw();
944 fvhHitLPChanEvo[uChan]->Draw("same");
945 fvhHitFAChanEvo[uChan]->Draw("same");
946 } // for( UInt_t uChan = 0; uChan < kuNbChanPsd; uChan ++)
947
948 AddCanvasToVector(fcZLevo, "canvases");
949 /*******************************************************************/
950
951 /*******************************************************************/
953 fcChargesFPGA = new TCanvas("cChargesFPGA", "Charges spectra in all channels calculated by FPGA", w, h);
954 fcChargesFPGA->DivideSquare(kuNbChanPsd);
955
956 for (UInt_t uChan = 0; uChan < kuNbChanPsd; uChan++) {
957 fcChargesFPGA->cd(1 + uChan)->SetLogy();
958 fvhHitChargeChan[uChan]->Draw();
959 } // for( UInt_t uChan = 0; uChan < kuNbChanPsd; uChan ++)
960
961 AddCanvasToVector(fcChargesFPGA, "canvases");
962 /*******************************************************************/
963
964 /*******************************************************************/
966 fcChargesWfm = new TCanvas("cChargesWfm", "Charges spectra in all channels calculated over waveform", w, h);
967 fcChargesWfm->DivideSquare(kuNbChanPsd);
968
969 for (UInt_t uChan = 0; uChan < kuNbChanPsd; uChan++) {
970 fcChargesWfm->cd(1 + uChan)->SetLogy();
971 fvhHitChargeByWfmChan[uChan]->Draw();
972 } // for( UInt_t uChan = 0; uChan < kuNbChanPsd; uChan ++)
973
974 AddCanvasToVector(fcChargesWfm, "canvases");
975 /*******************************************************************/
976
977 /*******************************************************************/
979 fcAmplitudes = new TCanvas("cAmplitudes", "Amplitude spectra in all channels", w, h);
980 fcAmplitudes->DivideSquare(kuNbChanPsd);
981
982 for (UInt_t uChan = 0; uChan < kuNbChanPsd; uChan++) {
983 fcAmplitudes->cd(1 + uChan)->SetLogy();
984 fvhHitAmplChan[uChan]->Draw();
985 } // for( UInt_t uChan = 0; uChan < kuNbChanPsd; uChan ++)
986
987 AddCanvasToVector(fcAmplitudes, "canvases");
988 /*******************************************************************/
989
990 /*******************************************************************/
992 fcZeroLevels = new TCanvas("cZeroLevels", "Zero Level spectra in all channels", w, h);
993 fcZeroLevels->DivideSquare(kuNbChanPsd);
994
995 for (UInt_t uChan = 0; uChan < kuNbChanPsd; uChan++) {
996 fcZeroLevels->cd(1 + uChan);
997 fvhHitZeroLevelChan[uChan]->Draw();
998 } // for( UInt_t uChan = 0; uChan < kuNbChanPsd; uChan ++)
999
1000 AddCanvasToVector(fcZeroLevels, "canvases");
1001 /*******************************************************************/
1002
1003 if (fbMonitorWfmMode) {
1004
1005 /*******************************************************************/
1007 fcWfmsAllChannels = new TCanvas("cWfmsAllChannels", "Last waveforms in PSD fired channels", w, h);
1008 fcWfmsAllChannels->DivideSquare(kuNbChanPsd);
1009
1010 for (UInt_t uChan = 0; uChan < kuNbChanPsd; uChan++) {
1011 fcWfmsAllChannels->cd(1 + uChan);
1012 if (!fbMonitorFitMode) fvhHitWfmChan[uChan]->Draw("HIST LP");
1013 if (fbMonitorFitMode) {
1014 fvhHitWfmChan[uChan]->Draw("HIST P");
1015 fvhHitFitWfmChan[uChan]->Draw("L SAME");
1016 }
1017 } // for( UInt_t uChan = 0; uChan < kuNbChanPsd; uChan ++)
1018
1020 /*******************************************************************/
1021
1022 /*******************************************************************/
1024 for (UInt_t uChan = 0; uChan < kuNbChanPsd; uChan++) {
1025 fvcWfmsChan[uChan] =
1026 new TCanvas(Form("cWfmsChan%03u", uChan), Form("Canvas with last waveforms in channel %03u", uChan), w, h);
1027 fvcWfmsChan[uChan]->Divide(kuNbWfmExamples, kuNbWfmRanges);
1028 UInt_t uHisto = 0;
1029
1030 for (UInt_t uWfmRangeIter = 0; uWfmRangeIter < kuNbWfmRanges; uWfmRangeIter++) {
1031 for (UInt_t uWfmExampleIter = 0; uWfmExampleIter < kuNbWfmExamples; uWfmExampleIter++) {
1032 fvcWfmsChan[uChan]->cd(1 + uHisto);
1033 UInt_t uFlatIndex = uWfmExampleIter * kuNbWfmRanges * kuNbChanPsd + uWfmRangeIter * kuNbChanPsd + uChan;
1034 fv3hHitWfmFlattenedChan[uFlatIndex]->Draw("HIST LP");
1035 uHisto++;
1036 } // for( UInt_t uWfmRangeIter = 0; uWfmRangeIter < kuNbWfmRanges; uWfmRangeIter ++)
1037 } // for( UInt_t uWfmExampleIter = 0; uWfmExampleIter < kuNbWfmExamples; uWfmExampleIter ++)
1038
1039 AddCanvasToVector(fvcWfmsChan[uChan], "waveforms");
1040 } // for( UInt_t uChan = 0; uChan < kuNbChanPsd; uChan ++)
1041
1042 /*******************************************************************/
1043
1044 } // if (fbMonitorWfmMode)
1045
1046 if (fbMonitorFitMode) {
1047
1048 fcPronyFit = new TCanvas("cPronyFit", "Prony wfm fitting", w, h);
1049 fcPronyFit->Divide(2);
1050
1051 fcPronyFit->cd(1);
1052 for (UInt_t uChan = 0; uChan < kuNbChanPsd; uChan++) {
1053 fvhFitHarmonic1Chan[uChan]->Draw("same");
1054 fvhFitHarmonic2Chan[uChan]->Draw("same");
1055 }
1056
1057 fcPronyFit->cd(2);
1058 for (UInt_t uChan = 0; uChan < kuNbChanPsd; uChan++) {
1059 fvhFitQaChan[uChan]->Draw("same");
1060 }
1061
1062 AddCanvasToVector(fcPronyFit, "PronyFit");
1063
1064 /*******************************************************************/
1065
1066 } // if (fbMonitorFitMode)
1067
1068 } // if (fbMonitorChanMode)
1069
1070 return kTRUE;
1071}
1072
1087
1089{
1090
1091 if (fbMonitorChanMode) {
1092 for (UInt_t uChan = 0; uChan < kuNbChanPsd; ++uChan) {
1093 fvhHitZLChanEvo[uChan]->Reset();
1094 fvhHitLPChanEvo[uChan]->Reset();
1095 fvhHitFAChanEvo[uChan]->Reset();
1096 fvhHitChargeChan[uChan]->Reset();
1097 fvhHitZeroLevelChan[uChan]->Reset();
1098 fvhHitAmplChan[uChan]->Reset();
1099 fvhHitChargeByWfmChan[uChan]->Reset();
1100 if (fbMonitorWfmMode) fvhHitWfmChan[uChan]->Reset();
1101
1102 if (fbMonitorFitMode) {
1103 fvhHitFitWfmChan[uChan]->Reset();
1104 fvhFitHarmonic1Chan[uChan]->Reset();
1105 fvhFitHarmonic2Chan[uChan]->Reset();
1106 fvhFitQaChan[uChan]->Reset();
1107 } // if (fbMonitorFitMode)
1108 } // for( UInt_t uChan = 0; uChan < kuNbChanPsd; ++uChan )
1109 } // if(fbMonitorChanMode)
1110
1111 if (fbMonitorWfmMode) {
1112 for (UInt_t uFlatIndex = 0; uFlatIndex < kuNbChanPsd * kuNbWfmRanges * kuNbWfmExamples; ++uFlatIndex)
1113 fv3hHitWfmFlattenedChan[uFlatIndex]->Reset();
1114 for (UInt_t uWfmIndex = 0; uWfmIndex < kuNbChanPsd * kuNbWfmRanges; ++uWfmIndex)
1115 kvuWfmInRangeToChangeChan[uWfmIndex] = 0;
1116 } // if(fbMonitorWfmMode)
1117
1119 fhHitChargeMap->Reset();
1120 fhHitMapEvo->Reset();
1121 fhChanHitMapEvo->Reset();
1122
1123 fhMissedData->Reset();
1124 fhAdcTime->Reset();
1125 fhMsLengthEvo->Reset();
1126
1127 fhMsgsCntPerMsEvo->Reset();
1128 fhReadMsgsCntPerMsEvo->Reset();
1129 fhLostMsgsCntPerMsEvo->Reset();
1130 fhReadEvtsCntPerMsEvo->Reset();
1131
1132 if (kTRUE == bResetTime) {
1134 fdStartTime = -1.0;
1135 } // if( kTRUE == bResetTime )
1136
1137
1138 return kTRUE;
1139}
1140
1141// -------------------------------------------------------------------------
std::vector< double > GenerateLogBinArray(uint32_t uNbDecadesLog, uint32_t uNbStepsDecade, uint32_t uNbSubStepsInStep, uint32_t &uNbBinsLog, int32_t iStartExp, bool bAddZeroStart)
std::string FormatMsHeaderPrintout(const fles::MicrosliceDescriptor &msDescriptor)
UInt_t fuNrOfChannelsPerGdpb
Number of channels in each FEE.
std::vector< UInt_t > kvuWfmInRangeToChangeChan
UInt_t fuMsIndex
Start Time in ns of previous MS from its index field in header.
std::vector< Bool_t > fvbMaskedComponents
static constexpr UInt_t GetNbChanPsd()
Double_t fdStartTime
Index of the DPB from which the MS currently unpacked is coming.
Bool_t fbFirstPackageError
Switch ON the filling of a additional set of histograms.
std::vector< TH1 * > fvhHitFitWfmChan
Waveform fitting.
UInt_t fuCurrDpbId
Current equipment ID, tells from which DPB the current MS is originating.
CbmMcbm2018PsdPar * fUnpackPar
Settings from parameter file.
UInt_t fuNrOfChannelsPerFee
Number of FEBs per GDPB.
std::vector< TH1 * > fvhHitChargeByWfmChan
Bool_t ProcessTs(const fles::Timeslice &ts)
UInt_t fuCurrentEquipmentId
Index of current MS within the TS.
UInt_t fuCurrDpbIdx
Temp holder until Current equipment ID is properly filled in MS.
std::vector< TH2 * > fvhHitZLChanEvo
Channel rate plots.
UInt_t fuHistoryHistoSize
Histograms related variables.
Bool_t ProcessMs(const fles::Timeslice &ts, size_t uMsCompIdx, size_t uMsIdx)
Double_t fdPrevMsTime
Start Time in ns of current MS from its index field in header.
Bool_t fbMonitorFitMode
Switch ON the filling waveforms histograms.
Bool_t fbMonitorChanMode
Switch ON the filling of a minimal set of histograms.
std::vector< TCanvas * > fvcWfmsChan
UInt_t fuNrOfFeePerGdpb
gDPB ID to index map
Double_t fdMsTime
Time in ns of current TS from the index of the first MS first component.
void AddMsComponentToList(size_t component, UShort_t usDetectorId)
std::map< UInt_t, UInt_t > fGdpbIdIndexMap
Total number of GDPBs in the system.
Bool_t fbMonitorWfmMode
Switch ON the filling channelwise histograms.
Bool_t ResetHistograms(Bool_t bResetTime=kTRUE)
std::vector< TH1 * > fv3hHitWfmFlattenedChan
Int_t GetGdpbId(Int_t i)
void AddHistoToVector(TNamed *pointer, std::string sFolder="")
std::vector< size_t > fvMsComponentsList
void AddCanvasToVector(TCanvas *pointer, std::string sFolder="")
Data class with information on a STS local track.
std::vector< PsdHitData > VectHitData
std::vector< PsdHitHeader > VectHitHdr
std::vector< struct PsdPackHeader > VectPackHdr
std::vector< struct PsdHitData > VectHitData
std::vector< struct PsdHitHeader > VectHitHdr
void SetSignalBegin(int SignalBeg)
void SetExternalHarmonics(std::complex< float > z1, std::complex< float > z2)
std::vector< uint16_t > GetFitWfm()
Definition PronyFitter.h:74
float GetRSquare(int gate_beg, int gate_end)
std::complex< float > * GetHarmonics()
int ChooseBestSignalBeginHarmonics(int first_sample, int last_sample)
void SetDebugMode(bool IsDebug)
Definition PronyFitter.h:50
void SetWaveform(std::vector< uint16_t > &uWfm, uint16_t uZeroLevel)
int ChooseBestSignalBegin(int first_sample, int last_sample)
float GetIntegral(int gate_beg, int gate_end)
uint64_t ulMicroSlice
Total number of hits.