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 if (pc.inputs().get<gsl::span<char>>("CTPLumi").size() == sizeof(o2::ctp::LumiInfo)) {
58 lumiPrev = lumiObj = pc.inputs().get<o2::ctp::LumiInfo>("CTPLumi");
59 } else {
60 if (dumRep < maxDumRep && lumiPrev.nHBFCounted == 0 && lumiPrev.nHBFCountedFV0 == 0) {
61 LOGP(alarm, "Previous TF lumi used to substitute dummy input is empty, warning {} of {}", ++dumRep, maxDumRep);
62 }
63 lumiObj = lumiPrev;
64 }
65 setInstLumiCTP(mInstLumiCTPFactor * (mLumiCTPSource == 0 ? lumiObj.getLumi() : lumiObj.getLumiAlt()));
66 if (getLumiScaleType() == 1) {
67 setInstLumi(getInstLumiCTP());
68 }
69 }
70 if (getLumiScaleType() == 2) {
71 float tpcScaler = pc.inputs().get<float>("tpcscaler");
72 setInstLumi(tpcScaler);
73 }
75 LOGP(info, "Setting M-Shape map");
76 const auto mapMShape = pc.inputs().get<o2::gpu::TPCFastTransform*>("mshape");
77 const_cast<o2::gpu::TPCFastTransform*>(mapMShape.get())->rectifyAfterReadingFromFile();
78 mCorrMapMShape = std::unique_ptr<TPCFastTransform>(new TPCFastTransform);
79 mCorrMapMShape->cloneFromObject(*(mapMShape.get()), nullptr);
82 }
83
84 // update inverse in case it is requested
85 if (!mScaleInverse) {
87 }
89}
90
91//________________________________________________________
92void CorrectionMapsLoader::requestCCDBInputs(std::vector<InputSpec>& inputs, std::vector<o2::framework::ConfigParamSpec>& options, const CorrectionMapsLoaderGloOpts& gloOpts)
93{
94 if (gloOpts.lumiMode == 0) {
95 addInput(inputs, {"tpcCorrMap", "TPC", "CorrMap", 0, Lifetime::Condition, ccdbParamSpec(CDBTypeMap.at(CDBType::CalCorrMap), {}, 1)}); // time-dependent
96 addInput(inputs, {"tpcCorrMapRef", "TPC", "CorrMapRef", 0, Lifetime::Condition, ccdbParamSpec(CDBTypeMap.at(CDBType::CalCorrMapRef), {}, 0)}); // load once
97 } else if (gloOpts.lumiMode == 1) {
98 addInput(inputs, {"tpcCorrMap", "TPC", "CorrMap", 0, Lifetime::Condition, ccdbParamSpec(CDBTypeMap.at(CDBType::CalCorrMap), {}, 1)}); // time-dependent
99 addInput(inputs, {"tpcCorrMapRef", "TPC", "CorrMapRef", 0, Lifetime::Condition, ccdbParamSpec(CDBTypeMap.at(CDBType::CalCorrDerivMap), {}, 1)}); // time-dependent
100 } else if (gloOpts.lumiMode == 2) {
101 // for MC corrections
102 addInput(inputs, {"tpcCorrMap", "TPC", "CorrMap", 0, Lifetime::Condition, ccdbParamSpec(CDBTypeMap.at(CDBType::CalCorrMapMC), {}, 1)}); // time-dependent
103 addInput(inputs, {"tpcCorrMapRef", "TPC", "CorrMapRef", 0, Lifetime::Condition, ccdbParamSpec(CDBTypeMap.at(CDBType::CalCorrDerivMapMC), {}, 1)}); // time-dependent
104 } else {
105 LOG(fatal) << "Correction mode unknown! Choose either 0 (default) or 1 (derivative map) for flag corrmap-lumi-mode.";
106 }
107
108 if (gloOpts.requestCTPLumi) {
109 addInput(inputs, {"CTPLumi", "CTP", "LUMI", 0, Lifetime::Timeframe});
110 }
111
112 if (gloOpts.lumiType == 2) {
113 addInput(inputs, {"tpcscaler", o2::header::gDataOriginTPC, "TPCSCALER", 0, Lifetime::Timeframe});
114 }
115
116 addInput(inputs, {"tpcCorrPar", "TPC", "CorrMapParam", 0, Lifetime::Condition, ccdbParamSpec(CDBTypeMap.at(CDBType::CorrMapParam), {}, 0)}); // load once
117
118 if (gloOpts.enableMShapeCorrection) {
119 addInput(inputs, {"mshape", o2::header::gDataOriginTPC, "TPCMSHAPE", 0, Lifetime::Timeframe});
120 }
121 addOptions(options);
122}
123
124//________________________________________________________
125void CorrectionMapsLoader::addOptions(std::vector<ConfigParamSpec>& options)
126{
127 // these are options which should be added at the level of device using TPC corrections
128 // At the moment - nothing, all options are moved to configurable param CorrMapParam
129 addOption(options, ConfigParamSpec{"recalculate-inverse-correction", o2::framework::VariantType::Bool, false, {"recalculate the inverse correction in case lumi mode 1 or 2 is used"}});
130 addOption(options, ConfigParamSpec{"nthreads-inverse-correction", o2::framework::VariantType::Int, 4, {"Number of threads used for calculating the inverse correction (-1=all threads)"}});
131}
132
133//________________________________________________________
134void CorrectionMapsLoader::addGlobalOptions(std::vector<ConfigParamSpec>& options)
135{
136 // these are options which should be added at the workflow level, since they modify the inputs of the devices
137 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"}});
138 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)"}});
139 addOption(options, ConfigParamSpec{"enable-M-shape-correction", o2::framework::VariantType::Bool, false, {"Enable M-shape distortion correction"}});
140 addOption(options, ConfigParamSpec{"disable-ctp-lumi-request", o2::framework::VariantType::Bool, false, {"do not request CTP lumi (regardless what is used for corrections)"}});
141}
142
143//________________________________________________________
145{
147 tpcopt.lumiType = opts.get<int>("lumi-type");
148 tpcopt.lumiMode = opts.get<int>("corrmap-lumi-mode");
149 tpcopt.enableMShapeCorrection = opts.get<bool>("enable-M-shape-correction");
150 tpcopt.requestCTPLumi = !opts.get<bool>("disable-ctp-lumi-request");
151 if (!tpcopt.requestCTPLumi && tpcopt.lumiType == 1) {
152 LOGP(fatal, "Scaling with CTP Lumi is requested but this input is disabled");
153 }
154 return tpcopt;
155}
156
157//________________________________________________________
158void CorrectionMapsLoader::addInput(std::vector<InputSpec>& inputs, InputSpec&& isp)
159{
160 if (std::find(inputs.begin(), inputs.end(), isp) == inputs.end()) {
161 inputs.emplace_back(isp);
162 }
163}
164
165//________________________________________________________
166void CorrectionMapsLoader::addOption(std::vector<ConfigParamSpec>& options, ConfigParamSpec&& osp)
167{
168 if (std::find(options.begin(), options.end(), osp) == options.end()) {
169 options.emplace_back(osp);
170 }
171}
172
173//________________________________________________________
175{
176 if (matcher == ConcreteDataMatcher("TPC", "CorrMap", 0)) {
179 if (getMeanLumiOverride() == 0 && mCorrMap->getLumi() > 0.) {
180 setMeanLumi(mCorrMap->getLumi(), false);
181 }
182 LOGP(debug, "MeanLumiOverride={} MeanLumiMap={} -> meanLumi = {}", getMeanLumiOverride(), mCorrMap->getLumi(), getMeanLumi());
184 return true;
185 }
186 if (matcher == ConcreteDataMatcher("TPC", "CorrMapRef", 0)) {
189 if (getMeanLumiRefOverride() == 0) {
190 setMeanLumiRef(mCorrMapRef->getLumi());
191 }
192 LOGP(debug, "MeanLumiRefOverride={} MeanLumiMap={} -> meanLumi = {}", getMeanLumiRefOverride(), mCorrMapRef->getLumi(), getMeanLumiRef());
194 return true;
195 }
196 if (matcher == ConcreteDataMatcher("TPC", "CorrMapParam", 0)) {
197 const auto& par = o2::tpc::CorrMapParam::Instance();
198 mMeanLumiOverride = par.lumiMean; // negative value switches off corrections !!!
199 mMeanLumiRefOverride = par.lumiMeanRef;
200 mInstCTPLumiOverride = par.lumiInst;
201 mInstLumiCTPFactor = par.lumiInstFactor;
202 mLumiCTPSource = par.ctpLumiSource;
203
204 if (mMeanLumiOverride != 0.) {
206 }
207 if (mMeanLumiRefOverride != 0.) {
209 }
210 if (mInstCTPLumiOverride != 0.) {
212 if (getLumiScaleType() == 1) {
213 setInstLumi(getInstLumiCTP(), false);
214 }
215 }
217 int scaleType = getLumiScaleType();
218 const std::array<std::string, 3> lumiS{"OFF", "CTP", "TPC scaler"};
219 if (scaleType >= lumiS.size()) {
220 LOGP(fatal, "Wrong lumi-scale-type provided!");
221 }
222
223 LOGP(info, "TPC correction map params updated: SP corrections: {} (corr.map scaling type={}, override values: lumiMean={} lumiRefMean={} lumiScaleMode={}), CTP Lumi: source={} lumiInstOverride={} , LumiInst scale={} ",
224 canUseCorrections() ? "ON" : "OFF",
226 }
227 return false;
228}
229
230//________________________________________________________
232{
233 if (getLumiScaleMode() < 0) {
234 LOGP(fatal, "TPC correction lumi scaling mode is not set");
235 }
236 const auto& inputRouts = ic.services().get<const o2::framework::DeviceSpec>().inputs;
237 bool foundCTP = false, foundTPCScl = false, foundMShape = false;
238 for (const auto& route : inputRouts) {
239 if (route.matcher == InputSpec{"CTPLumi", "CTP", "LUMI", 0, Lifetime::Timeframe}) {
240 foundCTP = true;
241 } else if (route.matcher == InputSpec{"tpcscaler", o2::header::gDataOriginTPC, "TPCSCALER", 0, Lifetime::Timeframe}) {
242 foundTPCScl = true;
243 } else if (route.matcher == InputSpec{"mshape", o2::header::gDataOriginTPC, "TPCMSHAPE", 0, Lifetime::Timeframe}) {
244 foundMShape = true;
245 }
246 }
247 setLumiCTPAvailable(foundCTP);
248 enableMShapeCorrection(foundMShape);
249 if ((getLumiScaleType() == 1 && !foundCTP) || (getLumiScaleType() == 2 && !foundTPCScl)) {
250 LOGP(fatal, "Lumi scaling source {}({}) is not available for TPC correction", getLumiScaleType(), getLumiScaleType() == 1 ? "CTP" : "TPCScaler");
251 }
252
253 if ((getLumiScaleMode() == 1) || (getLumiScaleMode() == 2)) {
254 mScaleInverse = !(ic.options().get<bool>("recalculate-inverse-correction"));
255 } else {
256 mScaleInverse = true;
257 }
258 const int nthreadsInv = (ic.options().get<int>("nthreads-inverse-correction"));
260}
261
262//________________________________________________________
264{
265 setInstLumi(src.getInstLumi(), false);
266 setInstLumiCTP(src.getInstLumiCTP());
267 setMeanLumi(src.getMeanLumi(), false);
268 setLumiCTPAvailable(src.getLumiCTPAvailable());
269 setMeanLumiRef(src.getMeanLumiRef());
270 setLumiScaleType(src.getLumiScaleType());
271 setMeanLumiOverride(src.getMeanLumiOverride());
272 setMeanLumiRefOverride(src.getMeanLumiRefOverride());
273 setInstCTPLumiOverride(src.getInstCTPLumiOverride());
274 setLumiScaleMode(src.getLumiScaleMode());
275 enableMShapeCorrection(src.getUseMShapeCorrection());
276 mInstLumiCTPFactor = src.mInstLumiCTPFactor;
277 mLumiCTPSource = src.mLumiCTPSource;
278 mLumiScaleMode = src.mLumiScaleMode;
279 mScaleInverse = src.getScaleInverse();
280}
281
283{
284 if (mLumiScaleMode == 1 || mLumiScaleMode == 2) {
285 LOGP(info, "Recalculating the inverse correction");
287 std::vector<float> scaling{1, mLumiScale};
288 std::vector<o2::gpu::TPCFastSpaceChargeCorrection*> corr{&(mCorrMap->getCorrection()), &(mCorrMapRef->getCorrection())};
289 if (mCorrMapMShape) {
290 scaling.emplace_back(1);
291 corr.emplace_back(&(mCorrMapMShape->getCorrection()));
292 }
294 } else {
295 LOGP(info, "Reinitializing inverse correction with lumi scale mode {} not supported for now", mLumiScaleMode);
296 }
297}
298
299#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.
class to create TPC fast space charge correction
class to create TPC fast transformation
std::ostringstream debug
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)
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)
recalculate inverse correction
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 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
GLsizeiptr size
Definition glcorearb.h:659
constexpr o2::header::DataOrigin gDataOriginTPC
Definition DataHeader.h:576
Defining PrimaryVertex explicitly as messageable.
Definition TFIDInfo.h:20
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:167
const std::unordered_map< CDBType, const std::string > CDBTypeMap
Storage name in CCDB for each calibration and parameter type.
Definition CDBTypes.h:94
@ 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"