Project
Loading...
Searching...
No Matches
EMCALTimeCalibData.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"
15#include "CCDB/CcdbApi.h"
17#include <boost/histogram.hpp>
18#include <boost/histogram/ostream.hpp>
19#include <boost/histogram/algorithm/sum.hpp>
20#include <boost/format.hpp>
21#include <cassert>
22#include <iostream>
23#include <sstream>
24#include <TStopwatch.h>
25
26namespace o2
27{
28namespace emcal
29{
30
33using boost::histogram::indexed;
34
35//===================================================================
36//_____________________________________________
37void EMCALTimeCalibData::PrintStream(std::ostream& stream) const
38{
39 stream << "EMCAL Cell ID: " << mTimeHisto[0] << "\n";
40}
41//_____________________________________________
43{
44 LOG(debug) << *this;
45}
46//_____________________________________________
47std::ostream& operator<<(std::ostream& stream, const EMCALTimeCalibData& emcdata)
48{
49 emcdata.PrintStream(stream);
50 return stream;
51}
52//_____________________________________________
54{
55 mEvents += prev->getNEvents();
56 mNEntriesInHisto += prev->getNEntriesInHisto();
57 mTimeHisto[0] += prev->getHisto();
58}
59//_____________________________________________
61{
62 bool enough = false;
63
64 LOG(debug) << "mNEntriesInHisto: " << mNEntriesInHisto << " needed: " << EMCALCalibParams::Instance().minNEntries_tc << " mEvents = " << mEvents;
65 // use enrties in histogram for calibration
66 if (!EMCALCalibParams::Instance().useNEventsForCalib_tc && mNEntriesInHisto > EMCALCalibParams::Instance().minNEntries_tc) {
67 enough = true;
68 }
69 // use number of events (from emcal trigger record) for calibration
70 if (EMCALCalibParams::Instance().useNEventsForCalib_tc && mEvents > EMCALCalibParams::Instance().minNEvents_tc) {
71 enough = true;
72 }
73
74 return enough;
75}
76
77//_____________________________________________
78void EMCALTimeCalibData::fill(const gsl::span<const o2::emcal::Cell> data)
79{
80 // the fill function is called once per event
81 mEvents++;
82
83 if (data.size() == 0) {
84 return;
85 }
86 auto fillfunction = [this](int thread, const gsl::span<const o2::emcal::Cell> data, double minCellEnergy) {
87 LOG(debug) << "filling in thread " << thread << " ncells = " << data.size();
88 auto& mCurrentHist = mTimeHisto[thread];
89 unsigned int nEntries = 0;
90 for (auto cell : data) {
91 double cellEnergy = cell.getEnergy();
92 double cellTime = cell.getTimeStamp();
93 int id = cell.getTower();
94 if (mApplyGainCalib) {
95 LOG(debug) << " gain calib factor for cell " << id << " = " << mGainCalibFactors->getGainCalibFactors(id);
96 cellEnergy *= mGainCalibFactors->getGainCalibFactors(id);
97 }
98 if (cellEnergy > minCellEnergy) {
99 LOG(debug) << "inserting in cell ID " << id << ": cellTime = " << cellTime;
100 mCurrentHist(cellTime, id);
101 nEntries++;
102 }
103 }
104 mVecNEntriesInHisto[thread] += nEntries;
105 };
106
107 std::vector<gsl::span<const o2::emcal::Cell>> ranges(mNThreads);
108 auto size_per_thread = static_cast<unsigned int>(std::ceil((static_cast<float>(data.size()) / mNThreads)));
109 unsigned int currentfirst = 0;
110 for (int ithread = 0; ithread < mNThreads; ithread++) {
111 unsigned int nelements = std::min(size_per_thread, static_cast<unsigned int>(data.size() - 1 - currentfirst));
112 ranges[ithread] = data.subspan(currentfirst, nelements);
113 currentfirst += nelements;
114 LOG(debug) << "currentfirst " << currentfirst << " nelements " << nelements;
115 }
116
117 double minCellEnergy = EMCALCalibParams::Instance().minCellEnergy_tc;
118
119#if (defined(WITH_OPENMP) && !defined(__CLING__))
120 LOG(debug) << "Number of threads that will be used = " << mNThreads;
121#pragma omp parallel for num_threads(mNThreads)
122#else
123 LOG(debug) << "OPEN MP will not be used for the bad channel calibration";
124#endif
125 for (int ithread = 0; ithread < mNThreads; ithread++) {
126 fillfunction(ithread, ranges[ithread], minCellEnergy);
127 }
128
129 // only sum up entries if needed
130 if (!o2::emcal::EMCALCalibParams::Instance().useNEventsForCalib_tc) {
131 for (auto& nEntr : mVecNEntriesInHisto) {
132 mNEntriesInHisto += nEntr;
133 nEntr = 0;
134 }
135 }
136}
137
138//_____________________________________________
140{
141 o2::emcal::TimeCalibrationParams TimeCalibParams;
142 return TimeCalibParams;
143}
144
145} // end namespace emcal
146} // end namespace o2
Utils and constants for calibration and related workflows.
std::ostringstream debug
long unsigned int getNEntriesInHisto() const
Get the number of entries in histogram.
void print()
Print a useful message about the container.
void PrintStream(std::ostream &stream) const
print stream
int getNEvents() const
Get number of events currently available for calibration.
void fill(const gsl::span< const o2::emcal::Cell > data)
Fill the container with the cell ID and amplitude and time information.
void merge(EMCALTimeCalibData *prev)
Merge the data of two slots.
const boostHisto & getHisto()
Get current histogram.
bool hasEnoughData() const
Check if enough data for calibration has been accumulated.
o2::emcal::TimeCalibrationParams process()
Actual function where calibration is done. Has to be called in has enough data when enough data is th...
float getGainCalibFactors(unsigned short iCell) const
Get the gain calibration factor for a certain cell.
GLboolean * data
Definition glcorearb.h:298
GLuint GLuint stream
Definition glcorearb.h:1806
std::ostream & operator<<(std::ostream &stream, const Cell &cell)
Stream operator for EMCAL cell.
Definition Cell.cxx:355
a couple of static helper functions to create timestamp values for CCDB queries or override obsolete ...
LOG(info)<< "Compressed in "<< sw.CpuTime()<< " s"