Project
Loading...
Searching...
No Matches
FV0CalibCollector.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
13#include "Framework/Logger.h"
14#include <cassert>
15#include <iostream>
16#include <sstream>
17#include <numeric>
18#include <TStopwatch.h>
19
20namespace o2
21{
22namespace fv0
23{
24
26
27//_____________________________________________
28void FV0CalibInfoSlot::fill(const gsl::span<const o2::fv0::FV0CalibrationInfoObject> data)
29{
30 // fill container
31 // we first order the data that arrived, to improve speed when filling
32 int nd = data.size();
33 LOG(info) << "FV0CalibInfoSlot::fill entries in incoming data = " << nd;
34 std::vector<int> ord(nd);
35 std::iota(ord.begin(), ord.end(), 0);
36 std::sort(ord.begin(), ord.end(), [&data](int i, int j) { return data[i].getChannelIndex() < data[j].getChannelIndex(); });
37 int chPrev = 0, offsPrev = 0;
38 for (int i = 0; i < nd; i++) {
39 if (std::abs(data[ord[i]].getTime()) > HISTO_RANGE) {
40 continue;
41 }
42 const auto& dti = data[ord[i]];
43 auto ch = dti.getChannelIndex();
44 auto offset = offsPrev;
45 if (ch > chPrev) {
46 offset += std::accumulate(mEntriesSlot.begin() + chPrev, mEntriesSlot.begin() + ch, 0);
47 }
48 offsPrev = offset;
49 chPrev = ch;
50 auto it = mFV0CollectedCalibInfoSlot.emplace(mFV0CollectedCalibInfoSlot.begin() + offset, data[ord[i]].getChannelIndex(), data[ord[i]].getTime(), data[ord[i]].getCharge(), data[ord[i]].getTimeStamp());
51 mEntriesSlot[ch]++;
52 }
53}
54//_____________________________________________
56{
57 // merge data of 2 slots
58
59 LOG(info) << "Merging two slots with entries: current slot -> " << mFV0CollectedCalibInfoSlot.size() << " , previous slot -> " << prev->mFV0CollectedCalibInfoSlot.size();
60
61 int offset = 0, offsetPrev = 0;
62 std::vector<o2::fv0::FV0CalibrationInfoObject> tmpVector;
63 for (int ch = 0; ch < NCHANNELS; ch++) {
64 if (mEntriesSlot[ch] != 0) {
65 for (int i = offset; i < offset + mEntriesSlot[ch]; i++) {
66 tmpVector.emplace_back(mFV0CollectedCalibInfoSlot[i]);
67 }
68 offset += mEntriesSlot[ch];
69 }
70 if (prev->mEntriesSlot[ch] != 0) {
71 for (int i = offsetPrev; i < offsetPrev + prev->mEntriesSlot[ch]; i++) {
72 tmpVector.emplace_back(prev->mFV0CollectedCalibInfoSlot[i]);
73 }
74 offsetPrev += prev->mEntriesSlot[ch];
75 mEntriesSlot[ch] += prev->mEntriesSlot[ch];
76 }
77 }
78 mFV0CollectedCalibInfoSlot.swap(tmpVector);
79 LOG(debug) << "After merging the size is " << mFV0CollectedCalibInfoSlot.size();
80 return;
81}
82//_____________________________________________
84{
85 // to print number of entries in the tree and the channel with the max number of entries
86
87 LOG(info) << "Total number of entries " << mFV0CollectedCalibInfoSlot.size();
88 auto maxElementIndex = std::max_element(mEntriesSlot.begin(), mEntriesSlot.end());
89 auto channelIndex = std::distance(mEntriesSlot.begin(), maxElementIndex);
90 LOG(info) << "The maximum number of entries per channel in the current mFV0CollectedCalibInfo is " << *maxElementIndex << " for channel " << channelIndex;
91 return;
92}
93
94//_____________________________________________
96{
97 // to print number of entries in the tree and per channel
98
99 LOG(debug) << "Total number of entries " << mFV0CollectedCalibInfoSlot.size();
100 for (int i = 0; i < mEntriesSlot.size(); ++i) {
101 if (mEntriesSlot[i] != 0) {
102 LOG(info) << "channel " << i << " has " << mEntriesSlot[i] << " entries";
103 }
104 }
105 return;
106}
107
108//===================================================================
109
110//_____________________________________________
112{
113 // emptying the vectors
114
115 mFV0CollectedCalibInfo.clear();
116 for (int ch = 0; ch < NCHANNELS; ch++) {
117 mEntries[ch] = 0;
118 }
119
120 return;
121}
122
123//_____________________________________________
125{
126 // We define that we have enough data if the tree is big enough.
127 // each FV0CalibrationInfoObject is composed of two int8 and one int16 --> 32 bytes
128 // E.g. supposing that we have 500000 entries per channel --> 500 eneries per one amplitude bin
129 // we can check if we have 500000*Constants::nFv0Channels entries in the vector
130
131 if (mTest) {
132 return true;
133 }
135 LOG(info) << "we have " << c->getCollectedCalibInfoSlot().size() << " entries";
136 int maxNumberOfHits = mAbsMaxNumOfHits ? mMaxNumOfHits : mMaxNumOfHits * NCHANNELS;
137 if (mTFsendingPolicy || c->getCollectedCalibInfoSlot().size() > maxNumberOfHits) {
138 return true;
139 }
140 return false;
141}
142
143//_____________________________________________
145{
146
148 mFV0CollectedCalibInfo = c->getCollectedCalibInfoSlot();
149 LOG(info) << "vector of received with size = " << mFV0CollectedCalibInfo.size();
150 mEntries = c->getEntriesPerChannel();
151 return;
152}
153
154//_____________________________________________
155Slot& FV0CalibCollector::emplaceNewSlot(bool front, TFType tstart, TFType tend)
156{
157
158 auto& cont = getSlots();
159 auto& slot = front ? cont.emplace_front(tstart, tend) : cont.emplace_back(tstart, tend);
160 slot.setContainer(std::make_unique<FV0CalibInfoSlot>());
161 return slot;
162}
163
164} // end namespace fv0
165} // end namespace o2
int32_t i
uint32_t j
Definition RawData.h:0
uint32_t c
Definition RawData.h:2
std::ostringstream debug
const Container * getContainer() const
Definition TimeSlot.h:53
void finalizeSlot(Slot &slot) final
bool hasEnoughData(const Slot &slot) const final
Slot & emplaceNewSlot(bool front, TFType tstart, TFType tend) final
void merge(const FV0CalibInfoSlot *prev)
static constexpr int NCHANNELS
void fill(const gsl::span< const o2::fv0::FV0CalibrationInfoObject > data)
static constexpr int HISTO_RANGE
GLboolean * data
Definition glcorearb.h:298
GLintptr offset
Definition glcorearb.h:660
a couple of static helper functions to create timestamp values for CCDB queries or override obsolete ...
LOG(info)<< "Compressed in "<< sw.CpuTime()<< " s"