Project
Loading...
Searching...
No Matches
DigitRecoSpec.cxx
Go to the documentation of this file.
1// Copyright 2019-2020 CERN and copyright holders of ALICE O2.
2// See https://alice-o2.web.cern.ch/copyright for details of the copyright holders.
3// All rights not expressly granted are reserved.
4//
5// This software is distributed under the terms of the GNU General Public
6// License v3 (GPL Version 3), copied verbatim in the file "COPYING".
7//
8// In applying this license CERN does not waive the privileges and immunities
9// granted to it by virtue of its status as an Intergovernmental Organization
10// or submit itself to any jurisdiction.
11
15
16#include <iostream>
17#include <vector>
18#include <string>
19#include <cstdlib>
22#include "Framework/Logger.h"
41
42using namespace o2::framework;
43
44namespace o2
45{
46namespace zdc
47{
48
50{
51 mTimer.Stop();
52 mTimer.Reset();
53}
54
55DigitRecoSpec::DigitRecoSpec(const int verbosity, const bool debugOut,
56 const bool enableZDCTDCCorr, const bool enableZDCEnergyParam, const bool enableZDCTowerParam, const bool enableBaselineParam)
57 : mVerbosity(verbosity), mDebugOut(debugOut), mEnableZDCTDCCorr(enableZDCTDCCorr), mEnableZDCEnergyParam(enableZDCEnergyParam), mEnableZDCTowerParam(enableZDCTowerParam), mEnableBaselineParam(enableBaselineParam)
58{
59 mTimer.Stop();
60 mTimer.Reset();
61}
62
64{
65 mMaxWave = ic.options().get<int>("max-wave");
66 if (mMaxWave > 0) {
67 LOG(warning) << "Limiting the number of waveforms in ourput to " << mMaxWave;
68 }
69 mRecoFraction = ic.options().get<double>("tf-fraction");
70 if (mRecoFraction < 0 || mRecoFraction > 1) {
71 LOG(error) << "Unphysical reconstructed fraction " << mRecoFraction << " set to 1.0";
72 mRecoFraction = 1.0;
73 }
74 if (mRecoFraction < 1) {
75 LOG(warning) << "Target fraction for reconstructed TFs = " << mRecoFraction;
76 }
77}
78
80{
81 // we call these methods just to trigger finaliseCCDB callback
82 pc.inputs().get<o2::zdc::ModuleConfig*>("moduleconfig");
83 pc.inputs().get<o2::zdc::RecoConfigZDC*>("recoconfig");
84 pc.inputs().get<o2::zdc::ZDCTDCParam*>("tdccalib");
85 if (mEnableZDCTDCCorr) {
86 pc.inputs().get<o2::zdc::ZDCTDCCorr*>("tdccorr");
87 }
88 if (mEnableZDCEnergyParam) {
89 pc.inputs().get<o2::zdc::ZDCEnergyParam*>("adccalib");
90 }
91 if (mEnableZDCTowerParam) {
92 pc.inputs().get<o2::zdc::ZDCTowerParam*>("towercalib");
93 }
94 if (mEnableBaselineParam) {
95 pc.inputs().get<o2::zdc::BaselineParam*>("basecalib");
96 }
97}
98
100{
101 if (matcher == ConcreteDataMatcher("ZDC", "MODULECONFIG", 0)) {
102 auto* config = (const o2::zdc::ModuleConfig*)obj;
103 if (mVerbosity > DbgZero) {
104 config->print();
105 }
106 mWorker.setModuleConfig(config);
107 }
108 if (matcher == ConcreteDataMatcher("ZDC", "RECOCONFIG", 0)) {
109 // Configuration parameters for ZDC reconstruction
110 auto* config = (const o2::zdc::RecoConfigZDC*)obj;
111 if (mVerbosity > DbgZero) {
112 config->print();
113 }
114 mWorker.setRecoConfigZDC(config);
115 }
116 if (matcher == ConcreteDataMatcher("ZDC", "TDCCALIB", 0)) {
117 // TDC centering
118 auto* config = (const o2::zdc::ZDCTDCParam*)obj;
119 if (mVerbosity > DbgZero) {
120 config->print();
121 }
122 mWorker.setTDCParam(config);
123 }
124 if (matcher == ConcreteDataMatcher("ZDC", "TDCCORR", 0)) {
125 // TDC correction parameters
126 auto* config = (const o2::zdc::ZDCTDCCorr*)obj;
127 if (mVerbosity > DbgZero) {
128 config->print();
129 }
130 mWorker.setTDCCorr(config);
131 }
132 if (matcher == ConcreteDataMatcher("ZDC", "ADCCALIB", 0)) {
133 // Energy calibration
134 auto* config = (const o2::zdc::ZDCEnergyParam*)obj;
135 if (mVerbosity > DbgZero) {
136 config->print();
137 }
138 mWorker.setEnergyParam(config);
139 }
140 if (matcher == ConcreteDataMatcher("ZDC", "TOWERCALIB", 0)) {
141 // Tower intercalibration
142 auto* config = (const o2::zdc::ZDCTowerParam*)obj;
143 if (mVerbosity > DbgZero) {
144 config->print();
145 }
146 mWorker.setTowerParam(config);
147 }
148 if (matcher == ConcreteDataMatcher("ZDC", "BASECALIB", 0)) {
149 // Average pedestals
150 auto* config = (const o2::zdc::BaselineParam*)obj;
151 if (mVerbosity > DbgZero) {
152 config->print();
153 }
154 mWorker.setBaselineParam(config);
155 }
156}
157
159{
160 if (!mInitialized) {
161 LOG(info) << "DigitRecoSpec::run initialization";
162 mInitialized = true;
164 if (mDebugOut) {
165 mWorker.setDebugOutput();
166 }
167 mWorker.setVerbosity(mVerbosity);
168 mWorker.init();
169 }
170 auto cput = mTimer.CpuTime();
171 mTimer.Start(false);
172
173 auto bcdata = pc.inputs().get<gsl::span<o2::zdc::BCData>>("trig");
174 auto chans = pc.inputs().get<gsl::span<o2::zdc::ChannelData>>("chan");
175 auto peds = pc.inputs().get<gsl::span<o2::zdc::OrbitData>>("peds");
176
177 // Reduce load by dropping random time frames
178 bool toProcess = true;
179 if (mRecoFraction <= 0.) {
180 toProcess = false;
181 } else if (mRecoFraction < 1.0) {
182 double frac = std::rand() / double(RAND_MAX);
183 if (frac > mRecoFraction) {
184 toProcess = false;
185 }
186 }
187
188 RecEvent recEvent;
189
190 if (toProcess) {
191 int rval = mWorker.process(peds, bcdata, chans);
192 if (rval != 0 || mWorker.inError()) {
193 LOG(warning) << bcdata.size() << " BC " << chans.size() << " CH " << peds.size() << " OD -> processing ended in ERROR @ line " << rval;
194 } else {
195 const std::vector<o2::zdc::RecEventAux>& recAux = mWorker.getReco();
196
197 // Transfer wafeform
198 bool fullinter = mWorker.getFullInterpolation();
199 if (mVerbosity > 0 || fullinter) {
200 LOGF(info, "BC processed during reconstruction %d%s", recAux.size(), (fullinter ? " FullInterpolation" : ""));
201 }
202 uint32_t nte = 0, ntt = 0, nti = 0, ntw = 0;
203 for (auto reca : recAux) {
204 bool toAddBC = true;
205 int32_t ne = reca.ezdc.size();
206 int32_t nt = 0;
207 // Store TDC hits
208 for (int32_t it = 0; it < o2::zdc::NTDCChannels; it++) {
209 for (int32_t ih = 0; ih < reca.ntdc[it]; ih++) {
210 if (toAddBC) {
211 recEvent.addBC(reca);
212 toAddBC = false;
213 }
214 nt++;
215 recEvent.addTDC(it, reca.TDCVal[it][ih], reca.TDCAmp[it][ih], reca.isBeg[it], reca.isEnd[it]);
216 }
217 }
218 // Add waveform information
219 if (fullinter) {
220 for (int32_t isig = 0; isig < o2::zdc::NChannels; isig++) {
221 if (reca.inter[isig].size() == NIS) {
222 if (toAddBC) {
223 recEvent.addBC(reca);
224 toAddBC = false;
225 }
226 recEvent.addWaveform(isig, reca.inter[isig]);
227 ntw++;
228 }
229 }
230 // Limit the number of waveforms in output message
231 if (mMaxWave > 0 && ntw >= mMaxWave) {
232 if (mVerbosity > DbgMinimal) {
233 LOG(warning) << "Maximum number of output waveforms per TF reached: " << mMaxWave;
234 }
235 break;
236 }
237 }
238 if (ne > 0) {
239 if (toAddBC) {
240 recEvent.addBC(reca);
241 toAddBC = false;
242 }
243 std::map<uint8_t, float>::iterator it;
244 for (it = reca.ezdc.begin(); it != reca.ezdc.end(); it++) {
245 recEvent.addEnergy(it->first, it->second);
246 }
247 }
248 nte += ne;
249 ntt += nt;
250 if (mVerbosity > 1 && (nt > 0 || ne > 0)) {
251 printf("Orbit %9u bc %4u ntdc %2d ne %2d channels=0x%08x\n", reca.ir.orbit, reca.ir.bc, nt, ne, reca.channels);
252 }
253 // Event information
254 nti += recEvent.addInfos(reca);
255 }
256 LOG(info) << bcdata.size() << " BC " << chans.size() << " CH " << peds.size() << " OD "
257 << "-> Reconstructed " << ntt << " signal TDCs and " << nte << " ZDC energies and "
258 << nti << " info messages in " << recEvent.mRecBC.size() << "/" << recAux.size() << " b.c. and "
259 << ntw << " waveform chunks";
260 }
261 } else {
262 LOG(info) << bcdata.size() << " BC " << chans.size() << " CH " << peds.size() << " OD "
263 << "-> SKIPPED because of requested reconstruction fraction = " << mRecoFraction;
264 }
265 // TODO: rate information for all channels
266 // TODO: summary of reconstruction to be collected by DQM?
267 pc.outputs().snapshot(Output{"ZDC", "BCREC", 0}, recEvent.mRecBC);
268 pc.outputs().snapshot(Output{"ZDC", "ENERGY", 0}, recEvent.mEnergy);
269 pc.outputs().snapshot(Output{"ZDC", "TDCDATA", 0}, recEvent.mTDCData);
270 pc.outputs().snapshot(Output{"ZDC", "INFO", 0}, recEvent.mInfo);
271 pc.outputs().snapshot(Output{"ZDC", "WAVE", 0}, recEvent.mWaveform);
272 mTimer.Stop();
273}
274
276{
277 mWorker.eor();
278 LOGF(info, "ZDC Reconstruction total timing: Cpu: %.3e Real: %.3e s in %d slots", mTimer.CpuTime(), mTimer.RealTime(), mTimer.Counter() - 1);
279}
280
281framework::DataProcessorSpec getDigitRecoSpec(const int verbosity = 0, const bool enableDebugOut = true, const bool enableZDCTDCCorr = true, const bool enableZDCEnergyParam = true, const bool enableZDCTowerParam = true, const bool enableBaselineParam = true)
282{
283 std::vector<InputSpec> inputs;
284 inputs.emplace_back("trig", "ZDC", "DIGITSBC", 0, Lifetime::Timeframe);
285 inputs.emplace_back("chan", "ZDC", "DIGITSCH", 0, Lifetime::Timeframe);
286 inputs.emplace_back("peds", "ZDC", "DIGITSPD", 0, Lifetime::Timeframe);
287 inputs.emplace_back("moduleconfig", "ZDC", "MODULECONFIG", 0, Lifetime::Condition, o2::framework::ccdbParamSpec(o2::zdc::CCDBPathConfigModule.data()));
288 inputs.emplace_back("recoconfig", "ZDC", "RECOCONFIG", 0, Lifetime::Condition, o2::framework::ccdbParamSpec(o2::zdc::CCDBPathRecoConfigZDC.data()));
289 inputs.emplace_back("tdccalib", "ZDC", "TDCCALIB", 0, Lifetime::Condition, o2::framework::ccdbParamSpec(o2::zdc::CCDBPathTDCCalib.data()));
290 if (enableZDCTDCCorr) {
291 inputs.emplace_back("tdccorr", "ZDC", "TDCCORR", 0, Lifetime::Condition, o2::framework::ccdbParamSpec(o2::zdc::CCDBPathTDCCorr.data()));
292 } else {
293 LOG(warning) << "ZDCTDCCorr has been disabled - no correction is applied";
294 }
295 if (enableZDCEnergyParam) {
296 inputs.emplace_back("adccalib", "ZDC", "ADCCALIB", 0, Lifetime::Condition, o2::framework::ccdbParamSpec(o2::zdc::CCDBPathEnergyCalib.data()));
297 } else {
298 LOG(warning) << "ZDCEnergyParam has been disabled - no energy calibration is applied";
299 }
300 if (enableZDCTowerParam) {
301 inputs.emplace_back("towercalib", "ZDC", "TOWERCALIB", 0, Lifetime::Condition, o2::framework::ccdbParamSpec(o2::zdc::CCDBPathTowerCalib.data()));
302 } else {
303 LOG(warning) << "ZDCTowerParam has been disabled - no tower intercalibration";
304 }
305 if (enableBaselineParam) {
306 inputs.emplace_back("basecalib", "ZDC", "BASECALIB", 0, Lifetime::Condition, o2::framework::ccdbParamSpec(o2::zdc::CCDBPathBaselineCalib.data()));
307 } else {
308 LOG(warning) << "BaselineParam has been disabled - no fallback in case orbit pedestals are missing";
309 }
310
311 std::vector<OutputSpec> outputs;
312 outputs.emplace_back("ZDC", "BCREC", 0, Lifetime::Timeframe);
313 outputs.emplace_back("ZDC", "ENERGY", 0, Lifetime::Timeframe);
314 outputs.emplace_back("ZDC", "TDCDATA", 0, Lifetime::Timeframe);
315 outputs.emplace_back("ZDC", "INFO", 0, Lifetime::Timeframe);
316 outputs.emplace_back("ZDC", "WAVE", 0, Lifetime::Timeframe);
317
318 return DataProcessorSpec{
319 "zdc-digi-reco",
320 inputs,
321 outputs,
322 AlgorithmSpec{adaptFromTask<DigitRecoSpec>(verbosity, enableDebugOut, enableZDCTDCCorr, enableZDCEnergyParam, enableZDCTowerParam, enableBaselineParam)},
323 o2::framework::Options{{"max-wave", o2::framework::VariantType::Int, 0, {"Maximum number of waveforms per TF in output"}},
324 {"tf-fraction", o2::framework::VariantType::Double, 1.0, {"Fraction of reconstructed TFs"}}}};
325}
326
327} // namespace zdc
328} // namespace o2
#define verbosity
Class to describe fired triggered and/or stored channels for the BC and to refer to channel data.
Baseline calibration data.
Run ZDC digits reconstruction.
Definition of the Names Generator class.
Class to describe pedestal data accumulated over the orbit.
Class to describe reconstructed ZDC event (single BC with signal in one of detectors)
ZDC reconstruction parameters.
ZDC TDC correction parameters.
Parameters to correct TDCs (produced by QA)
Container class to store NTimeBinsPerBC ADC values of single ZDC channel.
void snapshot(const Output &spec, T const &object)
ConfigParamRegistry const & options()
Definition InitContext.h:33
decltype(auto) get(R binding, int part=0) const
DataAllocator & outputs()
The data allocator is used to allocate memory for the output data.
InputRecord & inputs()
The inputs associated with this processing context.
void setDebugOutput(bool state=true)
Definition DigiReco.h:79
void setTowerParam(const ZDCTowerParam *param)
Definition DigiReco.h:100
bool getFullInterpolation()
Definition DigiReco.h:120
void setRecoConfigZDC(const RecoConfigZDC *cfg)
Definition DigiReco.h:104
void setModuleConfig(const ModuleConfig *moduleConfig)
Definition DigiReco.h:92
void setBaselineParam(const BaselineParam *param)
Definition DigiReco.h:102
const std::vector< o2::zdc::RecEventAux > & getReco()
Definition DigiReco.h:143
void setTDCParam(const ZDCTDCParam *param)
Definition DigiReco.h:94
void setVerbosity(int v)
Definition DigiReco.h:74
void setEnergyParam(const ZDCEnergyParam *param)
Definition DigiReco.h:98
void setTDCCorr(const ZDCTDCCorr *param)
Definition DigiReco.h:96
int process(const gsl::span< const o2::zdc::OrbitData > &orbitdata, const gsl::span< const o2::zdc::BCData > &bcdata, const gsl::span< const o2::zdc::ChannelData > &chdata)
Definition DigiReco.cxx:493
void endOfStream(o2::framework::EndOfStreamContext &ec) final
This is invoked whenever we have an EndOfStream event.
void init(o2::framework::InitContext &ic) final
void updateTimeDependentParams(o2::framework::ProcessingContext &pc)
void run(o2::framework::ProcessingContext &pc) final
void finaliseCCDB(o2::framework::ConcreteDataMatcher &matcher, void *obj) final
GLboolean * data
Definition glcorearb.h:298
Defining PrimaryVertex explicitly as messageable.
Definition TFIDInfo.h:20
std::vector< ConfigParamSpec > ccdbParamSpec(std::string const &path, int runDependent, std::vector< CCDBMetadata > metadata={}, int qrate=0)
std::vector< ConfigParamSpec > Options
struct o2::upgrades_utils::@463 zdc
structure to keep FT0 information
framework::DataProcessorSpec getDigitRecoSpec(const int verbosity, const bool enableDebugOut, const bool enableZDCTDCCorr, const bool enableZDCEnergyParam, const bool enableZDCTowerParam, const bool enableBaselineParam)
create a processor spec
const std::string CCDBPathTDCCorr
Definition Constants.h:223
const std::string CCDBPathBaselineCalib
Definition Constants.h:229
constexpr int NTDCChannels
Definition Constants.h:90
constexpr int NChannels
Definition Constants.h:65
const std::string CCDBPathRecoConfigZDC
Definition Constants.h:220
const std::string CCDBPathEnergyCalib
Definition Constants.h:224
const std::string CCDBPathTDCCalib
Definition Constants.h:221
constexpr int DbgMinimal
Definition Constants.h:208
constexpr int DbgZero
Definition Constants.h:207
constexpr int NIS
Definition Constants.h:97
const std::string CCDBPathTowerCalib
Definition Constants.h:225
const std::string CCDBPathConfigModule
Definition Constants.h:219
a couple of static helper functions to create timestamp values for CCDB queries or override obsolete ...
void print(bool printall=true) const
uint32_t addInfos(const RecEventAux &reca)
Definition RecEvent.cxx:94
void addTDC(uint8_t ch, int16_t val, int16_t amp, bool isbeg=false, bool isend=false)
Definition RecEvent.h:66
std::vector< o2::zdc::ZDCWaveform > mWaveform
Event quality information.
Definition RecEvent.h:42
void addEnergy(uint8_t ch, float energy)
Definition RecEvent.h:56
void addWaveform(uint8_t ch, std::vector< float > &wave)
Definition RecEvent.h:113
std::vector< o2::zdc::ZDCEnergy > mEnergy
Interaction record and references to data.
Definition RecEvent.h:39
std::vector< uint16_t > mInfo
ZDC TDC.
Definition RecEvent.h:41
void addBC(const RecEventAux &reca)
ZDC waveform.
Definition RecEvent.h:45
std::vector< o2::zdc::BCRecData > mRecBC
Definition RecEvent.h:38
std::vector< o2::zdc::ZDCTDCData > mTDCData
ZDC energy.
Definition RecEvent.h:40
LOG(info)<< "Compressed in "<< sw.CpuTime()<< " s"