Project
Loading...
Searching...
No Matches
CalibPadGainTracks.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 "TPCBase/PadPos.h"
19#include "TPCBase/ROC.h"
20#include "TPCBase/Mapper.h"
24#include "GPUO2InterfaceRefit.h"
25#include "GPUO2Interface.h"
29
30// root includes
31#include "TFile.h"
32#include <random>
33
34using namespace o2::tpc;
35
36void CalibPadGainTracks::processTracks(const int nMaxTracks)
37{
38 std::unique_ptr<o2::gpu::GPUO2InterfaceRefit> refit;
39 if (!mPropagateTrack) {
40 refit = std::make_unique<o2::gpu::GPUO2InterfaceRefit>(mClusterIndex, mTPCCorrMapsHelper, mFieldNominalGPUBz, mTPCTrackClIdxVecInput->data(), 0, mTPCRefitterShMap.data(), mTPCRefitterOccMap.data(), mTPCRefitterOccMap.size());
41 }
42
43 const size_t loopEnd = (nMaxTracks < 0) ? mTracks->size() : ((nMaxTracks > mTracks->size()) ? mTracks->size() : size_t(nMaxTracks));
44
45 if (loopEnd < mTracks->size()) {
46 // draw random tracks
47 std::vector<size_t> ind(mTracks->size());
48 std::iota(ind.begin(), ind.end(), 0);
49 std::minstd_rand rng(std::time(nullptr));
50 std::shuffle(ind.begin(), ind.end(), rng);
51 for (size_t i = 0; i < loopEnd; ++i) {
52 processTrack((*mTracks)[ind[i]], refit.get());
53 }
54 } else {
55 for (const auto& trk : *mTracks) {
56 processTrack(trk, refit.get());
57 }
58 }
59}
60
61void CalibPadGainTracks::processTrack(o2::tpc::TrackTPC track, o2::gpu::GPUO2InterfaceRefit* refit)
62{
63 // make momentum cut
64 const float mom = track.getP();
65 const int nClusters = track.getNClusterReferences();
66 if (mom < mMomMin || mom > mMomMax || std::abs(track.getEta()) > mEtaMax || nClusters < mMinClusters) {
67 return;
68 }
69
70 // clearing memory
71 for (auto& buffer : mDEdxBuffer) {
72 buffer.clear();
73 }
74 mClTrk.clear();
75
76 for (int iCl = 0; iCl < nClusters; iCl++) { // loop over cluster
77 const o2::tpc::ClusterNative& cl = track.getCluster(*mTPCTrackClIdxVecInput, iCl, *mClusterIndex);
78
79 const auto flagsCl = cl.getFlags();
81 continue;
82 }
83
84 unsigned char sectorIndex = 0;
85 unsigned char rowIndex = 0;
86 unsigned int clusterIndexNumb = 0;
87
88 // this function sets sectorIndex, rowIndex, clusterIndexNumb
89 track.getClusterReference(*mTPCTrackClIdxVecInput, iCl, sectorIndex, rowIndex, clusterIndexNumb);
90 const float xPosition = Mapper::instance().getPadCentre(PadPos(rowIndex, 0)).X();
91 if (!mPropagateTrack) {
92 refit->setTrackReferenceX(xPosition);
93 } else {
94 track.rotate(o2::math_utils::detail::sector2Angle<float>(sectorIndex));
95 }
96 const bool check = mPropagateTrack ? o2::base::Propagator::Instance()->PropagateToXBxByBz(track, xPosition, 0.9f, 2., o2::base::Propagator::MatCorrType::USEMatCorrLUT) : ((refit->RefitTrackAsGPU(track, false, true) < 0) ? false : true); // propagate this track to the plane X=xk (cm) in the field "b" (kG)
97
98 if (!check || std::isnan(track.getParam(1))) {
99 continue;
100 }
101
102 const int region = Mapper::REGION[rowIndex];
103 const float charge = (mChargeType == ChargeType::Max) ? cl.qMax : cl.qTot;
104 const float effectiveLength = mCalibTrackTopologyPol ? getTrackTopologyCorrectionPol(track, cl, region, charge) : getTrackTopologyCorrection(track, region);
105
106 const unsigned char pad = std::clamp(static_cast<unsigned int>(cl.getPad() + 0.5f), static_cast<unsigned int>(0), Mapper::PADSPERROW[region][Mapper::getLocalRowFromGlobalRow(rowIndex)] - 1); // the left side of the pad is defined at e.g. 3.5 and the right side at 4.5
107 const float gain = mGainMapRef ? mGainMapRef->getValue(sectorIndex, rowIndex, pad) : 1;
108 const float chargeNorm = charge / (effectiveLength * gain);
109
110 if (mMode == dedxTracking) {
111 const auto& dEdx = track.getdEdx();
112 float dedx = -1;
113 const CRU cru(Sector(sectorIndex), region);
114
115 if (mDedxRegion == stack) {
116 const auto stack = cru.gemStack();
117 if (stack == GEMstack::IROCgem && dEdx.NHitsIROC > mMinClusters) {
118 dedx = getdEdxIROC(dEdx);
119 } else if (stack == GEMstack::OROC1gem && dEdx.NHitsOROC1 > mMinClusters) {
120 dedx = getdEdxOROC1(dEdx);
121 } else if (stack == GEMstack::OROC2gem && dEdx.NHitsOROC2 > mMinClusters) {
122 dedx = getdEdxOROC2(dEdx);
123 } else if (stack == GEMstack::OROC3gem && dEdx.NHitsOROC3 > mMinClusters) {
124 dedx = getdEdxOROC3(dEdx);
125 }
126 } else if (mDedxRegion == chamber) {
127 if (cru.isIROC() && dEdx.NHitsIROC > mMinClusters) {
128 dedx = getdEdxIROC(dEdx);
129 } else {
130 int count = 0;
131 if (dEdx.NHitsOROC1 > mMinClusters) {
132 dedx += getdEdxOROC1(dEdx);
133 ++count;
134 }
135 if (dEdx.NHitsOROC2 > mMinClusters) {
136 dedx += getdEdxOROC2(dEdx);
137 ++count;
138 }
139 if (dEdx.NHitsOROC3 > mMinClusters) {
140 dedx += getdEdxOROC3(dEdx);
141 ++count;
142 }
143 if (count > 0) {
144 dedx /= count;
145 }
146 }
147 } else if (mDedxRegion == sector) {
148 dedx = getdEdxTPC(dEdx);
149 }
150
151 if (dedx <= 0) {
152 continue;
153 }
154
155 if (dedx < mDedxMin || (mDedxMax > 0 && dedx > mDedxMax)) {
156 continue;
157 }
158
159 const float fillVal = mDoNotNormCharge ? chargeNorm : chargeNorm / dedx;
161 if (cru.isOROC()) {
163 }
164
165 fillPadByPadHistogram(cru.roc().getRoc(), index, fillVal);
166 }
167
168 if (mMode == dedxTrack) {
169 const int indexBuffer = getdEdxBufferIndex(region);
170
171 const bool isEdge = (flagsCl & ClusterNative::flagEdge) == ClusterNative::flagEdge;
172 const int isSectorCentre = std::abs(static_cast<int>(pad) - static_cast<int>(Mapper::PADSPERROW[region][rowIndex] / 2)); // do not use clusters at the centre of a sector for dE/dx calculation
173 const int nPadsSector = 1;
174
175 if (!isEdge && (isSectorCentre > nPadsSector)) {
176 mDEdxBuffer[indexBuffer].emplace_back(chargeNorm);
177 }
178
179 mClTrk.emplace_back(std::make_tuple(sectorIndex, rowIndex, pad, chargeNorm)); // fill with dummy dedx value
180 }
181 }
182
183 if (mMode == dedxTrack) {
184 getTruncMean();
185
186 // set the dEdx
187 for (auto& x : mClTrk) {
188 const unsigned char globRow = std::get<1>(x);
189 const int region = Mapper::REGION[globRow];
190 const int indexBuffer = getdEdxBufferIndex(region);
191
192 const float dedxTmp = mDedxTmp[indexBuffer];
193 if (dedxTmp <= 0 || dedxTmp < mDedxMin || (mDedxMax > 0 && dedxTmp > mDedxMax)) {
194 continue;
195 }
196
197 const unsigned char pad = std::get<2>(x);
199 const ROC roc = CRU(Sector(std::get<0>(x)), region).roc();
200 if (roc.isOROC()) {
202 }
203
204 // fill the normalizes charge in pad histogram
205 float fillVal = mDoNotNormCharge ? std::get<3>(x) : std::get<3>(x) / dedxTmp;
206 if (getLogTransformQ()) {
207 fillVal = std::log(1 + fillVal);
208 }
209 fillPadByPadHistogram(roc.getRoc(), index, fillVal);
210 }
211 } else {
212 }
213}
214
215void CalibPadGainTracks::getTruncMean(float low, float high)
216{
217 mDedxTmp.clear();
218 mDedxTmp.reserve(mDEdxBuffer.size());
219 // returns the truncated mean for input vector
220 for (auto& charge : mDEdxBuffer) {
221 const int nClustersUsed = static_cast<int>(charge.size());
222 if (nClustersUsed < mMinClusters) {
223 mDedxTmp.emplace_back(-1);
224 continue;
225 }
226
227 std::sort(charge.begin(), charge.end()); // sort the vector for performing truncated mean
228
229 const int startInd = static_cast<int>(low * nClustersUsed);
230 const int endInd = static_cast<int>(high * nClustersUsed);
231
232 if (endInd <= startInd) {
233 mDedxTmp.emplace_back(-1);
234 continue;
235 }
236
237 const float dEdx = std::accumulate(charge.begin() + startInd, charge.begin() + endInd, 0.f);
238 const int nClustersTrunc = endInd - startInd; // count number of clusters
239 mDedxTmp.emplace_back(dEdx / nClustersTrunc);
240 }
241}
242
243float CalibPadGainTracks::getTrackTopologyCorrection(const o2::tpc::TrackTPC& track, const unsigned int region) const
244{
245 const float padLength = Mapper::instance().getPadRegionInfo(region).getPadHeight();
246 const float sinPhi = track.getSnp();
247 const float tgl = track.getTgl();
248 const float snp2 = sinPhi * sinPhi;
249 const float effectiveLength = padLength * std::sqrt((1 + tgl * tgl) / (1 - snp2)); // calculate the trace length of the track over the pad
250 return effectiveLength;
251}
252
253float CalibPadGainTracks::getTrackTopologyCorrectionPol(const o2::tpc::TrackTPC& track, const o2::tpc::ClusterNative& cl, const unsigned int region, const float charge) const
254{
255 const float trackSnp = std::abs(track.getSnp());
256 float snp2 = trackSnp * trackSnp;
257 const float sec2 = 1.f / (1.f - snp2);
258 const float trackTgl = track.getTgl();
259 const float tgl2 = trackTgl * trackTgl;
260 const float tanTheta = std::sqrt(tgl2 * sec2);
261
262 const float z = std::abs(track.getParam(1));
263 const float padTmp = cl.getPad();
264 const float absRelPad = std::abs(padTmp - int(padTmp + 0.5f));
265 const float relTime = cl.getTime() - int(cl.getTime() + 0.5f);
266
267 const float effectiveLength = (mChargeType == ChargeType::Max) ? mCalibTrackTopologyPol->getCorrectionqMax(region, tanTheta, trackSnp, z, absRelPad, relTime) : mCalibTrackTopologyPol->getCorrectionqTot(region, tanTheta, trackSnp, z, 3.5f /*dummy threshold for now*/, charge);
268 return effectiveLength;
269}
270
271void CalibPadGainTracks::reserveMemory()
272{
273 mClTrk.reserve(Mapper::PADROWS);
274 resizedEdxBuffer();
275}
276
277void CalibPadGainTracks::resizedEdxBuffer()
278{
279 if (mDedxRegion == stack) {
280 mDEdxBuffer.resize(4);
281 mDEdxBuffer[0].reserve(Mapper::getNumberOfRowsInIROC());
282 mDEdxBuffer[1].reserve(Mapper::getNumberOfRowsInOROC());
283 mDEdxBuffer[2].reserve(Mapper::getNumberOfRowsInOROC());
284 mDEdxBuffer[3].reserve(Mapper::getNumberOfRowsInOROC());
285 } else if (mDedxRegion == chamber) {
286 mDEdxBuffer.resize(2);
287 mDEdxBuffer[0].reserve(Mapper::getNumberOfRowsInIROC());
288 mDEdxBuffer[1].reserve(Mapper::getNumberOfRowsInOROC());
289 } else if (mDedxRegion == sector) {
290 mDEdxBuffer.resize(1);
291 mDEdxBuffer[0].reserve(Mapper::instance().getNumberOfRows());
292 } else {
293 LOGP(warning, "wrong dE/dx type");
294 }
295}
296
298{
299 mDedxRegion = dedx;
300 resizedEdxBuffer();
301}
302
303void CalibPadGainTracks::dumpToFile(const char* outFileName, const char* outName) const
304{
305 TFile fOut(outFileName, "RECREATE");
306 fOut.WriteObject(this, outName);
307 fOut.Close();
308}
309
310void CalibPadGainTracks::setMembers(gsl::span<const o2::tpc::TrackTPC>* vTPCTracksArrayInp, gsl::span<const o2::tpc::TPCClRefElem>* tpcTrackClIdxVecInput, const o2::tpc::ClusterNativeAccess& clIndex, gsl::span<const unsigned char> TPCRefitterShMap, gsl::span<const unsigned int> TPCRefitterOccMap)
311{
312 mTracks = vTPCTracksArrayInp;
313 mTPCTrackClIdxVecInput = tpcTrackClIdxVecInput;
314 mClusterIndex = &clIndex;
315 mTPCRefitterShMap = TPCRefitterShMap;
316 mTPCRefitterOccMap = TPCRefitterOccMap;
317}
318
319void CalibPadGainTracks::setMomentumRange(const float momMin, const float momMax)
320{
321 mMomMin = momMin;
322 mMomMax = momMax;
323}
324
325void CalibPadGainTracks::setRefGainMap(const char* inpFile, const char* mapName)
326{
327 TFile f(inpFile, "READ");
328 o2::tpc::CalPad* gainMap = nullptr;
329 f.GetObject(mapName, gainMap);
330
331 if (!gainMap) {
332 LOGP(info, "GainMap {} not found returning", mapName);
333 return;
334 }
335 setRefGainMap(*gainMap);
336 delete gainMap;
337}
338
339int CalibPadGainTracks::getdEdxBufferIndex(const int region) const
340{
341 if (mDedxRegion == stack) {
342 return static_cast<int>(CRU(region).gemStack());
343 } else if (mDedxRegion == chamber) {
344 return static_cast<int>(CRU(region).rocType());
345 } else if (mDedxRegion == sector) {
346 return 0;
347 } else {
348 LOGP(warning, "wrong dE/dx type");
349 return -1;
350 }
351}
352
354{
355 mCalibTrackTopologyPol = std::make_unique<CalibdEdxTrackTopologyPol>();
356 mCalibTrackTopologyPol->loadFromFile(fileName.data(), "CalibdEdxTrackTopologyPol");
357}
358
360{
361 mCalibTrackTopologyPol = std::make_unique<CalibdEdxTrackTopologyPol>();
362 mCalibTrackTopologyPol->setFromContainer(polynomials);
363}
364
365void CalibPadGainTracks::drawRefGainMapHelper(const bool type, const Sector sector, const std::string filename, const float minZ, const float maxZ) const
366{
367 if (!mGainMapRef) {
368 LOGP(error, "Map not set");
369 return;
370 }
371
372 std::function<float(const unsigned int, const unsigned int, const unsigned int, const unsigned int)> idcFunc = [mapTmp = mGainMapRef.get()](const unsigned int sector, const unsigned int region, const unsigned int lrow, const unsigned int pad) {
373 return mapTmp->getValue(sector, Mapper::getGlobalPadNumber(lrow, pad, region));
374 };
375
377 drawFun.mIDCFunc = idcFunc;
378 const std::string zAxisTitle = "rel. gain";
379 type ? IDCDrawHelper::drawSide(drawFun, sector.side(), zAxisTitle, filename, minZ, maxZ) : IDCDrawHelper::drawSector(drawFun, 0, Mapper::NREGIONS, sector, zAxisTitle, filename, minZ, maxZ);
380}
381
382void CalibPadGainTracks::dumpReferenceExtractedGainMap(const char* outFileName, const char* outName) const
383{
384 if (!mGainMapRef) {
385 LOGP(error, "Map not set");
386 return;
387 }
388 CalDet gainMapRef(*mGainMapRef.get());
389 gainMapRef *= getPadGainMap();
390
391 TFile f(outFileName, "RECREATE");
392 f.WriteObject(&gainMapRef, outName);
393}
394
395int CalibPadGainTracks::getIndex(o2::tpc::PadSubset padSub, int padSubsetNumber, const int row, const int pad)
396{
397 return Mapper::instance().getPadNumber(padSub, padSubsetNumber, row, pad);
398}
399
400//______________________________________________
402{
403 mTPCVDrift = v.getVDrift();
404 mTPCVDriftCorrFact = v.corrFact;
405 mTPCVDriftRef = v.refVDrift;
406 mTPCDriftTimeOffset = v.getTimeOffset();
407}
408
409//______________________________________________
411{
412 mTPCCorrMapsHelper = maph;
413}
Class of a TPC cluster in TPC-native coordinates (row, time)
Helper class to access correction maps.
int16_t charge
Definition RawEventData.h:5
int32_t i
helper class for drawing IDCs per region/side
uint32_t roc
Definition RawData.h:3
class to create TPC fast transformation
calibration data from laser track calibration
int nClusters
GPUd() value_type estimateLTFast(o2 static GPUd() float estimateLTIncrement(const o2 PropagatorImpl * Instance(bool uninitialized=false)
Definition Propagator.h:143
int32_t RefitTrackAsGPU(o2::tpc::TrackTPC &trk, bool outward=false, bool resetCov=false)
GEMstack gemStack() const
Definition CRU.h:82
RocType rocType() const
Definition CRU.h:66
const ROC roc() const
Definition CRU.h:61
void fillPadByPadHistogram(const size_t roc, const size_t padInROC, const float val)
void setTPCCorrMaps(o2::gpu::CorrectionMapsHelper *maph)
set cluster correction maps helper
void setTPCVDrift(const o2::tpc::VDriftCorrFact &v)
set VDrift correction
void setPolTopologyCorrectionFromContainer(const CalibdEdxTrackTopologyPolContainer &polynomials)
void setMomentumRange(const float momMin, const float momMax)
void setRefGainMap(const char *inpFile, const char *mapName)
void loadPolTopologyCorrectionFromFile(std::string_view fileName)
void processTracks(const int nMaxTracks=-1)
void setdEdxRegion(const DEdxRegion dedx)
set how the dedx is calculated which is used for normalizing the cluster charge
@ sector
use the dE/dx from the whole sector
@ chamber
use the dE/dx from IROC and OROC
@ stack
use the dE/dx from IROC, OROC1, OROC2, OROC3
@ dedxTrack
normalize qMax using the truncated mean from the track
@ dedxTracking
normalize qMax using the dEdx which was calculated during the tracking
void dumpToFile(const char *outFileName="calPadGainTracks.root", const char *outName="calPadGain") const
void setMembers(gsl::span< const o2::tpc::TrackTPC > *vTPCTracksArrayInp, gsl::span< const o2::tpc::TPCClRefElem > *tpcTrackClIdxVecInput, const o2::tpc::ClusterNativeAccess &clIndex, gsl::span< const unsigned char > TPCRefitterShMap, gsl::span< const unsigned int > TPCRefitterOccMap)
void dumpReferenceExtractedGainMap(const char *outFileName="GainMapRefExtracted.root", const char *outName="GainMap") const
static void drawSide(const IDCDraw &idc, const o2::tpc::Side side, const std::string zAxisTitle, const std::string filename, const float minZ=0, const float maxZ=-1)
static GlobalPadNumber getGlobalPadNumber(const unsigned int lrow, const unsigned int pad, const unsigned int region)
Definition Mapper.h:64
static constexpr auto getNumberOfRowsInIROC()
Definition Mapper.h:300
static const std::vector< unsigned int > PADSPERROW[NREGIONS]
number of pads per row in region
Definition Mapper.h:567
static constexpr unsigned int GLOBALPADOFFSET[NREGIONS]
offset of number of pads for region
Definition Mapper.h:531
static unsigned int getLocalRowFromGlobalRow(const unsigned int row)
Definition Mapper.h:72
static Mapper & instance(const std::string mappingDir="")
Definition Mapper.h:44
GlobalPadNumber getPadNumber(const PadSubset padSubset, const size_t padSubsetNumber, const int row, const int pad) const
Definition Mapper.h:134
static constexpr unsigned short getPadsInIROC()
Definition Mapper.h:409
const PadRegionInfo & getPadRegionInfo(const unsigned char region) const
Definition Mapper.h:385
static constexpr auto getNumberOfRowsInOROC()
Definition Mapper.h:303
static constexpr unsigned int PADROWS
total number of pad rows
Definition Mapper.h:528
static constexpr unsigned int OFFSETCRUGLOBAL[PADROWS]
row offset in cru for given global pad row
Definition Mapper.h:579
GlobalPosition2D getPadCentre(const PadSecPos &padSec) const
Definition Mapper.h:163
static constexpr unsigned REGION[PADROWS]
region for global pad row
Definition Mapper.h:537
float getPadHeight() const
GLint GLenum GLint x
Definition glcorearb.h:403
GLint GLsizei count
Definition glcorearb.h:399
GLuint buffer
Definition glcorearb.h:655
GLsizeiptr size
Definition glcorearb.h:659
const GLdouble * v
Definition glcorearb.h:832
GLuint index
Definition glcorearb.h:781
GLdouble f
Definition glcorearb.h:310
GLint GLint GLsizei GLint GLenum GLenum type
Definition glcorearb.h:275
GLfloat GLfloat GLfloat GLfloat GLfloat GLfloat maxZ
Definition glcorearb.h:2910
GLdouble GLdouble GLdouble z
Definition glcorearb.h:843
GLfloat GLfloat minZ
Definition glcorearb.h:2910
void check(const std::vector< std::string > &arguments, const std::vector< ConfigParamSpec > &workflowOptions, const std::vector< DeviceSpec > &deviceSpecs, CheckMatrix &matrix)
auto get(const std::byte *buffer, size_t=0)
Definition DataHeader.h:454
Global TPC definitions and constants.
Definition SimTraits.h:167
@ IROCgem
Definition Defs.h:53
@ OROC1gem
Definition Defs.h:54
@ OROC3gem
Definition Defs.h:56
@ OROC2gem
Definition Defs.h:55
@ Max
Definition Defs.h:71
PadSubset
Definition of the different pad subsets.
Definition Defs.h:63
Defining DataPointCompositeObject explicitly as copiable.
std::string filename()
simple struct to enable writing the MultivariatePolynomialCT to file
std::function< float(const unsigned int, const unsigned int, const unsigned int, const unsigned int)> mIDCFunc
function returning the value which will be drawn for sector, region, row, pad
std::vector< int > row