Project
Loading...
Searching...
No Matches
TPCResidualReaderSpec.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 <memory>
17#include <vector>
18#include <boost/algorithm/string/predicate.hpp>
19#include "TFile.h"
20#include "TTree.h"
21#include "TGrid.h"
22#include "Framework/Task.h"
25#include "DataFormatsTPC/Defs.h"
33
34using namespace o2::framework;
35
36namespace o2
37{
38namespace tpc
39{
40
41class TPCResidualReader : public Task
42{
43 public:
44 TPCResidualReader(bool doBinning, GID::mask_t src) : mDoBinning(doBinning), mSources(src) {}
45 ~TPCResidualReader() override = default;
46 void init(InitContext& ic) final;
47 void run(ProcessingContext& pc) final;
48
49 private:
50 void connectTree(const std::string& filename);
51
52 bool revalidateTrack(const TrackData& trk) const;
53
56 void fillResiduals(const int iSec);
57
58 std::unique_ptr<TFile> mFile;
59 std::unique_ptr<TTree> mTreeResiduals;
60 std::unique_ptr<TTree> mTreeStats;
61 std::unique_ptr<TTree> mTreeUnbinnedResiduals;
62 std::unique_ptr<TTree> mTreeTrackData;
64 bool mDoBinning{false};
65 bool mStoreBinnedResiduals{false};
66 GID::mask_t mSources;
67 TrackResiduals mTrackResiduals;
68 std::vector<std::string> mFileNames;
69 std::string mOutfile{"debugVoxRes.root"};
70 std::vector<TrackResiduals::LocalResid> mResiduals, *mResidualsPtr = &mResiduals;
71 std::array<std::vector<TrackResiduals::LocalResid>, SECTORSPERSIDE * SIDES> mResidualsSector;
72 std::array<std::vector<TrackResiduals::LocalResid>*, SECTORSPERSIDE * SIDES> mResidualsSectorPtr;
73 std::array<std::vector<TrackResiduals::VoxStats>, SECTORSPERSIDE * SIDES> mVoxStatsSector;
74 std::vector<UnbinnedResid> mUnbinnedResiduals, *mUnbinnedResidualsPtr = &mUnbinnedResiduals;
75 std::vector<TrackDataCompact> mTrackData, *mTrackDataPtr = &mTrackData;
76 std::vector<TrackData> mTrackDataExtra, *mTrackDataExtraPtr = &mTrackDataExtra;
77};
78
79void TPCResidualReader::fillResiduals(const int iSec)
80{
81 auto brStats = mTreeStats->GetBranch(Form("sec%d", iSec));
82 brStats->SetAddress(mTrackResiduals.getVoxStatPtr());
83 brStats->GetEntry(mTreeStats->GetEntries() - 1); // only the last entry is of interest
84 mTrackResiduals.fillStats(iSec);
85
86 // in case autosave was enabled, we have multiple entries in the statistics tree
87 auto brResid = mTreeResiduals->GetBranch(Form("sec%d", iSec));
88 brResid->SetAddress(&mResidualsPtr);
89
90 for (int iEntry = 0; iEntry < brResid->GetEntries(); ++iEntry) {
91 brResid->GetEntry(iEntry);
92 LOGF(debug, "Pushing %lu TPC residuals at entry %i for sector %i", mResiduals.size(), iEntry, iSec);
93 for (const auto& res : mResiduals) {
94 LOGF(debug, "Adding residual from Voxel %i-%i-%i. dy(%i), dz(%i), tg(%i)", res.bvox[0], res.bvox[1], res.bvox[2], res.dy, res.dz, res.tgSlp);
95 }
96 mTrackResiduals.getLocalResVec().insert(mTrackResiduals.getLocalResVec().end(), mResiduals.begin(), mResiduals.end());
97 }
98}
99
101{
102 const auto dontCheckFileAccess = ic.options().get<bool>("dont-check-file-access");
103
104 auto fileList = o2::RangeTokenizer::tokenize<std::string>(ic.options().get<std::string>("residuals-infiles"));
105 mOutfile = ic.options().get<std::string>("outfile");
106 mTrackResiduals.init();
107
108 // check if only one input file (a txt file contaning a list of files is provided)
109 if (fileList.size() == 1) {
110 if (boost::algorithm::ends_with(fileList.front(), "txt")) {
111 LOGP(info, "Reading files from input file {}", fileList.front());
112 std::ifstream is(fileList.front());
113 std::istream_iterator<std::string> start(is);
114 std::istream_iterator<std::string> end;
115 std::vector<std::string> fileNamesTmp(start, end);
116 fileList = fileNamesTmp;
117 }
118 }
119
120 const std::string inpDir = o2::utils::Str::rectifyDirectory(ic.options().get<std::string>("input-dir"));
121 for (const auto& file : fileList) {
122 if ((file.find("alien://") == 0) && !gGrid && !TGrid::Connect("alien://")) {
123 LOG(fatal) << "Failed to open alien connection";
124 }
125 const auto fileDir = o2::utils::Str::concat_string(inpDir, file);
126 if (!dontCheckFileAccess) {
127 std::unique_ptr<TFile> filePtr(TFile::Open(fileDir.data()));
128 if (!filePtr || !filePtr->IsOpen() || filePtr->IsZombie()) {
129 LOGP(warning, "Could not open file {}", fileDir);
130 continue;
131 }
132 }
133 mFileNames.emplace_back(fileDir);
134 }
135
136 if (mFileNames.size() == 0) {
137 LOGP(error, "No input files to process");
138 }
139
140 mStoreBinnedResiduals = ic.options().get<bool>("store-binned");
141}
142
144{
145 mTrackResiduals.createOutputFile(mOutfile.data()); // FIXME remove when map output is handled properly
146
147 if (mDoBinning) {
148 // initialize the trees for the binned residuals and the statistics
149 LOGP(info, "Preparing the binning of residuals. Storing the afterwards is set to {}", mStoreBinnedResiduals);
150 mTreeResiduals = std::make_unique<TTree>("resid", "TPC binned residuals");
151 if (!mStoreBinnedResiduals) {
152 // if not set to nullptr, the reisduals tree will be added to the output file of TrackResiduals created above
153 mTreeResiduals->SetDirectory(nullptr);
154 }
155 for (int iSec = 0; iSec < SECTORSPERSIDE * SIDES; ++iSec) {
156 mResidualsSectorPtr[iSec] = &mResidualsSector[iSec];
157 mVoxStatsSector[iSec].resize(mTrackResiduals.getNVoxelsPerSector());
158 for (int ix = 0; ix < mTrackResiduals.getNXBins(); ++ix) {
159 for (int ip = 0; ip < mTrackResiduals.getNY2XBins(); ++ip) {
160 for (int iz = 0; iz < mTrackResiduals.getNZ2XBins(); ++iz) {
161 auto& statsVoxel = mVoxStatsSector[iSec][mTrackResiduals.getGlbVoxBin(ix, ip, iz)];
162 // COG estimates are set to the bin center by default
163 mTrackResiduals.getVoxelCoordinates(iSec, ix, ip, iz, statsVoxel.meanPos[TrackResiduals::VoxX], statsVoxel.meanPos[TrackResiduals::VoxF], statsVoxel.meanPos[TrackResiduals::VoxZ]);
164 }
165 }
166 }
167 mTreeResiduals->Branch(Form("sec%d", iSec), &mResidualsSectorPtr[iSec]);
168 }
169 // now go through all input files, apply the track selection and fill the binned residuals
170 for (const auto& file : mFileNames) {
171 LOGP(info, "Processing residuals from file {}", file);
172 connectTree(file);
173 for (int iEntry = 0; iEntry < mTreeUnbinnedResiduals->GetEntries(); ++iEntry) {
174 mTreeUnbinnedResiduals->GetEntry(iEntry);
175 if (mParams.useTrackData) {
176 mTreeTrackData->GetEntry(iEntry);
177 }
178 auto nTracks = mTrackData.size();
179 for (size_t iTrack = 0; iTrack < nTracks; ++iTrack) {
180 const auto& trkInfo = mTrackData[iTrack];
181 if (!GID::includesSource(trkInfo.sourceId, mSources)) {
182 continue;
183 }
184 if (mParams.useTrackData) {
185 const auto& trkInfoExtra = mTrackDataExtra[iTrack];
186 if (!revalidateTrack(trkInfoExtra)) {
187 continue;
188 }
189 }
190 for (int i = trkInfo.idxFirstResidual; i < trkInfo.idxFirstResidual + trkInfo.nResiduals; ++i) {
191 const auto& residIn = mUnbinnedResiduals[i];
192 int sec = residIn.sec;
193 auto& residVecOut = mResidualsSector[sec];
194 auto& statVecOut = mVoxStatsSector[sec];
195 std::array<unsigned char, TrackResiduals::VoxDim> bvox;
196 float xPos = param::RowX[residIn.row];
197 float yPos = residIn.y * param::MaxY / 0x7fff + residIn.dy * param::MaxResid / 0x7fff;
198 float zPos = residIn.z * param::MaxZ / 0x7fff + residIn.dz * param::MaxResid / 0x7fff;
199 if (!mTrackResiduals.findVoxelBin(sec, xPos, yPos, zPos, bvox)) {
200 // we are not inside any voxel
201 LOGF(debug, "Dropping residual in sec(%i), x(%f), y(%f), z(%f)", sec, xPos, yPos, zPos);
202 continue;
203 }
204 residVecOut.emplace_back(residIn.dy, residIn.dz, residIn.tgSlp, bvox);
205 auto& stat = statVecOut[mTrackResiduals.getGlbVoxBin(bvox)];
206 float& binEntries = stat.nEntries;
207 float oldEntries = binEntries++;
208 float norm = 1.f / binEntries;
209 // update COG for voxel bvox (update for X only needed in case binning is not per pad row)
210 float xPosInv = 1.f / xPos;
211 stat.meanPos[TrackResiduals::VoxX] = (stat.meanPos[TrackResiduals::VoxX] * oldEntries + xPos) * norm;
212 stat.meanPos[TrackResiduals::VoxF] = (stat.meanPos[TrackResiduals::VoxF] * oldEntries + yPos * xPosInv) * norm;
213 stat.meanPos[TrackResiduals::VoxZ] = (stat.meanPos[TrackResiduals::VoxZ] * oldEntries + zPos * xPosInv) * norm;
214 }
215 }
216 }
217 mTreeResiduals->Fill();
218 for (auto& residuals : mResidualsSector) {
219 residuals.clear();
220 }
221 }
222 }
223
224 for (int iSec = 0; iSec < SECTORSPERSIDE * SIDES; ++iSec) {
225 if (mDoBinning) {
226 // for each sector fill the vector of local residuals from the respective branch
227 auto brResid = mTreeResiduals->GetBranch(Form("sec%d", iSec));
228 brResid->SetAddress(&mResidualsPtr);
229 for (int iEntry = 0; iEntry < brResid->GetEntries(); ++iEntry) {
230 brResid->GetEntry(iEntry);
231 mTrackResiduals.getLocalResVec().insert(mTrackResiduals.getLocalResVec().end(), mResiduals.begin(), mResiduals.end());
232 }
233 mTrackResiduals.setStats(mVoxStatsSector[iSec], iSec);
234 } else {
235 // we read the binned residuals directly from the input files
236 for (const auto& file : mFileNames) {
237 LOGP(info, "Processing residuals from file {}", file);
238 // set up the tree from the input file
239 connectTree(file);
240 // fill the residuals for one sector
241 fillResiduals(iSec);
242 }
243 }
244
245 // do processing
246 mTrackResiduals.processSectorResiduals(iSec);
247 // do cleanup
248 mTrackResiduals.clear();
249 }
250
251 mTrackResiduals.closeOutputFile(); // FIXME remove when map output is handled properly
252
253 // const auto& voxResArray = mTrackResiduals.getVoxelResults(); // array with one vector of results per sector
254 // pc.outputs().snapshot(Output{"GLO", "VOXELRESULTS", 0}, voxResArray); // send results as one large vector?
255
256 pc.services().get<ControlService>().endOfStream();
257 pc.services().get<ControlService>().readyToQuit(QuitRequest::Me);
258}
259
260void TPCResidualReader::connectTree(const std::string& filename)
261{
262 if (!mDoBinning) {
263 // in case we do the binning on-the-fly, we fill these trees manually
264 // and don't want to delete them in between
265 mTreeResiduals.reset(nullptr); // in case it was already loaded
266 mTreeStats.reset(nullptr); // in case it was already loaded
267 }
268 mTreeUnbinnedResiduals.reset(nullptr); // in case it was already loaded
269 mTreeTrackData.reset(nullptr); // in case it was already loaded
270 mFile.reset(TFile::Open(filename.c_str()));
271 assert(mFile && !mFile->IsZombie());
272 if (!mDoBinning) {
273 // we load the already binned residuals and the voxel statistics
274 mTreeResiduals.reset((TTree*)mFile->Get("resid"));
275 assert(mTreeResiduals);
276 mTreeStats.reset((TTree*)mFile->Get("stats"));
277 assert(mTreeStats);
278 LOG(info) << "Loaded tree from " << filename << " with " << mTreeResiduals->GetEntries() << " entries";
279 } else {
280 // we load the unbinned residuals
281 LOG(info) << "Loading the binned residuals";
282 mTreeUnbinnedResiduals.reset((TTree*)mFile->Get("unbinnedResid"));
283 assert(mTreeUnbinnedResiduals);
284 mTreeUnbinnedResiduals->SetBranchAddress("res", &mUnbinnedResidualsPtr);
285 mTreeUnbinnedResiduals->SetBranchAddress("trackInfo", &mTrackDataPtr);
286 LOG(info) << "Loaded tree from " << filename << " with " << mTreeUnbinnedResiduals->GetEntries() << " entries";
287 if (mParams.useTrackData) {
288 mTreeTrackData.reset((TTree*)mFile->Get("trackData"));
289 mTreeTrackData->SetBranchAddress("trk", &mTrackDataExtraPtr);
290 LOG(info) << "Loaded the trackData tree in addition with " << mTreeTrackData->GetEntries() << " entries";
291 assert(mTreeUnbinnedResiduals->GetEntries() == mTreeTrackData->GetEntries());
292 }
293 }
294}
295
296bool TPCResidualReader::revalidateTrack(const TrackData& trk) const
297{
298 if (trk.nClsITS < mParams.minITSNCls) {
299 return false;
300 }
301 if (trk.nClsTPC < mParams.minTPCNCls) {
302 return false;
303 }
304 if (trk.nTrkltsTRD < mParams.minTRDNTrklts) {
305 return false;
306 }
307 return true;
308}
309
311{
312 std::vector<InputSpec> inputs;
313 std::vector<OutputSpec> outputs;
314 // outputs.emplace_back("GLO", "VOXELRESULTS", 0, Lifetime::Timeframe);
315 return DataProcessorSpec{
316 "tpc-residual-reader",
317 inputs,
318 outputs,
319 AlgorithmSpec{adaptFromTask<TPCResidualReader>(doBinning, src)},
320 Options{
321 {"residuals-infiles", VariantType::String, "o2tpc_residuals.root", {"comma-separated list of input files or .txt file containing list of input files"}},
322 {"input-dir", VariantType::String, "none", {"Input directory"}},
323 {"outfile", VariantType::String, "debugVoxRes.root", {"Output file name"}},
324 {"store-binned", VariantType::Bool, false, {"Store the binned residuals together with the voxel results"}},
325 {"dont-check-file-access", VariantType::Bool, false, {"Deactivate check if all files are accessible before adding them to the list of files"}},
326 }};
327}
328
329} // namespace tpc
330} // namespace o2
int32_t i
Helper function to tokenize sequences and ranges of integral numbers.
uint32_t res
Definition RawData.h:0
Definition of the TrackInterpolation class.
Definition of the TrackResiduals class.
std::ostringstream debug
ConfigParamRegistry const & options()
Definition InitContext.h:33
ServiceRegistryRef services()
The services registry associated with this processing context.
void run(ProcessingContext &pc) final
~TPCResidualReader() override=default
void init(InitContext &ic) final
TPCResidualReader(bool doBinning, GID::mask_t src)
std::vector< LocalResid > & getLocalResVec()
void processSectorResiduals(Int_t iSec)
size_t getGlbVoxBin(int ix, int ip, int iz) const
void setStats(const std::vector< TrackResiduals::VoxStats > &statsIn, int iSec)
Set the voxel statistics directly from outside.
void closeOutputFile()
Closes the file with the debug output.
void clear()
clear member to be able to process new sector or new input files
void fillStats(int iSec)
Fill statistics from TTree.
int getNVoxelsPerSector() const
Get the total number of voxels per TPC sector (mNXBins * mNY2XBins * mNZ2XBins)
void createOutputFile(const char *filename="debugVoxRes.root")
Creates a file for the debug output.
void getVoxelCoordinates(int isec, int ix, int ip, int iz, float &x, float &p, float &z) const
bool findVoxelBin(int secID, float x, float y, float z, std::array< unsigned char, VoxDim > &bvox) const
Calculates the bin indices for given x, y, z in sector coordinates.
void init(bool doBinning=true)
std::vector< VoxStats > ** getVoxStatPtr()
GLenum src
Definition glcorearb.h:1767
GLuint GLuint end
Definition glcorearb.h:469
GLuint start
Definition glcorearb.h:469
Defining PrimaryVertex explicitly as messageable.
Definition TFIDInfo.h:20
std::vector< ConfigParamSpec > Options
constexpr unsigned char SECTORSPERSIDE
Definition Defs.h:40
constexpr unsigned char SIDES
Definition Defs.h:41
framework::DataProcessorSpec getTPCResidualReaderSpec(bool doBinning, GID::mask_t src)
create a processor spec
a couple of static helper functions to create timestamp values for CCDB queries or override obsolete ...
std::string filename()
int minTRDNTrklts
min number of TRD space points
int minTPCNCls
min number of TPC clusters
bool useTrackData
if we have the track data available, we can redefine the above cuts for the map creation,...
int minITSNCls
min number of ITS clusters
Structure filled for each track with track quality information and a vector with TPCClusterResiduals.
static std::string rectifyDirectory(const std::string_view p)
static std::string concat_string(Ts const &... ts)
LOG(info)<< "Compressed in "<< sw.CpuTime()<< " s"