Project
Loading...
Searching...
No Matches
TOFIntegrateClusterSpec.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
16
18#include "Framework/Task.h"
22#include "TOFBase/Geo.h"
24
25#include <algorithm>
26#include <numeric>
27
28using namespace o2::framework;
29
30namespace o2
31{
32namespace tof
33{
34
36{
37 public:
39 TOFIntegrateClusters(std::shared_ptr<o2::base::GRPGeomRequest> req, const bool disableWriter) : mCCDBRequest(req), mDisableWriter(disableWriter){};
40
42 {
44 mNSlicesTF = ic.options().get<int>("nSlicesTF");
45 mTimeCutNoisePS = ic.options().get<int>("time-cut-noise-PS");
46 mMinClustersNeighbours = ic.options().get<int>("min-clusters-neighbours");
47 mTagNoise = !(ic.options().get<bool>("do-not-tag-noise"));
48 }
49
50 void run(ProcessingContext& pc) final
51 {
53 const int nHBFPerTF = o2::base::GRPGeomHelper::instance().getNHBFPerTF();
55 const double maxTimeTF = orbitLength * nHBFPerTF; // maximum time if a TF in PS
56 const double sliceWidthPS = maxTimeTF / mNSlicesTF; // integration time
57 const double sliceWidthMS = sliceWidthPS * 1.E-9;
58 const int nSlices = maxTimeTF / sliceWidthPS; // number of slices a TF is divided into
59 const double sliceWidthPSinv = 1. / sliceWidthPS;
60 const float sliceWidthMSinv = 1. / float(sliceWidthMS);
61
62 // storage for integrated currents
63 o2::pmr::vector<float> iTOFCNCl(nSlices);
64 o2::pmr::vector<float> iTOFCqTot(nSlices);
65
66 const auto clusters = pc.inputs().get<gsl::span<o2::tof::Cluster>>("tofcluster");
67 if (mTagNoise) {
68 tagNoiseCluster(clusters);
69 }
70
71 for (size_t iCl = 0; iCl < clusters.size(); ++iCl) {
72 if (mTagNoise && (mCounterNeighbours[iCl] < mMinClustersNeighbours)) {
73 continue;
74 }
75
76 const double timePS = clusters[iCl].getTime();
77 const unsigned int sliceInTF = timePS * sliceWidthPSinv;
78 if (sliceInTF < mNSlicesTF) {
79 ++iTOFCNCl[sliceInTF]; // count number of cluster
80 iTOFCqTot[sliceInTF] += clusters[iCl].getTot(); // integrated charge
81 } else {
82 LOGP(info, "slice in TF of ICC {} is larger than max slice {} with nTSPerSlice {}", sliceInTF, mNSlicesTF, nSlices);
83 }
84 }
85
86 // normalize currents to integration time
87 std::transform(iTOFCNCl.begin(), iTOFCNCl.end(), iTOFCNCl.begin(), [sliceWidthMSinv](float& val) { return val * sliceWidthMSinv; });
88 std::transform(iTOFCqTot.begin(), iTOFCqTot.end(), iTOFCqTot.begin(), [sliceWidthMSinv](float& val) { return val * sliceWidthMSinv; });
89
90 sendOutput(pc, std::move(iTOFCNCl), std::move(iTOFCqTot));
91 }
92
94 {
95 LOGP(info, "Finalizing calibration");
96 }
97
99
100 private:
101 int mNSlicesTF = 11;
102 const bool mDisableWriter{false};
103 int mTimeCutNoisePS = 4000;
104 int mMinClustersNeighbours = 2;
105 bool mTagNoise{true};
106 std::vector<int> mCounterNeighbours;
107 std::shared_ptr<o2::base::GRPGeomRequest> mCCDBRequest;
108
109 void sendOutput(ProcessingContext& pc, o2::pmr::vector<float> iTOFCNCl, o2::pmr::vector<float> iTOFCqTot)
110 {
111 pc.outputs().adoptContainer(Output{header::gDataOriginTOF, "ITOFCN"}, std::move(iTOFCNCl));
112 pc.outputs().adoptContainer(Output{header::gDataOriginTOF, "ITOFCQ"}, std::move(iTOFCqTot));
113 // in case of ROOT output also store the TFinfo in the TTree
114 if (!mDisableWriter) {
117 pc.outputs().snapshot(Output{header::gDataOriginTOF, "ITOFTFID"}, tfinfo);
118 }
119 }
120
126 void tagNoiseCluster(const gsl::span<const o2::tof::Cluster>& tofArray)
127 {
128 const size_t nCl = tofArray.size();
129 // reset buffer
130 std::fill(mCounterNeighbours.begin(), mCounterNeighbours.end(), 0);
131 mCounterNeighbours.resize(nCl);
132
133 std::vector<size_t> tofTimeIndex(nCl);
134 std::iota(tofTimeIndex.begin(), tofTimeIndex.end(), 0);
135
136 // 0.) sort clusters in time: filling index to tofTimeIndex
137 std::stable_sort(tofTimeIndex.begin(), tofTimeIndex.end(), [&tofArray](size_t i1, size_t i2) { return tofArray[i1].getTime() < tofArray[i2].getTime(); });
138
139 // 1.) update counter - count timebin after the trigger bin in interval (-mTimeCutNoisePS,mTimeCutNoisePS) ~ (-4ns,4ns)
140 // Loop over all cluster and count for each cluster the neighouring clusters -defined by timeCut- in time
141 for (int i0 = 0; i0 < nCl; i0++) {
142 const double t0 = tofArray[tofTimeIndex[i0]].getTime();
143 // loop over clusters after current cluster which are close in time and count them
144 for (int idP = i0 + 1; idP < nCl; ++idP) {
145 if (std::abs(t0 - tofArray[tofTimeIndex[idP]].getTime()) > mTimeCutNoisePS) {
146 break;
147 }
148 // count clusters in local neighbourhood
149 ++mCounterNeighbours[tofTimeIndex[i0]];
150 }
151 // loop over clusters before current cluster which are close in time and count them
152 for (int idM = i0 - 1; idM >= 0; --idM) {
153 if (std::abs(t0 - tofArray[tofTimeIndex[idM]].getTime()) > mTimeCutNoisePS) {
154 break;
155 }
156 // count clusters in local neighbourhood
157 ++mCounterNeighbours[tofTimeIndex[i0]];
158 }
159 }
160 }
161};
162
164{
165 std::vector<InputSpec> inputs;
166 inputs.emplace_back("tofcluster", o2::header::gDataOriginTOF, "CLUSTERS", 0, Lifetime::Timeframe);
167
168 auto ccdbRequest = std::make_shared<o2::base::GRPGeomRequest>(false, // orbitResetTime
169 true, // GRPECS=true for nHBF per TF
170 false, // GRPLHCIF
171 false, // GRPMagField
172 false, // askMatLUT
174 inputs);
175
176 std::vector<OutputSpec> outputs;
177 outputs.emplace_back(o2::header::gDataOriginTOF, "ITOFCN", 0, Lifetime::Timeframe);
178 outputs.emplace_back(o2::header::gDataOriginTOF, "ITOFCQ", 0, Lifetime::Timeframe);
179 if (!disableWriter) {
180 outputs.emplace_back(o2::header::gDataOriginTOF, "ITOFTFID", 0, Lifetime::Timeframe);
181 }
182
183 return DataProcessorSpec{
184 "TOFIntegrateClusters",
185 inputs,
186 outputs,
187 AlgorithmSpec{adaptFromTask<TOFIntegrateClusters>(ccdbRequest, disableWriter)},
188 Options{{"nSlicesTF", VariantType::Int, 11, {"number of slices into which a TF is divided"}},
189 {"time-cut-noise-PS", VariantType::Int, 4000, {"time in ps for which neighbouring clusters are checked during the noise filtering"}},
190 {"min-clusters-neighbours", VariantType::Int, 2, {"minimum neighbouring clusters for the noise filtering"}},
191 {"do-not-tag-noise", VariantType::Bool, false, {"Do not reject noise"}}}};
192}
193
194} // end namespace tof
195} // end namespace o2
Definition of the TOF cluster.
Helper for geometry and GRP related CCDB requests.
void checkUpdates(o2::framework::ProcessingContext &pc)
bool finaliseCCDB(o2::framework::ConcreteDataMatcher &matcher, void *obj)
static GRPGeomHelper & instance()
void setRequest(std::shared_ptr< GRPGeomRequest > req)
void snapshot(const Output &spec, T const &object)
CacheId adoptContainer(const Output &, ContainerT &, CacheStrategy, o2::header::SerializationMethod)
DataAllocator & outputs()
The data allocator is used to allocate memory for the output data.
static constexpr Double_t BC_TIME_INPS
Definition Geo.h:103
void endOfStream(EndOfStreamContext &eos) final
This is invoked whenever we have an EndOfStream event.
void finaliseCCDB(o2::framework::ConcreteDataMatcher &matcher, void *obj) final
void run(ProcessingContext &pc) final
TOFIntegrateClusters(std::shared_ptr< o2::base::GRPGeomRequest > req, const bool disableWriter)
\constructor
void init(framework::InitContext &ic) final
GLuint GLfloat * val
Definition glcorearb.h:1582
GLuint GLfloat GLfloat GLfloat GLfloat GLfloat GLfloat GLfloat t0
Definition glcorearb.h:5034
constexpr o2::header::DataOrigin gDataOriginTOF
Definition DataHeader.h:575
constexpr int LHCMaxBunches
Defining PrimaryVertex explicitly as messageable.
Definition TFIDInfo.h:20
std::vector< ConfigParamSpec > Options
std::vector< T, o2::pmr::polymorphic_allocator< T > > vector
o2::framework::DataProcessorSpec getTOFIntegrateClusterSpec(const bool disableWriter)
a couple of static helper functions to create timestamp values for CCDB queries or override obsolete ...
static void fillTFIDInfo(o2::framework::ProcessingContext &pc, o2::dataformats::TFIDInfo &ti)
std::vector< Cluster > clusters