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"
30
31using namespace o2::framework;
32
33namespace o2
34{
35namespace tpc
36{
37
38class TPCScalerSpec : public Task
39{
40 public:
41 TPCScalerSpec(std::shared_ptr<o2::base::GRPGeomRequest> req, bool enableIDCs, bool enableMShape) : mCCDBRequest(req), mEnableIDCs(enableIDCs), mEnableMShape(enableMShape){};
42
44 {
46 mIonDriftTimeMS = ic.options().get<float>("ion-drift-time");
47 mMaxTimeWeightsMS = ic.options().get<float>("max-time-for-weights");
48 mMShapeScalingFac = ic.options().get<float>("m-shape-scaling-factor");
49 mEnableWeights = !(ic.options().get<bool>("disableWeights"));
50 const bool enableStreamer = ic.options().get<bool>("enableStreamer");
51 const int mshapeThreads = ic.options().get<int>("n-threads");
52 mKnotsYMshape = ic.options().get<int>("n-knots-y");
53 mKnotsZMshape = ic.options().get<int>("n-knots-z");
56
57 if (enableStreamer) {
58 mStreamer = std::make_unique<o2::utils::TreeStreamRedirector>("M_Shape.root", "recreate");
59 }
60 }
61
63 {
64 if (mStreamer) {
65 mStreamer->Close();
66 }
67 }
68
69 void run(ProcessingContext& pc) final
70 {
72 if (mEnableIDCs && pc.inputs().isValid("tpcscaler")) {
73 pc.inputs().get<TTree*>("tpcscaler");
74 }
75
76 if (mEnableWeights && mEnableIDCs) {
77 if (pc.inputs().isValid("tpcscalerw")) {
78 pc.inputs().get<TPCScalerWeights*>("tpcscalerw");
79 }
80 }
81 if (mEnableMShape && pc.inputs().isValid("mshape")) {
82 pc.inputs().get<TTree*>("mshape");
83 }
84
85 const auto orbitResetTimeMS = o2::base::GRPGeomHelper::instance().getOrbitResetTimeMS();
86 const auto firstTFOrbit = pc.services().get<o2::framework::TimingInfo>().firstTForbit;
87 const double timestamp = orbitResetTimeMS + firstTFOrbit * o2::constants::lhc::LHCOrbitMUS * 0.001;
88 int currRun = pc.services().get<o2::framework::TimingInfo>().runNumber;
89 if (mEnableMShape) {
90 static int runWarningMS = -1;
91 if ((mMShapeTPCScaler.getRun() != -1) && currRun != mMShapeTPCScaler.getRun() && runWarningMS != currRun) {
92 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());
93 runWarningMS = currRun;
94 }
95
96 const auto& boundaryPotential = mMShapeTPCScaler.getBoundaryPotential(timestamp);
97 if (!boundaryPotential.mPotential.empty()) {
98 LOGP(info, "Calculating M-shape correction from input boundary potential");
99 const Side side = Side::A;
100 o2::tpc::SpaceCharge<double> sc = o2::tpc::SpaceCharge<double>(mMShapeTPCScaler.getMShapes().mBField, mMShapeTPCScaler.getMShapes().mZ, mMShapeTPCScaler.getMShapes().mR, mMShapeTPCScaler.getMShapes().mPhi);
101 for (int iz = 0; iz < mMShapeTPCScaler.getMShapes().mZ; ++iz) {
102 for (int iphi = 0; iphi < mMShapeTPCScaler.getMShapes().mPhi; ++iphi) {
103 const float pot = mMShapeScalingFac * boundaryPotential.mPotential[iz];
104 sc.setPotential(iz, 0, iphi, side, pot);
105 }
106 }
107
109 sc.calcEField(side);
111
112 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) {
113 Side side = roc < 18 ? Side::A : Side::C;
114 if (side == Side::C) {
115 dx = 0;
116 dy = 0;
117 dz = 0;
118 } else {
119 sc.getCorrections(x, y, z, side, dx, dy, dz);
120 }
121 };
122
123 std::unique_ptr<TPCFastSpaceChargeCorrection> spCorrection = TPCFastSpaceChargeCorrectionHelper::instance()->createFromGlobalCorrection(getCorrections, mKnotsYMshape, mKnotsZMshape);
124 std::unique_ptr<TPCFastTransform> fastTransform(TPCFastTransformHelperO2::instance()->create(0, *spCorrection));
125 pc.outputs().snapshot(Output{header::gDataOriginTPC, "TPCMSHAPE"}, *fastTransform);
126 } else {
127 // send empty dummy object
128 LOGP(info, "Sending default (no) M-shape correction");
129 auto fastTransform = o2::tpc::TPCFastTransformHelperO2::instance()->create(0);
130 pc.outputs().snapshot(Output{header::gDataOriginTPC, "TPCMSHAPE"}, *fastTransform);
131 }
132
133 if (mStreamer) {
134 (*mStreamer) << "treeMShape"
135 << "firstTFOrbit=" << firstTFOrbit
136 << "timestamp=" << timestamp
137 << "boundaryPotential=" << boundaryPotential
138 << "mMShapeScalingFac=" << mMShapeScalingFac
139 << "\n";
140 }
141 }
142
143 if (mEnableIDCs) {
144 static int runWarningIDC = -1;
145 if (pc.services().get<o2::framework::TimingInfo>().runNumber != mTPCScaler.getRun() && runWarningIDC != currRun) {
146 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());
147 runWarningIDC = currRun;
148 }
149 float scalerA = mTPCScaler.getMeanScaler(timestamp, o2::tpc::Side::A);
150 float scalerC = mTPCScaler.getMeanScaler(timestamp, o2::tpc::Side::C);
151 float meanScaler = (scalerA + scalerC) / 2;
152 LOGP(info, "Publishing TPC scaler: {} for timestamp: {}, firstTFOrbit: {}", meanScaler, timestamp, firstTFOrbit);
153 pc.outputs().snapshot(Output{header::gDataOriginTPC, "TPCSCALER"}, meanScaler);
154 if (mStreamer) {
155 (*mStreamer) << "treeIDC"
156 << "scalerA=" << scalerA
157 << "scalerC=" << scalerC
158 << "firstTFOrbit=" << firstTFOrbit
159 << "timestamp=" << timestamp
160 << "\n";
161 }
162 }
163 }
164
165 void finaliseCCDB(o2::framework::ConcreteDataMatcher& matcher, void* obj) final
166 {
168 if (matcher == ConcreteDataMatcher(o2::header::gDataOriginTPC, "TPCSCALERCCDB", 0)) {
169 LOGP(info, "Updating TPC scaler");
170 mTPCScaler.setFromTree(*((TTree*)obj));
171 if (mIonDriftTimeMS > 0) {
172 LOGP(info, "Setting ion drift time to: {}", mIonDriftTimeMS);
173 mTPCScaler.setIonDriftTimeMS(mIonDriftTimeMS);
174 }
175 if (mScalerWeights.isValid()) {
176 LOGP(info, "Setting TPC scaler weights");
177 mTPCScaler.setScalerWeights(mScalerWeights);
178 mTPCScaler.useWeights(true);
179 if (mIonDriftTimeMS == -1) {
180 overWriteIntegrationTime();
181 }
182 }
183 }
184 if (matcher == ConcreteDataMatcher(o2::header::gDataOriginTPC, "TPCSCALERWCCDB", 0)) {
185 LOGP(info, "Updating TPC scaler weights");
186 mScalerWeights = *(TPCScalerWeights*)obj;
187 mTPCScaler.setScalerWeights(mScalerWeights);
188 mTPCScaler.useWeights(true);
189 // in case ion drift time is not specified it is overwritten by the value in the weight object
190 if (mIonDriftTimeMS == -1) {
191 overWriteIntegrationTime();
192 }
193 }
194 if (matcher == ConcreteDataMatcher(o2::header::gDataOriginTPC, "MSHAPEPOTCCDB", 0)) {
195 LOGP(info, "Updating M-shape TPC scaler");
196 mMShapeTPCScaler.setFromTree(*((TTree*)obj));
197 if (mMShapeTPCScaler.getRun() == -1) {
198 LOGP(info, "Loaded default M-Shape correction object from CCDB");
199 }
200 }
201 }
202
203 private:
204 std::shared_ptr<o2::base::GRPGeomRequest> mCCDBRequest;
205 const bool mEnableIDCs{true};
206 const bool mEnableMShape{false};
207 bool mEnableWeights{false};
208 TPCScalerWeights mScalerWeights{};
209 float mIonDriftTimeMS{-1};
210 float mMaxTimeWeightsMS{500};
211 TPCScaler mTPCScaler;
212 float mMShapeScalingFac{0};
213 TPCMShapeCorrection mMShapeTPCScaler;
214 int mKnotsYMshape{4};
215 int mKnotsZMshape{4};
216 std::unique_ptr<o2::utils::TreeStreamRedirector> mStreamer;
217
218 void overWriteIntegrationTime()
219 {
220 float integrationTime = std::abs(mScalerWeights.mFirstTimeStampMS);
221 if (integrationTime <= 0) {
222 return;
223 }
224 if (integrationTime > mMaxTimeWeightsMS) {
225 integrationTime = mMaxTimeWeightsMS;
226 }
227 LOGP(info, "Setting maximum integration time for weights to: {}", integrationTime);
228 mTPCScaler.setIonDriftTimeMS(integrationTime);
229 }
230};
231
232o2::framework::DataProcessorSpec getTPCScalerSpec(bool enableIDCs, bool enableMShape)
233{
234 std::vector<InputSpec> inputs;
235 if (enableIDCs) {
236 LOGP(info, "Publishing IDC scalers for space-charge distortion fluctuation correction");
237 inputs.emplace_back("tpcscaler", o2::header::gDataOriginTPC, "TPCSCALERCCDB", 0, Lifetime::Condition, ccdbParamSpec(o2::tpc::CDBTypeMap.at(o2::tpc::CDBType::CalScaler), {}, 1)); // time-dependent
238 inputs.emplace_back("tpcscalerw", o2::header::gDataOriginTPC, "TPCSCALERWCCDB", 0, Lifetime::Condition, ccdbParamSpec(o2::tpc::CDBTypeMap.at(o2::tpc::CDBType::CalScalerWeights), {}, 0)); // non time-dependent
239 }
240 if (enableMShape) {
241 LOGP(info, "Publishing M-shape correction map");
242 inputs.emplace_back("mshape", o2::header::gDataOriginTPC, "MSHAPEPOTCCDB", 0, Lifetime::Condition, ccdbParamSpec(o2::tpc::CDBTypeMap.at(o2::tpc::CDBType::CalMShape), {}, 1)); // time-dependent
243 }
244
245 auto ccdbRequest = std::make_shared<o2::base::GRPGeomRequest>(true, // orbitResetTime
246 false, // GRPECS=true for nHBF per TF
247 false, // GRPLHCIF
248 false, // GRPMagField
249 false, // askMatLUT
251 inputs);
252
253 std::vector<OutputSpec> outputs;
254 if (enableIDCs) {
255 outputs.emplace_back(o2::header::gDataOriginTPC, "TPCSCALER", 0, Lifetime::Timeframe);
256 }
257 if (enableMShape) {
258 outputs.emplace_back(o2::header::gDataOriginTPC, "TPCMSHAPE", 0, Lifetime::Timeframe);
259 }
260
261 return DataProcessorSpec{
262 "tpc-scaler",
263 inputs,
264 outputs,
265 AlgorithmSpec{adaptFromTask<TPCScalerSpec>(ccdbRequest, enableIDCs, enableMShape)},
266 Options{
267 {"ion-drift-time", VariantType::Float, -1.f, {"Overwrite ion drift time if a value >0 is provided"}},
268 {"max-time-for-weights", VariantType::Float, 500.f, {"Maximum possible integration time in ms when weights are used"}},
269 {"m-shape-scaling-factor", VariantType::Float, 1.f, {"Scale M-shape scaler with this value"}},
270 {"disableWeights", VariantType::Bool, false, {"Disable weights for TPC scalers"}},
271 {"enableStreamer", VariantType::Bool, false, {"Enable streaming of M-shape scalers"}},
272 {"n-threads", VariantType::Int, 4, {"Number of threads used for the M-shape correction"}},
273 {"n-knots-y", VariantType::Int, 4, {"Number of knots in y-direction used for the M-shape correction"}},
274 {"n-knots-z", VariantType::Int, 4, {"Number of knots in z-direction used for the M-shape correction"}}}};
275}
276
277} // namespace tpc
278} // namespace o2
Simple interface to the CDB manager.
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
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)
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 ________________________
static TPCFastSpaceChargeCorrectionHelper * instance()
Singleton.
std::unique_ptr< TPCFastTransform > create(Long_t TimeStamp)
_______________ Main functionality ________________________
static TPCFastTransformHelperO2 * instance()
Singleton.
BoundaryPotentialIFC getBoundaryPotential(const double timestamp) const
void setFromTree(TTree &tpcMShapeTree)
set this object from input tree
void endOfStream(EndOfStreamContext &eos) final
void finaliseCCDB(o2::framework::ConcreteDataMatcher &matcher, void *obj) final
TPCScalerSpec(std::shared_ptr< o2::base::GRPGeomRequest > req, bool enableIDCs, bool enableMShape)
void run(ProcessingContext &pc) final
void init(framework::InitContext &ic) final
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
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 gDataOriginTPC
Definition DataHeader.h:576
constexpr double LHCOrbitMUS
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)
std::vector< ConfigParamSpec > Options
o2::framework::DataProcessorSpec getTPCScalerSpec(bool enableIDCs, bool enableMShape)
const std::unordered_map< CDBType, const std::string > CDBTypeMap
Storage name in CCDB for each calibration and parameter type.
Definition CDBTypes.h:94
Side
TPC readout sidE.
Definition Defs.h:35
@ A
Definition Defs.h:35
@ C
Definition Defs.h:36
@ 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 ...
float mFirstTimeStampMS
first timestamp
Definition TPCScaler.h:35