Project
Loading...
Searching...
No Matches
CalibratorGain.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
15
22#include "TStopwatch.h"
23#include "CCDB/CcdbApi.h"
27#include <TFile.h>
28#include <TTree.h>
29#include <TF1Convolution.h>
30
31#include <string>
32#include <map>
33#include <memory>
34#include <ctime>
35
36using namespace o2::trd::constants;
37
38namespace o2::trd
39{
40
42
44{
45 // reset the CCDB output vectors
46 mInfoVector.clear();
47 mObjectVector.clear();
48}
49
51{
52 if (mInitDone) {
53 return;
54 }
55 LOG(info) << "Initializing the processing";
56 for (int iDet = 0; iDet < MAXCHAMBER; ++iDet) {
57 mdEdxhists[iDet] = std::make_unique<TH1F>(Form("hdEdx%d", iDet), "dEdx", NBINSGAINCALIB, 0., NBINSGAINCALIB);
58 }
59
60 // Init fitting functions
61 mFconv = std::make_unique<TF1Convolution>("[3] * TMath::Landau(x, [0], [1]) * exp(-[2]*x)", "TMath::Gaus(x, 0, [0])", 0., NBINSGAINCALIB, true);
62 mFconv->SetNofPointsFFT(2000);
63 mFitFunction = std::make_unique<TF1>("fitConvLandau", *mFconv, 0., NBINSGAINCALIB, mFconv->GetNpar());
64 mFitFunction->SetParameters(40, 15, 0.02, 1, 0.1);
65 mFitFunction->SetParLimits(0, 30., 100.0);
66 mFitFunction->SetParLimits(1, 5.0, 25.0);
67 mFitFunction->SetParLimits(2, -0.1, 0.5);
68 mFitFunction->SetParLimits(3, 1.0, 10.0);
69 mFitFunction->SetParLimits(4, -0.1, 0.5);
70
71 // set tree addresses
72 if (mEnableOutput) {
73 mOutTree->Branch("MPVdEdx", &mFitResults);
74 for (int iDet = 0; iDet < MAXCHAMBER; ++iDet) {
75 mOutTree->Branch(Form("dEdxHist_%d", iDet), mdEdxhists[iDet].get());
76 }
77 }
78
79 mInitDone = true;
80}
81
83{
84
85 mFitResults.fill(constants::MPVDEDXDEFAULT);
86 // We either get a pointer to a valid object from the last ~hour or to the default object
87 // which is always present. The first has precedence over the latter.
88 auto dataCalGain = pc.inputs().get<o2::trd::CalGain*>("calgain");
89 std::string msg = "Default Object";
90 // We check if the object we got is the default one by comparing it to the defaults.
91 for (int iDet = 0; iDet < MAXCHAMBER; ++iDet) {
92 if (dataCalGain->getMPVdEdx(iDet) != constants::MPVDEDXDEFAULT) {
93 msg = "Previous Object";
94 break;
95 }
96 }
97 LOG(info) << "Calibrator: From CCDB retrieved " << msg;
98
99 // Here we set each entry regardless if it is the default or not.
100 for (int iDet = 0; iDet < MAXCHAMBER; ++iDet) {
101 mFitResults[iDet] = dataCalGain->getMPVdEdx(iDet);
102 }
103}
104
106{
107 LOG(info) << "Finalizing gain calibration";
108 print(); // to see current number of slots and their entries
109 // do actual calibration for the data provided in the given slot
110 TStopwatch timer;
111 timer.Start();
113
114 auto dEdxHists = slot.getContainer();
115 LOGP(info, "Current slot has {} entries", dEdxHists->getNEntries());
116 for (int iDet = 0; iDet < MAXCHAMBER; ++iDet) {
117 mdEdxhists[iDet]->Reset();
118 int nEntries = 0;
119 for (int iBin = 0; iBin < NBINSGAINCALIB; ++iBin) {
120 nEntries += dEdxHists->getHistogramEntry(iDet * NBINSGAINCALIB + iBin);
121 mdEdxhists[iDet]->SetBinContent(iBin + 1, dEdxHists->getHistogramEntry(iDet * NBINSGAINCALIB + iBin));
122 mdEdxhists[iDet]->SetBinError(iBin + 1, sqrt(dEdxHists->getHistogramEntry(iDet * NBINSGAINCALIB + iBin)));
123 }
124 // Check if we have the minimum amount of entries
125 if (nEntries < mMinEntriesChamber) {
126 LOGF(debug, "Chamber %d did not reach minimum amount of %d entries for refit", iDet, mMinEntriesChamber);
127 continue;
128 }
129 mdEdxhists[iDet]->Scale(1. / nEntries);
130
131 // Fitting histogram
132 mFitFunction->SetParameter(0, mdEdxhists[iDet]->GetMean() / 1.25);
133 int fitStatus = mdEdxhists[iDet]->Fit("fitConvLandau", "LQB0", "", 1, NBINSGAINCALIB - 4);
134
135 if (fitStatus != 0) {
136 LOGF(warn, "Fit for chamber %i failed, nEntries: %d", iDet, nEntries);
137 continue;
138 }
139
140 mFitResults[iDet] = mFitFunction->GetMaximumX(0., NBINSGAINCALIB);
141 LOGF(debug, "Fit result for chamber %i: dEdx MPV = ", iDet, mFitResults[iDet]);
142 }
143 timer.Stop();
144 LOGF(info, "Done fitting dEdx histograms. CPU time: %f, real time: %f", timer.CpuTime(), timer.RealTime());
145
146 // Fill Tree and log to debug
147 if (mEnableOutput) {
148 mOutTree->Fill();
149 for (int iDet = 0; iDet < MAXCHAMBER; ++iDet) {
150 LOGF(debug, "Fit result for chamber %i: dEdx MPV = ", iDet, mFitResults[iDet]);
151 }
152 }
153
154 // assemble CCDB object
155 CalGain calObject;
156 for (int iDet = 0; iDet < MAXCHAMBER; ++iDet) {
157 // For Chambers which did not have the minimum amount of entries in this slot e.g. missing, empty chambers.
158 // We will reuse the previous one. This would have been read either from the ccdb or come from a previous successful fit.
159 calObject.setMPVdEdx(iDet, mFitResults[iDet]);
160 }
161 auto clName = o2::utils::MemFileHelper::getClassName(calObject);
162 auto flName = o2::ccdb::CcdbApi::generateFileName(clName);
163 std::map<std::string, std::string> metadata; // TODO: do we want to store any meta data?
164 long startValidity = slot.getStartTimeMS();
165 mInfoVector.emplace_back("TRD/Calib/CalGain", clName, flName, metadata, startValidity, startValidity + 3 * o2::ccdb::CcdbObjectInfo::DAY);
166 mObjectVector.push_back(calObject);
167}
168
170{
171 auto& container = getSlots();
172 auto& slot = front ? container.emplace_front(tStart, tEnd) : container.emplace_back(tStart, tEnd);
173 slot.setContainer(std::make_unique<GainCalibHistos>());
174 slot.getContainer()->init();
175 return slot;
176}
177
179{
180 mEnableOutput = true;
181 mOutFile = std::make_unique<TFile>("trd_calibgain.root", "RECREATE");
182 if (mOutFile->IsZombie()) {
183 LOG(error) << "Failed to create output file!";
184 mEnableOutput = false;
185 return;
186 }
187
188 mOutTree = std::make_unique<TTree>("calib", "Gain calibration");
189
190 LOG(info) << "Created output file trd_calibgain.root";
191}
192
194{
195 if (!mEnableOutput) {
196 return;
197 }
198
199 try {
200 mOutFile->cd();
201 mOutTree->Write();
202 mOutTree.reset();
203 mOutFile->Close();
204 mOutFile.reset();
205 } catch (std::exception const& e) {
206 LOG(error) << "Failed to write calibration data file, reason: " << e.what();
207 }
208 mEnableOutput = false;
209}
210
211} // namespace o2::trd
TimeSlot-based calibration of gain.
Global TRD definitions and constants.
Definition of the Names Generator class.
std::ostringstream debug
long getStartTimeMS() const
Definition TimeSlot.h:50
const Container * getContainer() const
Definition TimeSlot.h:53
static std::string generateFileName(const std::string &inp)
Definition CcdbApi.cxx:798
static constexpr long DAY
decltype(auto) get(R binding, int part=0) const
InputRecord & inputs()
The inputs associated with this processing context.
void setMPVdEdx(int iDet, float mpv)
Definition CalGain.h:34
void finalizeSlot(Slot &slot) final
Slot & emplaceNewSlot(bool front, TFType tStart, TFType tEnd) final
void retrievePrev(o2::framework::ProcessingContext &pc)
constexpr float MPVDEDXDEFAULT
default Most Probable Value of TRD dEdx
Definition Constants.h:80
constexpr int MAXCHAMBER
the maximum number of installed chambers
Definition Constants.h:30
constexpr int NBINSGAINCALIB
number of bins in the charge (Q0+Q1+Q2) histogram for gain calibration
Definition Constants.h:79
static std::string getClassName(const T &obj)
get the class name of the object
LOG(info)<< "Compressed in "<< sw.CpuTime()<< " s"
uint64_t const void const *restrict const msg
Definition x9.h:153