Project
Loading...
Searching...
No Matches
CorrectionMapsLoader.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#include "Framework/Logger.h"
27
28using namespace o2::tpc;
29using namespace o2::framework;
30
31#ifndef GPUCA_GPUCODE_DEVICE
32
33//________________________________________________________
34void CorrectionMapsLoader::updateVDrift(float vdriftCorr, float vdrifRef, float driftTimeOffset)
35{
36 o2::tpc::TPCFastTransformHelperO2::instance()->updateCalibration(*mCorrMap, 0, vdriftCorr, vdrifRef, driftTimeOffset);
37 if (mCorrMapRef) {
38 o2::tpc::TPCFastTransformHelperO2::instance()->updateCalibration(*mCorrMapRef, 0, vdriftCorr, vdrifRef, driftTimeOffset);
39 }
40 if (mCorrMapMShape) {
41 o2::tpc::TPCFastTransformHelperO2::instance()->updateCalibration(*mCorrMapMShape, 0, vdriftCorr, vdrifRef, driftTimeOffset);
42 }
43}
44
45//________________________________________________________
47{
48 pc.inputs().get<o2::tpc::CorrMapParam*>("tpcCorrPar");
49 pc.inputs().get<o2::gpu::TPCFastTransform*>("tpcCorrMap");
50 pc.inputs().get<o2::gpu::TPCFastTransform*>("tpcCorrMapRef");
51 const int maxDumRep = 5;
52 int dumRep = 0;
53 o2::ctp::LumiInfo lumiObj;
54 static o2::ctp::LumiInfo lumiPrev;
55
57 float tpcScaler = pc.inputs().get<float>("tpcscaler");
58 // check if tpcScaler is valid and CTP fallback is allowed
59 if (tpcScaler == -1.f) {
60 const bool canUseCTPScaling = mCorrMap && mCorrMapRef && mCorrMap->isIDCSet() && mCorrMapRef->isIDCSet() && mCorrMap->isLumiSet() && mCorrMapRef->isLumiSet();
61 if (canUseCTPScaling) {
62 LOGP(info, "Invalid TPC scaler value {} received for IDC-based scaling! Using CTP fallback", tpcScaler);
64 setMeanLumi(mCorrMap->getLumi(), false);
65 setMeanLumiRef(mCorrMapRef->getLumi());
67 } else if (mCorrMap) {
68 // CTP scaling is not possible, dont do any scaling to avoid applying wrong corrections
69 const float storedIDC = mCorrMap->getIDC();
70 LOGP(warning, "Invalid TPC scaler value {} received for IDC-based scaling! CTP fallback not possible, using stored IDC of {} from the map to avoid applying wrong corrections", tpcScaler, storedIDC);
71 setInstLumi(storedIDC);
72 }
73 } else {
75 // reset back to normal operation
76 LOGP(info, "Valid TPC scaler value {} received, switching back to IDC-based scaling", tpcScaler);
78 setMeanLumi(mCorrMap->getIDC(), false);
79 setMeanLumiRef(mCorrMapRef->getIDC());
81 }
82 // correct IDC received
83 setInstLumi(tpcScaler);
84 }
85 }
86
88 if (pc.inputs().get<gsl::span<char>>("CTPLumi").size() == sizeof(o2::ctp::LumiInfo)) {
89 lumiPrev = lumiObj = pc.inputs().get<o2::ctp::LumiInfo>("CTPLumi");
90 } else {
91 if (dumRep < maxDumRep && lumiPrev.nHBFCounted == 0 && lumiPrev.nHBFCountedFV0 == 0) {
92 LOGP(alarm, "Previous TF lumi used to substitute dummy input is empty, warning {} of {}", ++dumRep, maxDumRep);
93 }
94 lumiObj = lumiPrev;
95 }
96 setInstLumiCTP(mInstLumiCTPFactor * (mLumiCTPSource == 0 ? lumiObj.getLumi() : lumiObj.getLumiAlt()));
97 if (getLumiScaleType() == 1) {
98 setInstLumi(getInstLumiCTP());
99 }
100 }
101
103 LOGP(info, "Setting M-Shape map");
104 const auto mapMShape = pc.inputs().get<o2::gpu::TPCFastTransform*>("mshape");
105 const_cast<o2::gpu::TPCFastTransform*>(mapMShape.get())->rectifyAfterReadingFromFile();
106 mCorrMapMShape = std::unique_ptr<TPCFastTransform>(new TPCFastTransform);
107 mCorrMapMShape->cloneFromObject(*(mapMShape.get()), nullptr);
110 }
111
112 // update inverse in case it is requested
113 if (!mScaleInverse) {
115 }
117}
118
119//________________________________________________________
120void CorrectionMapsLoader::requestCCDBInputs(std::vector<InputSpec>& inputs, std::vector<o2::framework::ConfigParamSpec>& options, const CorrectionMapsLoaderGloOpts& gloOpts)
121{
122 if (gloOpts.lumiMode == 0) {
123 addInput(inputs, {"tpcCorrMap", "TPC", "CorrMap", 0, Lifetime::Condition, ccdbParamSpec(CDBTypeMap.at(CDBType::CalCorrMap), {}, 1)}); // time-dependent
124 addInput(inputs, {"tpcCorrMapRef", "TPC", "CorrMapRef", 0, Lifetime::Condition, ccdbParamSpec(CDBTypeMap.at(CDBType::CalCorrMapRef), {}, 0)}); // load once
125 } else if (gloOpts.lumiMode == 1) {
126 addInput(inputs, {"tpcCorrMap", "TPC", "CorrMap", 0, Lifetime::Condition, ccdbParamSpec(CDBTypeMap.at(CDBType::CalCorrMap), {}, 1)}); // time-dependent
127 addInput(inputs, {"tpcCorrMapRef", "TPC", "CorrMapRef", 0, Lifetime::Condition, ccdbParamSpec(CDBTypeMap.at(CDBType::CalCorrDerivMap), {}, 1)}); // time-dependent
128 } else if (gloOpts.lumiMode == 2) {
129 // for MC corrections
130 addInput(inputs, {"tpcCorrMap", "TPC", "CorrMap", 0, Lifetime::Condition, ccdbParamSpec(CDBTypeMap.at(CDBType::CalCorrMapMC), {}, 1)}); // time-dependent
131 addInput(inputs, {"tpcCorrMapRef", "TPC", "CorrMapRef", 0, Lifetime::Condition, ccdbParamSpec(CDBTypeMap.at(CDBType::CalCorrDerivMapMC), {}, 1)}); // time-dependent
132 } else {
133 LOG(fatal) << "Correction mode unknown! Choose either 0 (default) or 1 (derivative map) for flag corrmap-lumi-mode.";
134 }
135
136 if (gloOpts.requestCTPLumi) {
137 addInput(inputs, {"CTPLumi", "CTP", "LUMI", 0, Lifetime::Timeframe});
138 }
139
140 if (gloOpts.lumiType == 2) {
141 addInput(inputs, {"tpcscaler", o2::header::gDataOriginTPC, "TPCSCALER", 0, Lifetime::Timeframe});
142 }
143
144 addInput(inputs, {"tpcCorrPar", "TPC", "CorrMapParam", 0, Lifetime::Condition, ccdbParamSpec(CDBTypeMap.at(CDBType::CorrMapParam), {}, 0)}); // load once
145
146 if (gloOpts.enableMShapeCorrection) {
147 addInput(inputs, {"mshape", o2::header::gDataOriginTPC, "TPCMSHAPE", 0, Lifetime::Timeframe});
148 }
149 addOptions(options);
150}
151
152//________________________________________________________
153void CorrectionMapsLoader::addOptions(std::vector<ConfigParamSpec>& options)
154{
155 // these are options which should be added at the level of device using TPC corrections
156 // At the moment - nothing, all options are moved to configurable param CorrMapParam
157 addOption(options, ConfigParamSpec{"recalculate-inverse-correction", o2::framework::VariantType::Bool, false, {"recalculate the inverse correction in case lumi mode 1 or 2 is used"}});
158 addOption(options, ConfigParamSpec{"nthreads-inverse-correction", o2::framework::VariantType::Int, 4, {"Number of threads used for calculating the inverse correction (-1=all threads)"}});
159}
160
161//________________________________________________________
162void CorrectionMapsLoader::addGlobalOptions(std::vector<ConfigParamSpec>& options)
163{
164 // these are options which should be added at the workflow level, since they modify the inputs of the devices
165 addOption(options, ConfigParamSpec{"lumi-type", o2::framework::VariantType::Int, 0, {"1 = use CTP lumi for TPC correction scaling, 2 = use TPC scalers for TPC correction scaling"}});
166 addOption(options, ConfigParamSpec{"corrmap-lumi-mode", o2::framework::VariantType::Int, 0, {"scaling mode: (default) 0 = static + scale * full; 1 = full + scale * derivative; 2 = full + scale * derivative (for MC)"}});
167 addOption(options, ConfigParamSpec{"enable-M-shape-correction", o2::framework::VariantType::Bool, false, {"Enable M-shape distortion correction"}});
168 addOption(options, ConfigParamSpec{"disable-ctp-lumi-request", o2::framework::VariantType::Bool, false, {"do not request CTP lumi (regardless what is used for corrections)"}});
169 addOption(options, ConfigParamSpec{"disable-lumi-type-consistency-check", o2::framework::VariantType::Bool, false, {"disable check of selected CTP or IDC scaling source being consistent with the map"}});
170}
171
172//________________________________________________________
174{
176 tpcopt.lumiType = opts.get<int>("lumi-type");
177 tpcopt.lumiMode = opts.get<int>("corrmap-lumi-mode");
178 tpcopt.enableMShapeCorrection = opts.get<bool>("enable-M-shape-correction");
179 tpcopt.requestCTPLumi = !opts.get<bool>("disable-ctp-lumi-request");
180 tpcopt.checkCTPIDCconsistency = !opts.get<bool>("disable-lumi-type-consistency-check");
181 if (!tpcopt.requestCTPLumi && tpcopt.lumiType == 1) {
182 LOGP(fatal, "Scaling with CTP Lumi is requested but this input is disabled");
183 }
184 return tpcopt;
185}
186
187//________________________________________________________
188void CorrectionMapsLoader::addInput(std::vector<InputSpec>& inputs, InputSpec&& isp)
189{
190 if (std::find(inputs.begin(), inputs.end(), isp) == inputs.end()) {
191 inputs.emplace_back(isp);
192 }
193}
194
195//________________________________________________________
196void CorrectionMapsLoader::addOption(std::vector<ConfigParamSpec>& options, ConfigParamSpec&& osp)
197{
198 if (std::find(options.begin(), options.end(), osp) == options.end()) {
199 options.emplace_back(osp);
200 }
201}
202
203//________________________________________________________
205{
206 if (matcher == ConcreteDataMatcher("TPC", "CorrMap", 0)) {
210 if (getMeanLumiOverride() != 0) {
211 if (getLumiScaleType() == 1) {
213 LOGP(info, "CorrMap mean lumi rate is overridden to {}", mCorrMap->getLumi());
214 } else if (getLumiScaleType() == 2) {
216 LOGP(info, "CorrMap mean IDC rate is overridden to {}", mCorrMap->getIDC());
217 }
218 }
219 float mapMeanRate = 0;
220 if (getLumiScaleType() == 1) {
221 mapMeanRate = mCorrMap->getLumi();
222 } else if (getLumiScaleType() == 2) {
223 mapMeanRate = mCorrMap->getIDC();
224 }
226 checkMeanScaleConsistency(mapMeanRate, mCorrMap->getCTP2IDCFallBackThreshold());
227 }
228 if (getMeanLumiOverride() == 0 && mapMeanRate > 0.) {
229 setMeanLumi(mapMeanRate, false);
230 }
231 LOGP(debug, "MeanLumiOverride={} MeanLumiMap={} -> meanLumi = {}", getMeanLumiOverride(), mapMeanRate, getMeanLumi());
233 return true;
234 }
235 if (matcher == ConcreteDataMatcher("TPC", "CorrMapRef", 0)) {
239 if (getMeanLumiRefOverride() != 0) {
240 if (getLumiScaleType() == 1) {
242 LOGP(info, "CorrMapRef mean lumi rate is overridden to {}", mCorrMapRef->getLumi());
243 } else if (getLumiScaleType() == 2) {
245 LOGP(info, "CorrMapRef mean IDC rate is overridden to {}", mCorrMapRef->getIDC());
246 }
247 }
248 float mapRefMeanRate = 0;
249 if (getLumiScaleType() == 1) {
250 mapRefMeanRate = mCorrMapRef->getLumi();
251 } else if (getLumiScaleType() == 2) {
252 mapRefMeanRate = mCorrMapRef->getIDC();
253 }
255 checkMeanScaleConsistency(mapRefMeanRate, mCorrMapRef->getCTP2IDCFallBackThreshold());
256 }
257 if (getMeanLumiRefOverride() == 0) {
258 setMeanLumiRef(mapRefMeanRate);
259 }
260 LOGP(debug, "MeanLumiRefOverride={} MeanLumiMap={} -> meanLumi = {}", getMeanLumiRefOverride(), mapRefMeanRate, getMeanLumiRef());
262 return true;
263 }
264 if (matcher == ConcreteDataMatcher("TPC", "CorrMapParam", 0)) {
265 const auto& par = o2::tpc::CorrMapParam::Instance();
266 mMeanLumiOverride = par.lumiMean; // negative value switches off corrections !!!
267 mMeanLumiRefOverride = par.lumiMeanRef;
268 mInstCTPLumiOverride = par.lumiInst;
269 mInstLumiCTPFactor = par.lumiInstFactor;
270 mLumiCTPSource = par.ctpLumiSource;
271
272 if (mMeanLumiOverride != 0.) {
274 }
275 if (mMeanLumiRefOverride != 0.) {
277 }
278 if (mInstCTPLumiOverride != 0.) {
280 if (getLumiScaleType() == 1) {
281 setInstLumi(getInstLumiCTP(), false);
282 }
283 }
285 int scaleType = getLumiScaleType();
286 const std::array<std::string, 3> lumiS{"OFF", "CTP", "TPC scaler"};
287 if (scaleType >= lumiS.size()) {
288 LOGP(fatal, "Wrong corrmap-lumi-mode provided!");
289 }
290
291 LOGP(info, "TPC correction map params updated: SP corrections: {} (corr.map scaling type={}, override values: lumiMean={} lumiRefMean={} lumiScaleMode={}), CTP Lumi: source={} lumiInstOverride={} , LumiInst scale={} ",
292 canUseCorrections() ? "ON" : "OFF",
294 }
295 return false;
296}
297
298//________________________________________________________
300{
301 if (getLumiScaleMode() < 0) {
302 LOGP(fatal, "TPC correction lumi scaling mode is not set");
303 }
304 const auto& inputRouts = ic.services().get<const o2::framework::DeviceSpec>().inputs;
305 bool foundCTP = false, foundTPCScl = false, foundMShape = false;
306 for (const auto& route : inputRouts) {
307 if (route.matcher == InputSpec{"CTPLumi", "CTP", "LUMI", 0, Lifetime::Timeframe}) {
308 foundCTP = true;
309 } else if (route.matcher == InputSpec{"tpcscaler", o2::header::gDataOriginTPC, "TPCSCALER", 0, Lifetime::Timeframe}) {
310 foundTPCScl = true;
311 } else if (route.matcher == InputSpec{"mshape", o2::header::gDataOriginTPC, "TPCMSHAPE", 0, Lifetime::Timeframe}) {
312 foundMShape = true;
313 }
314 }
315 setLumiCTPAvailable(foundCTP);
316 enableMShapeCorrection(foundMShape);
317 if ((getLumiScaleType() == 1 && !foundCTP) || (getLumiScaleType() == 2 && !foundTPCScl)) {
318 LOGP(fatal, "Lumi scaling source {}({}) is not available for TPC correction", getLumiScaleType(), getLumiScaleType() == 1 ? "CTP" : "TPCScaler");
319 }
320
321 if ((getLumiScaleMode() == 1) || (getLumiScaleMode() == 2)) {
322 mScaleInverse = !(ic.options().get<bool>("recalculate-inverse-correction"));
323 } else {
324 mScaleInverse = true;
325 }
326 const int nthreadsInv = (ic.options().get<int>("nthreads-inverse-correction"));
328}
329
330//________________________________________________________
332{
333 setInstLumi(src.getInstLumi(), false);
334 setInstLumiCTP(src.getInstLumiCTP());
335 setMeanLumi(src.getMeanLumi(), false);
336 setLumiCTPAvailable(src.getLumiCTPAvailable());
337 setMeanLumiRef(src.getMeanLumiRef());
338 setLumiScaleType(src.getLumiScaleType());
339 setMeanLumiOverride(src.getMeanLumiOverride());
340 setMeanLumiRefOverride(src.getMeanLumiRefOverride());
341 setInstCTPLumiOverride(src.getInstCTPLumiOverride());
342 setLumiScaleMode(src.getLumiScaleMode());
343 enableMShapeCorrection(src.getUseMShapeCorrection());
344 mInstLumiCTPFactor = src.mInstLumiCTPFactor;
345 mLumiCTPSource = src.mLumiCTPSource;
346 mLumiScaleMode = src.mLumiScaleMode;
347 mScaleInverse = src.getScaleInverse();
348 mIDC2CTPFallbackActive = src.mIDC2CTPFallbackActive;
349}
350
352{
353 if (mLumiScaleMode == 1 || mLumiScaleMode == 2) {
354 LOGP(info, "Recalculating the inverse correction");
356 std::vector<float> scaling{1, mLumiScale};
357 std::vector<o2::gpu::TPCFastSpaceChargeCorrection*> corr{&(mCorrMap->getCorrection()), &(mCorrMapRef->getCorrection())};
358 if (mCorrMapMShape) {
359 scaling.emplace_back(1);
360 corr.emplace_back(&(mCorrMapMShape->getCorrection()));
361 }
363 } else {
364 LOGP(info, "Reinitializing inverse correction with lumi scale mode {} not supported for now", mLumiScaleMode);
365 }
366}
367
368void CorrectionMapsLoader::checkMeanScaleConsistency(float meanLumi, float threshold) const
369{
370 if (getLumiScaleType() == 1) {
371 if (meanLumi < threshold) {
372 LOGP(fatal, "CTP Lumi scaling source is requested, but the map mean scale {} is below the threshold {}", meanLumi, threshold);
373 }
374 } else if (getLumiScaleType() == 2) {
375 if (meanLumi > threshold) {
376 LOGP(fatal, "IDC scaling source is requested, but the map mean scale {} is above the threshold {}", meanLumi, threshold);
377 }
378 }
379}
380
381#endif // #ifndef GPUCA_GPUCODE_DEVICE
Simple interface to the CDB manager.
Implementation of the parameter class for the CorrectionMapsLoader options.
Helper class to access load maps from CCDB.
std::ostringstream debug
class to create TPC fast space charge correction
class to create TPC fast transformation
ServiceRegistryRef services()
Definition InitContext.h:34
decltype(auto) get(R binding, int part=0) const
InputRecord & inputs()
The inputs associated with this processing context.
void setCorrMapMShape(o2::gpu::TPCFastTransform *m)
void setCorrMapRef(o2::gpu::TPCFastTransform *m)
o2::gpu::TPCFastTransform * mCorrMap
void setMeanLumiRef(float v, bool report=false)
void setInstLumi(float v, bool report=false)
void setMeanLumi(float v, bool report=false)
o2::gpu::TPCFastTransform * mCorrMapRef
void setCorrMap(std::unique_ptr< o2::gpu::TPCFastTransform > &&m)
void setLumi(float l)
Set Lumi info.
void setCTP2IDCFallBackThreshold(float v)
TPCFastSpaceChargeCorrection & getCorrection()
Gives a reference for external initialization of TPC corrections.
void extractCCDBInputs(o2::framework::ProcessingContext &pc)
void updateVDrift(float vdriftCorr, float vdrifRef, float driftTimeOffset=0)
static void addGlobalOptions(std::vector< o2::framework::ConfigParamSpec > &options)
bool accountCCDBInputs(const o2::framework::ConcreteDataMatcher &matcher, void *obj)
static void addOptions(std::vector< o2::framework::ConfigParamSpec > &options)
static void requestCCDBInputs(std::vector< o2::framework::InputSpec > &inputs, std::vector< o2::framework::ConfigParamSpec > &options, const CorrectionMapsLoaderGloOpts &gloOpts)
void init(o2::framework::InitContext &ic)
static void addOption(std::vector< o2::framework::ConfigParamSpec > &options, o2::framework::ConfigParamSpec &&osp)
std::unique_ptr< o2::gpu::TPCFastTransform > mCorrMapMShape
void checkMeanScaleConsistency(float meanLumi, float threshold) const
recalculate inverse correction
void copySettings(const CorrectionMapsLoader &src)
static CorrectionMapsLoaderGloOpts parseGlobalOptions(const o2::framework::ConfigParamRegistry &opts)
static void addInput(std::vector< o2::framework::InputSpec > &inputs, o2::framework::InputSpec &&isp)
void initInverse(o2::gpu::TPCFastSpaceChargeCorrection &correction, bool prn)
initialise inverse transformation
void setNthreads(int n)
_______________ Settings ________________________
void setNthreadsToMaximum()
sets number of threads to N cpu cores
static TPCFastSpaceChargeCorrectionHelper * instance()
Singleton.
int updateCalibration(TPCFastTransform &transform, Long_t TimeStamp, float vDriftFactor=1.f, float vDriftRef=0.f, float driftTimeOffset=0.f)
Updates the transformation with the new time stamp.
static TPCFastTransformHelperO2 * instance()
Singleton.
GLenum src
Definition glcorearb.h:1767
constexpr o2::header::DataOrigin gDataOriginTPC
Definition DataHeader.h:576
Defining PrimaryVertex explicitly as messageable.
std::vector< ConfigParamSpec > ccdbParamSpec(std::string const &path, int runDependent, std::vector< CCDBMetadata > metadata={}, int qrate=0)
Global TPC definitions and constants.
Definition SimTraits.h:168
const std::unordered_map< CDBType, const std::string > CDBTypeMap
Storage name in CCDB for each calibration and parameter type.
Definition CDBTypes.h:96
@ CalCorrMapMC
Cluster correction map (high IR rate distortions) for MC.
@ CalCorrMapRef
Cluster correction reference map (static distortions)
@ CalCorrMap
Cluster correction map (high IR rate distortions)
float getLumiAlt() const
Definition LumiInfo.h:34
uint32_t nHBFCounted
Definition LumiInfo.h:26
uint32_t nHBFCountedFV0
Definition LumiInfo.h:27
float getLumi() const
Definition LumiInfo.h:32
int lumiType
what estimator to used for corrections scaling: 0: no scaling, 1: CTP, 2: IDC
int lumiMode
what corrections method to use: 0: classical scaling, 1: Using of the derivative map,...
LOG(info)<< "Compressed in "<< sw.CpuTime()<< " s"