Project
Loading...
Searching...
No Matches
TPCVDriftTglCalibration.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
13#include "Framework/Logger.h"
14#include "MathUtils/fit.h"
16#include "CCDB/CcdbApi.h"
18#include <Math/SMatrix.h>
19#include <TH2F.h>
20
21namespace o2
22{
23namespace tpc
24{
25
28
31
32//_____________________________________________
34{
35 // Here we initialize the vector of our output objects
36 mVDPerSlot.clear();
37 mCCDBInfoPerSlot.clear();
38 return;
39}
40
41//_____________________________________________
43{
44 // Extract results for the single slot
45 const auto* cont = slot.getContainer();
46 const auto* h = cont->histo.get();
47 auto nx = h->getNBinsX(), ny = h->getNBinsY();
48 std::array<double, 3> parg;
50 double sS = 0, sX = 0, sY = 0, sXX = 0, sXY = 0;
51 int npval = 0;
52 for (auto ix = nx; ix--;) {
53 const auto sliceY = h->getSliceY(ix);
54 auto chi2 = o2::math_utils::fitGaus(ny, sliceY.data(), h->getYMin(), h->getYMax(), parg, &mat);
55 if (chi2 < -2) { // failed
56 continue;
57 }
58 double err2 = mat(1, 1);
59 if (err2 <= 0.f) {
60 continue;
61 }
62 double err2i = 1. / err2, xb = h->getBinXCenter(ix), xbw = xb * err2i, ybw = parg[1] * err2i;
63 sS += err2i;
64 sX += xbw;
65 sXX += xb * xbw;
66 sY += ybw;
67 sXY += xb * ybw;
68 npval++;
69 }
70 if (!mSaveHistosFile.empty()) {
71 TFile savf(mSaveHistosFile.c_str(), "update");
72 auto th2f = h->createTH2F(fmt::format("vdtgl{}_{}", slot.getTFStart(), slot.getTFEnd()));
73 th2f->Write();
74 LOGP(info, "Saved histo for slot {}-{} to {}", slot.getTFStart(), slot.getTFEnd(), mSaveHistosFile);
75 }
76 double det = sS * sXX - sX * sX;
77 if (!det || npval < 2) {
78 LOGP(alarm, "VDrift fit failed for slot {}<=TF<={} with {} entries: det={} npoints={}", slot.getTFStart(), slot.getTFEnd(), cont->entries, det, npval);
79 } else {
80 det = 1. / det;
81 double offs = (sXX * sY - sX * sXY) * det, slope = (sS * sXY - sX * sY) * det;
82 double offsErr = sXX * det, slopErr = sS * det;
83 offsErr = offsErr > 0. ? std::sqrt(offsErr) : 0.;
84 slopErr = slopErr > 0. ? std::sqrt(slopErr) : 0.;
85 float corrFact = 1. / (1. - slope);
86 float corrFactErr = corrFact * corrFact * slopErr;
87 auto& vd = mVDPerSlot.emplace_back(o2::tpc::VDriftCorrFact{slot.getStartTimeMS(),
88 slot.getEndTimeMS(),
89 std::chrono::duration_cast<std::chrono::milliseconds>(std::chrono::system_clock::now().time_since_epoch()).count(),
90 corrFact,
91 corrFactErr,
92 float(cont->driftVFullMean),
93 cont->tOffsetRef, 0.f, cont->tp});
94 // at this stage the correction object is defined wrt average corrected drift used for the slot processing, we want to redefine it to run-constant reference vdrift
95 vd.normalize(cont->driftVRef);
96
97 auto clName = o2::utils::MemFileHelper::getClassName(mVDPerSlot.back());
98 auto flName = o2::ccdb::CcdbApi::generateFileName(clName);
99 std::map<std::string, std::string> md;
100 mCCDBInfoPerSlot.emplace_back("TPC/Calib/VDriftTgl", clName, flName, md,
101 vd.firstTime - 10 * o2::ccdb::CcdbObjectInfo::SECOND, vd.lastTime + o2::ccdb::CcdbObjectInfo::MONTH);
102 LOGP(info, "Finalize slot {}({})<=TF<={}({}) with {} entries | dTgl vs Tgl_ITS offset: {:.4f}+-{:.4f} Slope: {:.4f}+-{:.4f} -> Corr factor = {:.4f}+-{:.4f} wrt <VD>={:.4f} -> {:.4f} wrt {:.4f}",
103 slot.getTFStart(), slot.getStartTimeMS(), slot.getTFEnd(), slot.getEndTimeMS(), cont->entries, offs, offsErr, slope, slopErr, corrFact, corrFactErr, cont->driftVFullMean, vd.corrFact, vd.refVDrift);
104 }
105}
106
107//_____________________________________________
109{
110 auto& cont = getSlots();
111 auto& slot = front ? cont.emplace_front(tstart, tend) : cont.emplace_back(tstart, tend);
112 slot.setContainer(std::make_unique<TPCVDTglContainer>(mNBinsTgl, mMaxTgl, mNBinsDTgl, mMaxDTgl));
113 return slot;
114}
115
116} // namespace tpc
117} // end namespace o2
Utils and constants for calibration and related workflows.
uint16_t slope
Definition RawData.h:1
Class for time synchronization of RawReader instances.
TFType getTFEnd() const
Definition TimeSlot.h:47
long getStartTimeMS() const
Definition TimeSlot.h:50
long getEndTimeMS() const
Definition TimeSlot.h:51
const Container * getContainer() const
Definition TimeSlot.h:53
TFType getTFStart() const
Definition TimeSlot.h:46
static std::string generateFileName(const std::string &inp)
Definition CcdbApi.cxx:816
static constexpr long MONTH
static constexpr long SECOND
Slot & emplaceNewSlot(bool front, TFType tstart, TFType tend) final
Double_t fitGaus(const size_t nBins, const T *arr, const T xMin, const T xMax, std::vector< T > &param)
Definition fit.h:232
a couple of static helper functions to create timestamp values for CCDB queries or override obsolete ...
static std::string getClassName(const T &obj)
get the class name of the object