Project
Loading...
Searching...
No Matches
TPCFourierTransformAggregatorSpec.h
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
16
17#ifndef O2_TPCFOURIERTRANSFORMAGGREGATORSPEC_H
18#define O2_TPCFOURIERTRANSFORMAGGREGATORSPEC_H
19
20#include <vector>
21#include <fmt/format.h>
22#include "Framework/Task.h"
24#include "Framework/Logger.h"
26#include "Headers/DataHeader.h"
27#include "CCDB/CcdbApi.h"
31#include "TPCBase/CRU.h"
34
35using namespace o2::framework;
37using namespace o2::tpc;
38
39namespace o2::tpc
40{
41
43 std::array<long, 2> timestamp;
44 std::array<std::vector<float>, 2> idc1;
45};
46
48 std::unordered_map<long, TPCScalerProc> idcs;
49};
50
52{
53 public:
54 // Fourier type
56
57 TPCFourierTransformAggregatorSpec(const unsigned int nFourierCoefficientsStore, const unsigned int rangeIDC, const bool senddebug = false, const bool processSACs = false, const int inputLanes = 1)
58 : mIDCFourierTransform{IDCFType(rangeIDC, nFourierCoefficientsStore), IDCFType(rangeIDC, nFourierCoefficientsStore)}, mSendOutDebug{senddebug}, mProcessSACs{processSACs}, mInputLanes{inputLanes} {};
59
61 {
62 mDumpFFT = ic.options().get<bool>("dump-coefficients-agg");
63 mIntervalsSACs = ic.options().get<int>("intervalsSACs");
64 mLengthIDCScalerSeconds = ic.options().get<float>("tpcScalerLengthS");
65 mDisableScaler = ic.options().get<bool>("disable-scaler");
66 mEnableFFTCCDB = ic.options().get<bool>("enable-fft-CCDB");
67 int nthreads = ic.options().get<int>("nthreads");
69 resizeBuffer(mInputLanes);
70 }
71
73 {
74 const int lane = pc.inputs().get<int>("lane");
75 if (lane >= mInputLanes) {
76 LOGP(error, "Received data from lane {} which is >= than the specified number of expected lanes of {}!", lane, mInputLanes);
77 return;
78 }
79
80 const auto tsTmp = pc.inputs().get<std::vector<long>>("tsccdb");
81 if (tsTmp.front() == 0) {
82 LOGP(warning, "Received dummy data with empty timestamp");
83 return;
84 }
85 mCCDBBuffer[lane] = tsTmp;
86 if (mProcessedTimeStamp > mCCDBBuffer[lane].front()) {
87 LOGP(warning, "Already received data from a later time stamp {} then the currently received time stamp {}! (This might not be an issue)", mProcessedTimeStamp, mCCDBBuffer[lane].front());
88 } else {
89 mProcessedTimeStamp = mCCDBBuffer[lane].front();
90 }
91
92 if (!mProcessSACs) {
93 mIntervalsBuffer[lane] = pc.inputs().get<std::vector<unsigned int>>("intervals");
94 }
95
96 for (auto& ref : InputRecordWalker(pc.inputs(), mFilter[mProcessSACs])) {
97 auto const* dataHeader = o2::framework::DataRefUtils::getHeader<o2::header::DataHeader*>(ref);
98 const int side = dataHeader->subSpecification;
99 mIDCOneBuffer[lane][side].mIDCOne = pc.inputs().get<std::vector<float>>(ref);
100 LOGP(info, "Received {} 1D-IDCs for side {}", mIDCOneBuffer[lane][side].mIDCOne.size(), side);
101
102 if (mProcessSACs && mIntervalsBuffer[lane].empty()) {
103 const auto nValues = mIDCOneBuffer[lane][side].mIDCOne.size();
104 const int nIntervals = nValues / mIntervalsSACs;
105 const int nFirstInterval = nValues % mIntervalsSACs;
106 if (nFirstInterval == 0) {
107 mIntervalsBuffer[lane] = std::vector<unsigned int>(nIntervals, mIntervalsSACs);
108 } else {
109 mIntervalsBuffer[lane] = std::vector<unsigned int>(nIntervals + 1, mIntervalsSACs);
110 mIntervalsBuffer[lane].front() = nFirstInterval;
111 }
112 }
113 }
114
115 // buffer IDCs for TPC scaler
116 if (!mProcessSACs && !mDisableScaler) {
117 const long startTS = mCCDBBuffer[lane].front();
118 TPCScalerProc& scaler = mTPCScalerCont.idcs[startTS];
119 scaler.timestamp[0] = startTS;
120 scaler.timestamp[1] = mCCDBBuffer[lane].back();
121 for (auto& ref : InputRecordWalker(pc.inputs(), mFilterI0)) {
122 auto const* dataHeader = o2::framework::DataRefUtils::getHeader<o2::header::DataHeader*>(ref);
123 const int side = dataHeader->subSpecification;
124 const float idc0mean = pc.inputs().get<float>(ref);
125 LOGP(info, "Received {} IDC0 mean for side {}", idc0mean, side);
126 scaler.idc1[side] = mIDCOneBuffer[lane][side].mIDCOne;
127 auto& vecIDC = scaler.idc1[side];
128
129 // normalize 1D-IDCs to unity as it should be
130 if (vecIDC.size() > 0) {
131 const float mean = std::reduce(vecIDC.begin(), vecIDC.end()) / static_cast<float>(vecIDC.size());
132 LOGP(info, "normalizing by {}", mean);
133 if (std::abs(mean) > 0.001) {
134 std::transform(vecIDC.begin(), vecIDC.end(), vecIDC.begin(), [&mean](auto val) { return val / mean; });
135 }
136 }
137
138 // scale IDC1 with IDC0Mean
139 std::transform(scaler.idc1[side].begin(), scaler.idc1[side].end(), scaler.idc1[side].begin(), [&idc0mean](auto idc) { return idc * idc0mean; });
140 }
141 // check if A- and C-side has the same length!
142 const int lenA = scaler.idc1[0].size();
143 const int lenC = scaler.idc1[1].size();
144 if (lenA != lenC) {
145 // This should never happen
146 LOGP(warning, "Received IDCs have different length! A-side length: {} and C-side length: {}", lenA, lenC);
147 // add dummy to shorter vector
148 const int maxLen = std::max(lenA, lenC);
149 scaler.idc1[0].resize(maxLen);
150 scaler.idc1[1].resize(maxLen);
151 }
152
154 makeTPCScaler(pc.outputs(), false);
155 }
156
157 FourierCoeffSAC coeffSAC;
158 if (lane == mExpectedInputLane) {
159 const int nSides = mIDCOneBuffer[lane][Side::A].mIDCOne.empty() + mIDCOneBuffer[lane][Side::C].mIDCOne.empty();
160 // int iProcessLane = lane;
161 for (int iProcessLaneTmp = 0; iProcessLaneTmp < mInputLanes; ++iProcessLaneTmp) {
162 const int nSidesCurrLane = mIDCOneBuffer[mExpectedInputLane][Side::A].mIDCOne.empty() + mIDCOneBuffer[mExpectedInputLane][Side::C].mIDCOne.empty();
163 if (nSidesCurrLane != nSides) {
164 break;
165 }
166
167 for (int iSide = 0; iSide < SIDES; ++iSide) {
168 const Side side = (iSide == 0) ? A : C;
169 if (mIDCOneBuffer[mExpectedInputLane][side].mIDCOne.empty()) {
170 continue;
171 }
172 LOGP(info, "Processing input lane: {} for Side: {}", mExpectedInputLane, iSide);
173
174 // perform fourier transform of 1D-IDCs
175 mIDCFourierTransform[side].setIDCs(std::move(mIDCOneBuffer[mExpectedInputLane][side]), mIntervalsBuffer[mExpectedInputLane]);
176 mIDCFourierTransform[side].calcFourierCoefficients(mIntervalsBuffer[mExpectedInputLane].size());
177
178 if (!mProcessSACs) {
179 if (mEnableFFTCCDB) {
180 o2::ccdb::CcdbObjectInfo ccdbInfo(CDBTypeMap.at(((side == 0) ? CDBType::CalIDCFourierA : CDBType::CalIDCFourierC)), std::string{}, std::string{}, std::map<std::string, std::string>{}, mCCDBBuffer[mExpectedInputLane].front(), mCCDBBuffer[mExpectedInputLane].back());
181 auto imageFFT = o2::ccdb::CcdbApi::createObjectImage(&mIDCFourierTransform[side].getFourierCoefficients(), &ccdbInfo);
182 LOGP(info, "Sending object {} / {} of size {} bytes, valid for {} : {} ", ccdbInfo.getPath(), ccdbInfo.getFileName(), imageFFT->size(), ccdbInfo.getStartValidityTimestamp(), ccdbInfo.getEndValidityTimestamp());
183 pc.outputs().snapshot(Output{o2::calibration::Utils::gDataOriginCDBPayload, getDataDescriptionCCDBFourier(), 0}, *imageFFT.get());
185 }
186 } else {
187 coeffSAC.mCoeff[side] = mIDCFourierTransform[side].getFourierCoefficients();
188 }
189
190 if (mDumpFFT) {
191 LOGP(info, "dumping FT to file");
192 mIDCFourierTransform[side].dumpToFile(fmt::format("FourierAGG_{:02}_side{}.root", processing_helpers::getCurrentTF(pc), (int)side).data());
193 }
194
195 if (mSendOutDebug) {
196 sendOutput(pc.outputs(), side);
197 }
198 }
199
200 if (mProcessSACs && mEnableFFTCCDB) {
201 o2::ccdb::CcdbObjectInfo ccdbInfo(CDBTypeMap.at(CDBType::CalSACFourier), std::string{}, std::string{}, std::map<std::string, std::string>{}, mCCDBBuffer[mExpectedInputLane].front(), mCCDBBuffer[mExpectedInputLane].back());
202 auto imageFFT = o2::ccdb::CcdbApi::createObjectImage(&coeffSAC, &ccdbInfo);
203 LOGP(info, "Sending object {} / {} of size {} bytes, valid for {} : {} ", ccdbInfo.getPath(), ccdbInfo.getFileName(), imageFFT->size(), ccdbInfo.getStartValidityTimestamp(), ccdbInfo.getEndValidityTimestamp());
204 pc.outputs().snapshot(Output{o2::calibration::Utils::gDataOriginCDBPayload, getDataDescriptionCCDBFourier(), 0}, *imageFFT.get());
206 }
207 mExpectedInputLane = ++mExpectedInputLane % mInputLanes;
208 }
209 }
210 }
211
213 {
214 if (!mDisableScaler && !mProcessSACs) {
215 makeTPCScaler(ec.outputs(), true);
216 }
217 ec.services().get<ControlService>().readyToQuit(QuitRequest::Me);
218 }
219
223
224 private:
225 std::array<IDCFType, SIDES> mIDCFourierTransform{};
226 const bool mSendOutDebug{false};
227 const bool mProcessSACs{false};
228 const int mInputLanes{1};
229 bool mDumpFFT{false};
230 uint64_t mProcessedTimeStamp{0};
231 std::vector<std::vector<long>> mCCDBBuffer{};
232 std::vector<std::vector<unsigned int>> mIntervalsBuffer{};
233 std::vector<std::array<o2::tpc::IDCOne, SIDES>> mIDCOneBuffer{};
234 unsigned int mIntervalsSACs{12};
235 int mExpectedInputLane{0};
236 TPCScalerProcContainer mTPCScalerCont;
237 float mLengthIDCScalerSeconds = 300;
238 long mIDCSCalerEndTSLast = 0;
239 o2::tpc::TPCScaler mScalerLast;
240 bool mDisableScaler{false};
241 bool mEnableFFTCCDB{false};
242 int mRun{};
243 const std::array<std::vector<InputSpec>, 2> mFilter = {std::vector<InputSpec>{{"idcone", ConcreteDataTypeMatcher{o2::header::gDataOriginTPC, TPCFactorizeIDCSpec::getDataDescriptionIDC1()}, Lifetime::Sporadic}},
244 std::vector<InputSpec>{{"sacone", ConcreteDataTypeMatcher{o2::header::gDataOriginTPC, TPCFactorizeSACSpec::getDataDescriptionSAC1()}, Lifetime::Sporadic}}};
245 const std::vector<InputSpec> mFilterI0 = std::vector<InputSpec>{{"idczeromean", ConcreteDataTypeMatcher{o2::header::gDataOriginTPC, TPCFactorizeIDCSpec::getDataDescriptionIDC0Mean()}, Lifetime::Sporadic}};
246
247 void sendOutput(DataAllocator& output, const int side)
248 {
249 output.snapshot(Output{gDataOriginTPC, TPCFourierTransformAggregatorSpec::getDataDescriptionFourier()}, mIDCFourierTransform[side].getFourierCoefficients());
250 }
251
252 void resizeBuffer(const int expectedLanes)
253 {
254 mCCDBBuffer.resize(expectedLanes);
255 mIntervalsBuffer.resize(expectedLanes);
256 mIDCOneBuffer.resize(expectedLanes);
257 }
258
259 void makeTPCScaler(DataAllocator& output, const bool eos)
260 {
261 LOGP(info, "Making TPC scalers");
262 if (mTPCScalerCont.idcs.empty()) {
263 LOGP(warning, "No IDCs received for TPC scaler creation");
264 return;
265 }
266
267 // check if IDC scalers can be created - check length of continous received IDCs
268 std::vector<std::pair<long, long>> times;
269 times.reserve(mTPCScalerCont.idcs.size());
270 for (const auto& idc : mTPCScalerCont.idcs) {
271 times.emplace_back(idc.second.timestamp[0], idc.second.timestamp[1]);
272 }
273
274 // sort received times of the IDCs
275 std::sort(times.begin(), times.end());
276
277 // loop over times and make checks
278 const int checkGapp = 10;
279 // if the diff between end of data[i] and start of data[i+1] is smaller than this, the data is contigous
280 long timesDuration = (times.front().second - times.front().first);
281
282 // make check and store lastValid index in case IDC scalers can be created
283 int lastValidIdx = -1;
284 for (int i = 1; i < times.size(); ++i) {
285 // check time diff between start of current IDCs and end time of last IDCs
286 const auto deltaTime = times[i].first - times[i - 1].second;
287 // check if IDCs are contigous
288 if (deltaTime > (timesDuration / checkGapp)) {
289 // check if the gap is very large - in this case the gapp might be lost, so just write out the TPC scaler until the gap
290 if (deltaTime > (checkGapp * timesDuration)) {
291 lastValidIdx = i - 1;
292 }
293 LOGP(info, "Breaking as big gap between IDCs of {} detected", deltaTime);
294 break;
295 }
296
297 // check if time length is >= than mLengthIDCScalerSeconds
298 if ((times[i].first - times.front().first) / 1000 >= mLengthIDCScalerSeconds) {
299 lastValidIdx = i;
300 }
301 }
302
303 LOGP(info, "Creating IDC scalers with {} IDC objects", lastValidIdx);
304
305 if (eos) {
306 // in case of eos write out everything
307 lastValidIdx = times.empty() ? -1 : times.size() - 1;
308 LOGP(info, "End of stream detected: Creating IDC scalers with {} IDC objects", lastValidIdx);
309 }
310
311 // create IDC scaler in case index is valid
312 if (lastValidIdx >= 0) {
313 o2::tpc::TPCScaler scaler;
314 scaler.setIonDriftTimeMS(170);
315 scaler.setRun(mRun);
316 scaler.setStartTimeStampMS(times.front().first);
317 const auto idcIntegrationTime = 12 /*12 orbits integration interval per IDC*/ * o2::constants::lhc::LHCOrbitMUS / 1000;
318 scaler.setIntegrationTimeMS(idcIntegrationTime);
319
320 std::vector<float> idc1A;
321 std::vector<float> idc1C;
322 long idc1ASize = 0;
323 long idc1CSize = 0;
324
325 // in case already one object is stored add internal overlap
326 if (mIDCSCalerEndTSLast != 0) {
327 const int nOverlap = 500;
328 idc1ASize += nOverlap;
329 idc1CSize += nOverlap;
330 const auto& scalerALast = mScalerLast.getScalers(o2::tpc::Side::A);
331 const auto& scalerCLast = mScalerLast.getScalers(o2::tpc::Side::C);
332 if (scalerALast.size() > nOverlap) {
333 idc1A.insert(idc1A.end(), scalerALast.end() - nOverlap, scalerALast.end());
334 idc1C.insert(idc1C.end(), scalerCLast.end() - nOverlap, scalerCLast.end());
335 // adjust start time
336 scaler.setStartTimeStampMS(scaler.getStartTimeStampMS() - nOverlap * idcIntegrationTime);
337 }
338 } else {
339 // store end timestamp as start time stamp for first object for correct time stamp in CCDB
340 mIDCSCalerEndTSLast = scaler.getStartTimeStampMS();
341 }
342
343 for (int iter = 0; iter < 2; ++iter) {
344 if (iter == 1) {
345 idc1A.reserve(idc1ASize);
346 idc1C.reserve(idc1CSize);
347 }
348 for (int i = 0; i <= lastValidIdx; ++i) {
349 const auto& time = times[i];
350 const auto& idc = mTPCScalerCont.idcs[time.first];
351 if (iter == 0) {
352 idc1ASize += idc.idc1[0].size();
353 idc1CSize += idc.idc1[1].size();
354 } else {
355 idc1A.insert(idc1A.end(), idc.idc1[0].begin(), idc.idc1[0].end());
356 idc1C.insert(idc1C.end(), idc.idc1[1].begin(), idc.idc1[1].end());
357 // in case of eos check if the IDCs are contigous and add dummy values!
358 if (eos && (i < lastValidIdx)) {
359 const float deltaTime = times[i + 1].first - time.second;
360 // if delta time is too large add dummy values
361 if (deltaTime > (timesDuration / checkGapp)) {
362 int nDummyValues = deltaTime / idcIntegrationTime + 0.5;
363 // restrict dummy values
364 const int nMaxDummyValues = checkGapp * timesDuration / idcIntegrationTime;
365 if (nDummyValues > nMaxDummyValues) {
366 nDummyValues = nMaxDummyValues;
367 }
368
369 // add dummy to A
370 if (idc.idc1[0].size() > 0) {
371 float meanA = std::reduce(idc.idc1[0].begin(), idc.idc1[0].end()) / static_cast<float>(idc.idc1[0].size());
372 idc1A.insert(idc1A.end(), nDummyValues, meanA);
373 }
374
375 if (idc.idc1[1].size() > 0) {
376 // add dummy to C
377 float meanC = std::reduce(idc.idc1[1].begin(), idc.idc1[1].end()) / static_cast<float>(idc.idc1[1].size());
378 idc1C.insert(idc1C.end(), nDummyValues, meanC);
379 }
380 }
381 }
382 mTPCScalerCont.idcs.erase(time.first);
383 }
384 }
385 }
386 scaler.setScaler(idc1A, o2::tpc::Side::A);
387 scaler.setScaler(idc1C, o2::tpc::Side::C);
388
389 // store in CCDB
390 TTree tree("ccdb_object", "ccdb_object");
391 tree.Branch("TPCScaler", &scaler);
392 tree.Fill();
393
394 o2::ccdb::CcdbObjectInfo ccdbInfoIDC(CDBTypeMap.at(CDBType::CalScaler), std::string{}, std::string{}, std::map<std::string, std::string>{}, mIDCSCalerEndTSLast, scaler.getEndTimeStampMS(o2::tpc::Side::A));
395 auto imageIDC = o2::ccdb::CcdbApi::createObjectImage(&tree, &ccdbInfoIDC);
396 LOGP(info, "Sending object {} / {} of size {} bytes, valid for {} : {} ", ccdbInfoIDC.getPath(), ccdbInfoIDC.getFileName(), imageIDC->size(), ccdbInfoIDC.getStartValidityTimestamp(), ccdbInfoIDC.getEndValidityTimestamp());
399
400 // store end timestamp
401 mIDCSCalerEndTSLast = scaler.getEndTimeStampMS(o2::tpc::Side::A);
402
403 // for debugging
404 if (mDumpFFT) {
405 static int countwrite = 0;
406 scaler.dumpToFile(fmt::format("TPCScaler_snapshot_{}.root", countwrite++).data(), "ccdb_object");
407 }
408
409 // buffer current scaler object
410 mScalerLast = std::move(scaler);
411 }
412 }
413};
414DataProcessorSpec getTPCFourierTransformAggregatorSpec(const unsigned int rangeIDC, const unsigned int nFourierCoefficientsStore, const bool senddebug, const bool processSACs, const int inputLanes)
415{
416 std::vector<OutputSpec> outputSpecs;
421
422 if (senddebug) {
423 outputSpecs.emplace_back(ConcreteDataTypeMatcher{gDataOriginTPC, TPCFourierTransformAggregatorSpec::getDataDescriptionFourier()}, Lifetime::Sporadic);
424 }
425
426 std::vector<InputSpec> inputSpecs;
427 if (!processSACs) {
428 inputSpecs.emplace_back(InputSpec{"idczeromean", ConcreteDataTypeMatcher{gDataOriginTPC, TPCFactorizeIDCSpec::getDataDescriptionIDC0Mean()}, Lifetime::Sporadic});
429 inputSpecs.emplace_back(InputSpec{"idcone", ConcreteDataTypeMatcher{gDataOriginTPC, TPCFactorizeIDCSpec::getDataDescriptionIDC1()}, Lifetime::Sporadic});
430 inputSpecs.emplace_back(InputSpec{"tsccdb", gDataOriginTPC, TPCFactorizeIDCSpec::getDataDescriptionTimeStamp(), Lifetime::Sporadic});
431 inputSpecs.emplace_back(InputSpec{"intervals", gDataOriginTPC, TPCFactorizeIDCSpec::getDataDescriptionIntervals(), Lifetime::Sporadic});
432 inputSpecs.emplace_back(InputSpec{"lane", gDataOriginTPC, TPCFactorizeIDCSpec::getDataDescriptionLane(), Lifetime::Sporadic});
433 } else {
434 inputSpecs.emplace_back(InputSpec{"sacone", ConcreteDataTypeMatcher{gDataOriginTPC, TPCFactorizeSACSpec::getDataDescriptionSAC1()}, Lifetime::Sporadic});
435 inputSpecs.emplace_back(InputSpec{"tsccdb", gDataOriginTPC, TPCFactorizeSACSpec::getDataDescriptionTimeStamp(), Lifetime::Sporadic});
436 inputSpecs.emplace_back(InputSpec{"lane", gDataOriginTPC, TPCFactorizeSACSpec::getDataDescriptionLane(), Lifetime::Sporadic});
437 }
438
439 std::string processorName = "tpc-aggregator-ft";
440 if (processSACs) {
441 processorName = "tpc-aggregator-ft-sac";
442 }
443
444 return DataProcessorSpec{
445 processorName,
446 inputSpecs,
447 outputSpecs,
448 AlgorithmSpec{adaptFromTask<TPCFourierTransformAggregatorSpec>(nFourierCoefficientsStore, rangeIDC, senddebug, processSACs, inputLanes)},
449 Options{{"intervalsSACs", VariantType::Int, 11, {"Number of integration intervals which will be sampled for the fourier coefficients"}},
450 {"dump-coefficients-agg", VariantType::Bool, false, {"Dump fourier coefficients to file"}},
451 {"tpcScalerLengthS", VariantType::Float, 300.f, {"Length of the TPC scalers in seconds"}},
452 {"disable-scaler", VariantType::Bool, false, {"Disable creation of IDC scaler"}},
453 {"enable-fft-CCDB", VariantType::Bool, false, {"Enable writing of FFT coefficients to CCDB"}},
454 {"nthreads", VariantType::Int, 1, {"Number of threads which will be used during the calculation of the fourier coefficients."}}}};
455}
456
457} // namespace o2::tpc
458
459#endif
std::vector< unsigned long > times
int16_t time
Definition RawEventData.h:4
int32_t i
class for calculating the fourier coefficients from 1D-IDCs
void output(const std::map< std::string, ChannelStat > &channels)
Definition rawdump.cxx:197
Definition of the Names Generator class.
uint32_t side
Definition RawData.h:0
TPC factorization of SACs.
Definition A.h:16
static std::unique_ptr< std::vector< char > > createObjectImage(const T *obj, CcdbObjectInfo *info=nullptr)
Definition CcdbApi.h:103
long getEndValidityTimestamp() const
const std::string & getPath() const
long getStartValidityTimestamp() const
const std::string & getFileName() const
A helper class to iteratate over all parts of all input routes.
static void setNThreads(const int nThreads)
static constexpr header::DataDescription getDataDescriptionTimeStamp()
static constexpr header::DataDescription getDataDescriptionLane()
static constexpr header::DataDescription getDataDescriptionIDC0Mean()
static constexpr header::DataDescription getDataDescriptionIDC1()
static constexpr header::DataDescription getDataDescriptionIntervals()
static constexpr header::DataDescription getDataDescriptionSAC1()
static constexpr header::DataDescription getDataDescriptionLane()
static constexpr header::DataDescription getDataDescriptionTimeStamp()
void init(o2::framework::InitContext &ic) final
void endOfStream(o2::framework::EndOfStreamContext &ec) final
This is invoked whenever we have an EndOfStream event.
static constexpr header::DataDescription getDataDescriptionCCDBTPCScaler()
static constexpr header::DataDescription getDataDescriptionFourier()
static constexpr header::DataDescription getDataDescriptionCCDBFourier()
TPCFourierTransformAggregatorSpec(const unsigned int nFourierCoefficientsStore, const unsigned int rangeIDC, const bool senddebug=false, const bool processSACs=false, const int inputLanes=1)
void run(o2::framework::ProcessingContext &pc) final
float getScalers(unsigned int idx, o2::tpc::Side side) const
Definition TPCScaler.h:97
void setScaler(const std::vector< float > &values, const o2::tpc::Side side)
Definition TPCScaler.h:55
void setIonDriftTimeMS(float ionDriftTimeMS)
Definition TPCScaler.h:58
double getStartTimeStampMS() const
Definition TPCScaler.h:112
double getEndTimeStampMS(o2::tpc::Side side) const
Definition TPCScaler.h:115
void dumpToFile(const char *file, const char *name)
Definition TPCScaler.cxx:25
void setStartTimeStampMS(double timeStampMS)
Definition TPCScaler.h:67
void setRun(int run)
Definition TPCScaler.h:61
void setIntegrationTimeMS(float integrationTimeMS)
Definition TPCScaler.h:70
GLsizeiptr size
Definition glcorearb.h:659
GLboolean * data
Definition glcorearb.h:298
GLuint GLfloat * val
Definition glcorearb.h:1582
GLint ref
Definition glcorearb.h:291
constexpr o2::header::DataOrigin gDataOriginTPC
Definition DataHeader.h:576
constexpr double LHCOrbitMUS
Defining ITS Vertex explicitly as messageable.
Definition Cartesian.h:288
std::vector< ConfigParamSpec > Options
uint32_t getCurrentTF(o2::framework::ProcessingContext &pc)
uint64_t getRunNumber(o2::framework::ProcessingContext &pc)
Global TPC definitions and constants.
Definition SimTraits.h:168
const std::unordered_map< CDBType, const std::string > CDBTypeMap
Storage name in CCDB for each calibration and parameter type.
Definition CDBTypes.h:96
DataProcessorSpec getTPCFourierTransformAggregatorSpec(const unsigned int rangeIDC, const unsigned int nFourierCoefficientsStore, const bool senddebug, const bool processSACs, const int inputLanes)
constexpr unsigned char SIDES
Definition Defs.h:41
Side
TPC readout sidE.
Definition Defs.h:35
@ A
Definition Defs.h:35
@ C
Definition Defs.h:36
@ CalSACFourier
Fourier coefficients of CalSAC1.
@ CalScaler
Scaler from IDCs or combined estimator.
@ CalIDCFourierC
Fourier coefficients of CalIDC1.
@ CalIDCFourierA
Fourier coefficients of CalIDC1.
void empty(int)
static constexpr o2::header::DataOrigin gDataOriginCDBWrapper
Definition Utils.h:44
static constexpr o2::header::DataOrigin gDataOriginCDBPayload
Definition Utils.h:43
std::array< FourierCoeff, SIDES > mCoeff
std::unordered_map< long, TPCScalerProc > idcs
std::array< long, 2 > timestamp
start -> end timestamp
std::array< std::vector< float >, 2 > idc1
IDC1.
std::unique_ptr< TTree > tree((TTree *) flIn.Get(std::string(o2::base::NameConf::CTFTREENAME).c_str()))