Project
Loading...
Searching...
No Matches
TPCScalerSpec.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
17#include "Framework/Task.h"
18#include "Framework/Logger.h"
26#include "TTree.h"
32
33using namespace o2::framework;
34
35namespace o2
36{
37namespace tpc
38{
39
40class TPCScalerSpec : public Task
41{
42 public:
43 TPCScalerSpec(std::shared_ptr<o2::base::GRPGeomRequest> req, const o2::tpc::CorrectionMapsGloOpts& sclOpts, bool enableIDCs, bool enableMShape) : mCCDBRequest(req), mEnableIDCs(enableIDCs), mEnableMShape(enableMShape), mGlobOpts(sclOpts)
44 {
45 mTPCCorrMapsLoader.setLumiScaleType(sclOpts.lumiType);
46 mTPCCorrMapsLoader.setLumiScaleMode(sclOpts.lumiMode);
47 mTPCCorrMapsLoader.setCheckCTPIDCConsistency(sclOpts.checkCTPIDCconsistency);
48 };
49
51 {
53 mIonDriftTimeMS = ic.options().get<float>("ion-drift-time");
54 mMaxTimeWeightsMS = ic.options().get<float>("max-time-for-weights");
55 mMShapeScalingFac = ic.options().get<float>("m-shape-scaling-factor");
56 mEnableWeights = !(ic.options().get<bool>("disableWeights"));
57 const bool enableStreamer = ic.options().get<bool>("enableStreamer");
58 const int mshapeThreads = ic.options().get<int>("n-threads");
59 mKnotsYMshape = ic.options().get<int>("n-knots-y");
60 mKnotsZMshape = ic.options().get<int>("n-knots-z");
63
64 if (enableStreamer) {
65 mStreamer = std::make_unique<o2::utils::TreeStreamRedirector>("M_Shape.root", "recreate");
66 }
67 mTPCCorrMapsLoader.init(ic, mEnableIDCs);
68 }
69
71 {
72 if (mStreamer) {
73 mStreamer->Close();
74 }
75 }
76
77 void run(ProcessingContext& pc) final
78 {
80 mTPCVDriftHelper.extractCCDBInputs(pc);
81 if (mTPCVDriftHelper.isUpdated()) {
82 mTPCVDriftHelper.acknowledgeUpdate();
83 }
84
85 if (mEnableIDCs && pc.inputs().isValid("tpcscaler")) {
86 pc.inputs().get<TTree*>("tpcscaler");
87 }
88
89 if (mEnableWeights && mEnableIDCs) {
90 if (pc.inputs().isValid("tpcscalerw")) {
91 pc.inputs().get<TPCScalerWeights*>("tpcscalerw");
92 }
93 }
94 if (mEnableMShape && pc.inputs().isValid("mshape")) {
95 pc.inputs().get<TTree*>("mshape");
96 }
97
98 const auto orbitResetTimeMS = o2::base::GRPGeomHelper::instance().getOrbitResetTimeMS();
99 const auto firstTFOrbit = pc.services().get<o2::framework::TimingInfo>().firstTForbit;
100 const double timestamp = orbitResetTimeMS + firstTFOrbit * o2::constants::lhc::LHCOrbitMUS * 0.001;
101 int currRun = pc.services().get<o2::framework::TimingInfo>().runNumber;
102 if (mEnableMShape) {
103 static int runWarningMS = -1;
104 if ((mMShapeTPCScaler.getRun() != -1) && currRun != mMShapeTPCScaler.getRun() && runWarningMS != currRun) {
105 LOGP(alarm, "Run number {} of processed data and run number {} of loaded TPC M-shape scaler doesnt match!", pc.services().get<o2::framework::TimingInfo>().runNumber, mMShapeTPCScaler.getRun());
106 runWarningMS = currRun;
107 }
108
109 const auto& boundaryPotential = mMShapeTPCScaler.getBoundaryPotential(timestamp);
110 if (!boundaryPotential.mPotential.empty()) {
111 LOGP(info, "Calculating M-shape correction from input boundary potential");
112 const Side side = Side::A;
113 o2::tpc::SpaceCharge<double> sc = o2::tpc::SpaceCharge<double>(mMShapeTPCScaler.getMShapes().mBField, mMShapeTPCScaler.getMShapes().mZ, mMShapeTPCScaler.getMShapes().mR, mMShapeTPCScaler.getMShapes().mPhi);
114 for (int iz = 0; iz < mMShapeTPCScaler.getMShapes().mZ; ++iz) {
115 for (int iphi = 0; iphi < mMShapeTPCScaler.getMShapes().mPhi; ++iphi) {
116 const float pot = mMShapeScalingFac * boundaryPotential.mPotential[iz];
117 sc.setPotential(iz, 0, iphi, side, pot);
118 }
119 }
120
122 sc.calcEField(side);
124
125 std::function<void(int, double, double, double, double&, double&, double&)> getCorrections = [&sc = sc](const int roc, double x, double y, double z, double& dx, double& dy, double& dz) {
126 Side side = roc < 18 ? Side::A : Side::C;
127 if (side == Side::C) {
128 dx = 0;
129 dy = 0;
130 dz = 0;
131 } else {
132 sc.getCorrections(x, y, z, side, dx, dy, dz);
133 }
134 };
135
136 std::unique_ptr<TPCFastSpaceChargeCorrection> spCorrection = TPCFastSpaceChargeCorrectionHelper::instance()->createFromGlobalCorrection(getCorrections, mKnotsYMshape, mKnotsZMshape);
137 std::unique_ptr<TPCFastTransform> fastTransform(TPCFastTransformHelperO2::instance()->create(0, *spCorrection));
138 mTPCCorrMapsLoader.setCorrMapMShape(std::move(fastTransform));
139 }
140
141 if (mStreamer) {
142 (*mStreamer) << "treeMShape"
143 << "firstTFOrbit=" << firstTFOrbit
144 << "timestamp=" << timestamp
145 << "boundaryPotential=" << boundaryPotential
146 << "mMShapeScalingFac=" << mMShapeScalingFac
147 << "\n";
148 }
149 }
150
151 float tpcScaler = -1.f;
152 if (mEnableIDCs) {
153 static int runWarningIDC = -1;
154 if (pc.services().get<o2::framework::TimingInfo>().runNumber != mTPCScaler.getRun() && runWarningIDC != currRun) {
155 LOGP(alarm, "Run number {} of processed data and run number {} of loaded TPC scaler doesnt match!", pc.services().get<o2::framework::TimingInfo>().runNumber, mTPCScaler.getRun());
156 runWarningIDC = currRun;
157 }
158 float scalerA = mTPCScaler.getMeanScaler(timestamp, o2::tpc::Side::A);
159 float scalerC = mTPCScaler.getMeanScaler(timestamp, o2::tpc::Side::C);
160 float meanScaler = (scalerA + scalerC) / 2;
161 tpcScaler = meanScaler;
162 if (mStreamer) {
163 (*mStreamer) << "treeIDC"
164 << "scalerA=" << scalerA
165 << "scalerC=" << scalerC
166 << "firstTFOrbit=" << firstTFOrbit
167 << "timestamp=" << timestamp
168 << "\n";
169 }
170 }
171 // check for Maps update
172 mTPCCorrMapsLoader.extractCCDBInputs(pc, tpcScaler);
173
174 if (mGlobOpts.requestCTPLumi) {
175 const float lumiCTP = mTPCCorrMapsLoader.getInstLumiCTP();
176 // if CTP lumi was notrequest - defualt of 0 is published, otherwise the value is scaled with the provided factor
177 LOGP(info, "Publishing CTP Lumi: {} for timestamp: {}, firstTFOrbit: {}", lumiCTP, timestamp, firstTFOrbit);
178 pc.outputs().snapshot(Output{header::gDataOriginCTP, "LUMICTP"}, lumiCTP);
179 }
180
181 buildMap(pc);
182 }
183
185 {
186 // reference map
187 auto* corrMap = mTPCCorrMapsLoader.getCorrMap();
188
189 // // new correction map
191 finalMap.cloneFromObject(*corrMap, nullptr);
192 finalMap.setApplyCorrectionOn();
193
194 const auto* corrMapRef = mTPCCorrMapsLoader.getCorrMapRef();
195 const float lumiScale = mTPCCorrMapsLoader.getLumiScale();
196 std::vector<std::pair<const o2::gpu::TPCFastSpaceChargeCorrection*, float>> additionalCorrections;
197
198 // if standard scaling is used: map(lumi) = (mean_map - ref_map) * lumiScale + ref_map
199 if (mTPCCorrMapsLoader.getLumiScaleMode() == LumiScaleMode::Linear) {
200 const std::vector<std::pair<const o2::gpu::TPCFastSpaceChargeCorrection*, float>> step0{{&(corrMapRef->getCorrection()), -1.f}};
201 // finalMap = (mean_map - finalMap)
203
204 // finalMap = finalMap * lumiScale + ref_map
205 const std::vector<std::pair<const o2::gpu::TPCFastSpaceChargeCorrection*, float>> step1{{&(corrMapRef->getCorrection()), 1.f}};
207
208 } else if (mTPCCorrMapsLoader.getLumiScaleMode() == LumiScaleMode::DerivativeMap || mTPCCorrMapsLoader.getLumiScaleMode() == LumiScaleMode::DerivativeMapMC) {
209 additionalCorrections.emplace_back(&(corrMapRef->getCorrection()), lumiScale);
210 }
211
212 // if mshape map valid
213 if (!mTPCCorrMapsLoader.isCorrMapMShapeDummy()) {
214 LOGP(info, "Adding M-shape correction to the final map with scaling factor {}", mMShapeScalingFac);
215 additionalCorrections.emplace_back(&(mTPCCorrMapsLoader.getCorrMapMShape()->getCorrection()), 1.f);
216 }
217
218 if (!additionalCorrections.empty()) {
219 TPCFastSpaceChargeCorrectionHelper::instance()->mergeCorrections(finalMap.getCorrection(), 1, additionalCorrections, true);
220 }
221
222 Output corrMapOutput{header::gDataOriginTPC, "TPCCORRMAP", 0};
223 auto outputBuffer = o2::pmr::vector<char>(pc.outputs().getMemoryResource(corrMapOutput));
224 outputBuffer.resize(TPCFastTransformPOD::estimateSize(finalMap.getCorrection()));
225 auto* pod = TPCFastTransformPOD::create(outputBuffer.data(), outputBuffer.size(), finalMap);
226 const auto& vd = mTPCVDriftHelper.getVDriftObject();
227 o2::tpc::TPCFastTransformHelperO2::instance()->updateCalibration(*pod, 0, vd.corrFact, vd.refVDrift, vd.getTimeOffset());
228 pc.outputs().adoptContainer(corrMapOutput, std::move(outputBuffer));
229 }
230
231 void finaliseCCDB(o2::framework::ConcreteDataMatcher& matcher, void* obj) final
232 {
234 mTPCVDriftHelper.accountCCDBInputs(matcher, obj);
235 if (matcher == ConcreteDataMatcher(o2::header::gDataOriginTPC, "TPCSCALERCCDB", 0)) {
236 LOGP(info, "Updating TPC scaler");
237 mTPCScaler.setFromTree(*((TTree*)obj));
238 if (mIonDriftTimeMS > 0) {
239 LOGP(info, "Setting ion drift time to: {}", mIonDriftTimeMS);
240 mTPCScaler.setIonDriftTimeMS(mIonDriftTimeMS);
241 }
242 if (mScalerWeights.isValid()) {
243 LOGP(info, "Setting TPC scaler weights");
244 mTPCScaler.setScalerWeights(mScalerWeights);
245 mTPCScaler.useWeights(true);
246 if (mIonDriftTimeMS == -1) {
247 overWriteIntegrationTime();
248 }
249 }
250 }
251 if (matcher == ConcreteDataMatcher(o2::header::gDataOriginTPC, "TPCSCALERWCCDB", 0)) {
252 LOGP(info, "Updating TPC scaler weights");
253 mScalerWeights = *(TPCScalerWeights*)obj;
254 mTPCScaler.setScalerWeights(mScalerWeights);
255 mTPCScaler.useWeights(true);
256 // in case ion drift time is not specified it is overwritten by the value in the weight object
257 if (mIonDriftTimeMS == -1) {
258 overWriteIntegrationTime();
259 }
260 }
261 if (matcher == ConcreteDataMatcher(o2::header::gDataOriginTPC, "MSHAPEPOTCCDB", 0)) {
262 LOGP(info, "Updating M-shape TPC scaler");
263 mMShapeTPCScaler.setFromTree(*((TTree*)obj));
264 if (mMShapeTPCScaler.getRun() == -1) {
265 LOGP(info, "Loaded default M-Shape correction object from CCDB");
266 }
267 }
268 if (mTPCCorrMapsLoader.accountCCDBInputs(matcher, obj)) {
269 return;
270 }
271 }
272
273 private:
274 std::shared_ptr<o2::base::GRPGeomRequest> mCCDBRequest;
275 const bool mEnableIDCs{true};
276 const bool mEnableMShape{false};
277 const o2::tpc::CorrectionMapsGloOpts mGlobOpts;
278 bool mEnableWeights{false};
279 TPCScalerWeights mScalerWeights{};
280 float mIonDriftTimeMS{-1};
281 float mMaxTimeWeightsMS{500};
282 TPCScaler mTPCScaler;
283 float mMShapeScalingFac{0};
284 TPCMShapeCorrection mMShapeTPCScaler;
285 int mKnotsYMshape{4};
286 int mKnotsZMshape{4};
287 std::unique_ptr<o2::utils::TreeStreamRedirector> mStreamer;
288 o2::tpc::CorrectionMapsLoader mTPCCorrMapsLoader{};
289 o2::tpc::VDriftHelper mTPCVDriftHelper{};
290
291 void overWriteIntegrationTime()
292 {
293 float integrationTime = std::abs(mScalerWeights.mFirstTimeStampMS);
294 if (integrationTime <= 0) {
295 return;
296 }
297 if (integrationTime > mMaxTimeWeightsMS) {
298 integrationTime = mMaxTimeWeightsMS;
299 }
300 LOGP(info, "Setting maximum integration time for weights to: {}", integrationTime);
301 mTPCScaler.setIonDriftTimeMS(integrationTime);
302 }
303};
304
306{
307 std::vector<InputSpec> inputs;
308 if (enableIDCs) {
309 LOGP(info, "Publishing IDC scalers for space-charge distortion fluctuation correction");
310 inputs.emplace_back("tpcscaler", o2::header::gDataOriginTPC, "TPCSCALERCCDB", 0, Lifetime::Condition, ccdbParamSpec(o2::tpc::CDBTypeMap.at(o2::tpc::CDBType::CalScaler), {}, 1)); // time-dependent
311 inputs.emplace_back("tpcscalerw", o2::header::gDataOriginTPC, "TPCSCALERWCCDB", 0, Lifetime::Condition, ccdbParamSpec(o2::tpc::CDBTypeMap.at(o2::tpc::CDBType::CalScalerWeights), {}, 0)); // non time-dependent
312 }
313 if (enableMShape) {
314 LOGP(info, "Publishing M-shape correction map");
315 inputs.emplace_back("mshape", o2::header::gDataOriginTPC, "MSHAPEPOTCCDB", 0, Lifetime::Condition, ccdbParamSpec(o2::tpc::CDBTypeMap.at(o2::tpc::CDBType::CalMShape), {}, 1)); // time-dependent
316 }
317
318 auto ccdbRequest = std::make_shared<o2::base::GRPGeomRequest>(true, // orbitResetTime
319 false, // GRPECS=true for nHBF per TF
320 false, // GRPLHCIF
321 false, // GRPMagField
322 false, // askMatLUT
324 inputs);
325
326 std::vector<OutputSpec> outputs;
327 outputs.emplace_back(o2::header::gDataOriginTPC, "TPCCORRMAP", 0, Lifetime::Timeframe);
328 if (sclOpts.requestCTPLumi) {
329 outputs.emplace_back(o2::header::gDataOriginCTP, "LUMICTP", 0, Lifetime::Timeframe);
330 }
333
334 return DataProcessorSpec{
335 "tpc-scaler",
336 inputs,
337 outputs,
338 AlgorithmSpec{adaptFromTask<TPCScalerSpec>(ccdbRequest, sclOpts, enableIDCs, enableMShape)},
339 Options{
340 {"ion-drift-time", VariantType::Float, -1.f, {"Overwrite ion drift time if a value >0 is provided"}},
341 {"max-time-for-weights", VariantType::Float, 500.f, {"Maximum possible integration time in ms when weights are used"}},
342 {"m-shape-scaling-factor", VariantType::Float, 1.f, {"Scale M-shape scaler with this value"}},
343 {"disableWeights", VariantType::Bool, false, {"Disable weights for TPC scalers"}},
344 {"enableStreamer", VariantType::Bool, false, {"Enable streaming of M-shape scalers"}},
345 {"n-threads", VariantType::Int, 4, {"Number of threads used for the M-shape correction"}},
346 {"n-knots-y", VariantType::Int, 4, {"Number of knots in y-direction used for the M-shape correction"}},
347 {"n-knots-z", VariantType::Int, 4, {"Number of knots in z-direction used for the M-shape correction"}}}};
348}
349
350} // namespace tpc
351} // namespace o2
Simple interface to the CDB manager.
Helper class to access load maps from CCDB.
Helper for geometry and GRP related CCDB requests.
uint32_t roc
Definition RawData.h:3
uint32_t side
Definition RawData.h:0
This class contains the algorithms for calculation the distortions and corrections.
class to create TPC fast space charge correction
Helper class to extract VDrift from different sources.
auto getOrbitResetTimeMS() const
void checkUpdates(o2::framework::ProcessingContext &pc)
bool finaliseCCDB(o2::framework::ConcreteDataMatcher &matcher, void *obj)
static GRPGeomHelper & instance()
void setRequest(std::shared_ptr< GRPGeomRequest > req)
o2::pmr::FairMQMemoryResource * getMemoryResource(const Output &spec)
CacheId adoptContainer(const Output &, ContainerT &, CacheStrategy, o2::header::SerializationMethod)
DataAllocator & outputs()
The data allocator is used to allocate memory for the output data.
void setLumiScaleType(tpc::LumiScaleType v)
void setLumiScaleMode(tpc::LumiScaleMode v)
const o2::gpu::TPCFastTransform * getCorrMapMShape() const
const o2::gpu::TPCFastTransform * getCorrMap() const
void setCorrMapMShape(std::unique_ptr< o2::gpu::TPCFastTransform > &&m)
bool isCorrMapMShapeDummy() const
return returns if the correction map for the M-shape correction is a dummy spline object
const o2::gpu::TPCFastTransform * getCorrMapRef() const
tpc::LumiScaleMode getLumiScaleMode() const
static TPCFastTransformPOD * create(aligned_unique_buffer_ptr< TPCFastTransformPOD > &destVector, const TPCFastTransform &src)
Create POD transform from old flat-buffer one. Provided vector will serve as a buffer.
static size_t estimateSize(const TPCFastTransform &src)
void cloneFromObject(const TPCFastTransform &obj, char *newFlatBufferPtr)
Construction interface.
TPCFastSpaceChargeCorrection & getCorrection()
Gives a reference for external initialization of TPC corrections.
bool accountCCDBInputs(const o2::framework::ConcreteDataMatcher &matcher, void *obj)
void extractCCDBInputs(o2::framework::ProcessingContext &pc, float tpcScaler=-1.f)
void init(o2::framework::InitContext &ic, bool idcsAvailable)
static void requestCCDBInputs(std::vector< o2::framework::InputSpec > &inputs, const o2::tpc::CorrectionMapsGloOpts &gloOpts)
static void setNThreads(const int nThreads)
set the number of threads used for some of the calculations
void getCorrections(const DataT x, const DataT y, const DataT z, const Side side, DataT &corrX, DataT &corrY, DataT &corrZ) const
NumericalFields< DataT > getElectricFieldsInterpolator(const Side side) const
void calcEField(const Side side)
void setPotential(int iz, int ir, int iphi, Side side, float val)
setting the potential directly for given vertex
void calcGlobalCorrections(const Fields &formulaStruct, const int type=3)
void poissonSolver(const Side side, const DataT stoppingConvergence=1e-6, const int symmetry=0)
std::unique_ptr< TPCFastSpaceChargeCorrection > createFromGlobalCorrection(std::function< void(int roc, double gx, double gy, double gz, double &dgx, double &dgy, double &dgz)> correctionGlobal, const int nKnotsY=10, const int nKnotsZ=20)
creates TPCFastSpaceChargeCorrection object from a continious space charge correction in global coord...
void setNthreads(int n)
_______________ Settings ________________________
void mergeCorrections(o2::gpu::TPCFastSpaceChargeCorrection &mainCorrection, float scale, const std::vector< std::pair< const o2::gpu::TPCFastSpaceChargeCorrection *, float > > &additionalCorrections, bool prn)
static TPCFastSpaceChargeCorrectionHelper * instance()
Singleton.
static TPCFastTransformHelperO2 * instance()
Singleton.
int updateCalibration(TPCFastTransform &fastTransform, int64_t TimeStamp, float vDriftFactor=1.f, float vDriftRef=0.f, float driftTimeOffset=0.f)
Updates the transformation with the new time stamp.
BoundaryPotentialIFC getBoundaryPotential(const double timestamp) const
void setFromTree(TTree &tpcMShapeTree)
set this object from input tree
void endOfStream(EndOfStreamContext &eos) final
TPCScalerSpec(std::shared_ptr< o2::base::GRPGeomRequest > req, const o2::tpc::CorrectionMapsGloOpts &sclOpts, bool enableIDCs, bool enableMShape)
void finaliseCCDB(o2::framework::ConcreteDataMatcher &matcher, void *obj) final
void run(ProcessingContext &pc) final
void init(framework::InitContext &ic) final
void buildMap(ProcessingContext &pc)
float getMeanScaler(double timestamp, o2::tpc::Side side) const
void setIonDriftTimeMS(float ionDriftTimeMS)
Definition TPCScaler.h:58
int getRun() const
Definition TPCScaler.h:106
void useWeights(bool useWeights)
enable usage of weights
Definition TPCScaler.h:142
void setFromTree(TTree &tpcScalerTree)
set this object from input tree
void setScalerWeights(const TPCScalerWeights &weights)
setting the weights for the scalers
Definition TPCScaler.h:136
static void requestCCDBInputs(std::vector< o2::framework::InputSpec > &inputs, bool laser=true, bool itstpcTgl=true)
void extractCCDBInputs(o2::framework::ProcessingContext &pc, bool laser=true, bool itstpcTgl=true)
const VDriftCorrFact & getVDriftObject() const
bool accountCCDBInputs(const o2::framework::ConcreteDataMatcher &matcher, void *obj)
bool isUpdated() const
GLint GLenum GLint x
Definition glcorearb.h:403
GLint y
Definition glcorearb.h:270
typedef void(APIENTRYP PFNGLCULLFACEPROC)(GLenum mode)
GLdouble GLdouble GLdouble z
Definition glcorearb.h:843
constexpr o2::header::DataOrigin gDataOriginCTP
Definition DataHeader.h:564
constexpr o2::header::DataOrigin gDataOriginTPC
Definition DataHeader.h:576
constexpr double LHCOrbitMUS
Defining ITS Vertex explicitly as messageable.
Definition Cartesian.h:288
std::vector< ConfigParamSpec > ccdbParamSpec(std::string const &path, int runDependent, std::vector< CCDBMetadata > metadata={}, int qrate=0)
std::vector< ConfigParamSpec > Options
std::vector< T, fair::mq::pmr::polymorphic_allocator< T > > vector
@ DerivativeMap
map(lumi) = mean_map + lumiScale * (derivativeMap) where derivativeMap = (mean_map_A - mean_map_B)
@ Linear
map(lumi) = (mean_map - referenceMap) * lumiScale + referenceMap
@ DerivativeMapMC
same DerivativeMap, but for MC
const std::unordered_map< CDBType, const std::string > CDBTypeMap
Storage name in CCDB for each calibration and parameter type.
Definition CDBTypes.h:96
@ TPCScaler
use TPC scaler for scaling
Side
TPC readout sidE.
Definition Defs.h:35
@ A
Definition Defs.h:35
@ C
Definition Defs.h:36
o2::framework::DataProcessorSpec getTPCScalerSpec(bool enableIDCs, bool enableMShape, const o2::tpc::CorrectionMapsGloOpts &sclOpts)
@ CalScaler
Scaler from IDCs or combined estimator.
@ CalScalerWeights
Weights for scalers.
a couple of static helper functions to create timestamp values for CCDB queries or override obsolete ...
LumiScaleType lumiType
what estimator to used for corrections scaling: 0: no scaling, 1: CTP, 2: IDC
bool checkCTPIDCconsistency
check the selected CTP or IDC scaling source being consistent with mean scaler of the map
bool requestCTPLumi
request CTP Lumi regardless of what is used for corrections scaling
LumiScaleMode lumiMode
what corrections method to use: 0: classical scaling, 1: Using of the derivative map,...
float mFirstTimeStampMS
first timestamp
Definition TPCScaler.h:35