Project
Loading...
Searching...
No Matches
CalibratorVdExB.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"
28#include <TFile.h>
29#include <TTree.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
41double FitFunctor::calculateDeltaAlphaSim(double vdFit, double laFit, double impactAng) const
42{
43 auto xDir = TMath::Cos(impactAng);
44 auto yDir = TMath::Sin(impactAng);
45 double slope = (TMath::Abs(xDir) < 1e-7) ? 1e7 : yDir / xDir;
46 auto laTan = TMath::Tan(laFit);
47 double lorentzSlope = (TMath::Abs(laTan) < 1e-7) ? 1e7 : 1. / laTan;
48
49 // hit point of incoming track with anode plane
50 double xAnodeHit = mAnodePlane / slope;
51 double yAnodeHit = mAnodePlane;
52
53 // hit point at anode plane of Lorentz angle shifted cluster from the entrance -> independent of true drift velocity
54 double xLorentzAnodeHit = mAnodePlane / lorentzSlope;
55 double yLorentzAnodeHit = mAnodePlane;
56
57 // cluster location within drift cell of cluster from entrance after drift velocity ratio is applied
58 double xLorentzDriftHit = xLorentzAnodeHit;
59 double yLorentzDriftHit = mAnodePlane - mAnodePlane * (vdPreCorr[currDet] / vdFit);
60
61 // reconstructed hit of first cluster at chamber entrance after pre Lorentz angle correction
62 double xLorentzDriftHitPreCorr = xLorentzAnodeHit - (mAnodePlane - yLorentzDriftHit) * TMath::Tan(laPreCorr[currDet]);
63 double yLorentzDriftHitPreCorr = mAnodePlane - mAnodePlane * (vdPreCorr[currDet] / vdFit);
64
65 double impactAngleSim = TMath::ATan2(yAnodeHit, xAnodeHit);
66
67 double deltaXLorentzDriftHit = xAnodeHit - xLorentzDriftHitPreCorr;
68 double deltaYLorentzDriftHit = yAnodeHit - yLorentzDriftHitPreCorr;
69 double impactAngleRec = TMath::ATan2(deltaYLorentzDriftHit, deltaXLorentzDriftHit);
70
71 double deltaAngle = (impactAngleRec - impactAngleSim); // * TMath::RadToDeg());
72
73 return deltaAngle;
74}
75
76double FitFunctor::operator()(const double* par) const
77{
78 double sum = 0; // this value is minimized
79 for (int iBin = 1; iBin <= profiles[currDet]->GetNbinsX(); ++iBin) {
80 auto impactAngle = (profiles[currDet]->GetBinCenter(iBin) + 90) * TMath::DegToRad();
81 auto deltaAlpha = profiles[currDet]->GetBinContent(iBin) * TMath::DegToRad();
82 if (TMath::Abs(deltaAlpha) < 1e-7) {
83 continue;
84 }
85 if (impactAngle < lowerBoundAngleFit || impactAngle > upperBoundAngleFit) {
86 continue;
87 }
89 sum += TMath::Power(deltaAlphaSim - deltaAlpha, 2);
90 }
91 return sum;
92}
93
95
97{
98 // reset the CCDB output vectors
99 mInfoVector.clear();
100 mObjectVector.clear();
101}
102
104{
105 if (mInitDone) {
106 return;
107 }
108
109 // fit is done in region where ion tails are small, close to lorentz angle
110 // we want an approximate value of the lorentz angle in order to define better fit boundaries, coinciding with profile bin edges
111 float bz = o2::base::Propagator::Instance()->getNominalBz();
112 // default angle with zero field is slightly shifted
113 float lorentzAngleAvg = -1.f;
114 if (TMath::Abs(bz - 2) < 0.1f) {
115 lorentzAngleAvg = 3.f;
116 }
117 if (TMath::Abs(bz + 2) < 0.1f) {
118 lorentzAngleAvg = -5.f;
119 }
120 if (TMath::Abs(bz - 5) < 0.1f) {
121 lorentzAngleAvg = 7.f;
122 }
123 if (TMath::Abs(bz + 5) < 0.1f) {
124 lorentzAngleAvg = -9.f;
125 }
126
127 LOGP(info, "b field: {} lorentz angle start: {}", bz, lorentzAngleAvg);
128
129 mFitFunctor.lowerBoundAngleFit = (80 + lorentzAngleAvg) * TMath::DegToRad();
130 mFitFunctor.upperBoundAngleFit = (100 + lorentzAngleAvg) * TMath::DegToRad();
131 mFitFunctor.mAnodePlane = GeometryBase::camHght() / (2.f * 100.f);
132 for (int iDet = 0; iDet < MAXCHAMBER; ++iDet) {
133 mFitFunctor.profiles[iDet] = std::make_unique<TProfile>(Form("profAngleDiff_%i", iDet), Form("profAngleDiff_%i", iDet), NBINSANGLEDIFF, -MAXIMPACTANGLE, MAXIMPACTANGLE);
134 }
135
136 mFitter.SetFCN<FitFunctor>(2, mFitFunctor, mParamsStart);
137 mFitter.Config().ParSettings(ParamIndex::LA).SetLimits(constants::EXBMIN, constants::EXBMAX);
138 mFitter.Config().ParSettings(ParamIndex::LA).SetStepSize(.01);
139 mFitter.Config().ParSettings(ParamIndex::VD).SetLimits(constants::VDRIFTMIN, constants::VDRIFTMAX);
140 mFitter.Config().ParSettings(ParamIndex::VD).SetStepSize(.01);
141 ROOT::Math::MinimizerOptions opt;
142 opt.SetMinimizerType("Minuit2");
143 opt.SetMinimizerAlgorithm("Migrad");
144 opt.SetPrintLevel(0);
145 opt.SetMaxFunctionCalls(1'000);
146 opt.SetTolerance(.001);
147 mFitter.Config().SetMinimizerOptions(opt);
148
149 // set tree addresses
150 if (mEnableOutput) {
151 mOutTree->Branch("lorentzAngle", &mFitFunctor.laPreCorr);
152 mOutTree->Branch("vDrift", &mFitFunctor.vdPreCorr);
153 for (int iDet = 0; iDet < MAXCHAMBER; ++iDet) {
154 mOutTree->Branch(fmt::format("residuals_{:d}", iDet).c_str(), mFitFunctor.profiles[iDet].get());
155 }
156 }
157
158 mInitDone = true;
159}
160
162{
163 // We either get a pointer to a valid object from the last ~hour or to the default object
164 // which is always present. The first has precedence over the latter.
165 auto dataCalVdriftExB = pc.inputs().get<o2::trd::CalVdriftExB*>("calvdexb");
166 std::string msg = "Default Object";
167 // We check if the object we got is the default one by comparing it to the defaults.
168 for (int iDet = 0; iDet < MAXCHAMBER; ++iDet) {
169 if (dataCalVdriftExB->getVdrift(iDet) != VDRIFTDEFAULT ||
170 dataCalVdriftExB->getExB(iDet) != EXBDEFAULT) {
171 msg = "Previous Object";
172 break;
173 }
174 }
175 LOG(info) << "Calibrator: From CCDB retrieved " << msg;
176
177 // Here we set each entry regardless if it is the default or not.
178 for (int iDet = 0; iDet < MAXCHAMBER; ++iDet) {
179 mFitFunctor.laPreCorr[iDet] = dataCalVdriftExB->getExB(iDet);
180 mFitFunctor.vdPreCorr[iDet] = dataCalVdriftExB->getVdrift(iDet);
181 }
182}
183
185{
186 // do actual calibration for the data provided in the given slot
187 TStopwatch timer;
188 timer.Start();
190 auto laFitResults = mFitFunctor.laPreCorr;
191 auto vdFitResults = mFitFunctor.vdPreCorr;
192 auto residHists = slot.getContainer();
193 for (int iDet = 0; iDet < MAXCHAMBER; ++iDet) {
194 int sumEntries = 0;
195 mFitFunctor.profiles[iDet]->Reset();
196 mFitFunctor.currDet = iDet;
197 for (int iBin = 0; iBin < NBINSANGLEDIFF; ++iBin) {
198 // fill profiles
199 auto angleDiffSum = residHists->getHistogramEntry(iDet * NBINSANGLEDIFF + iBin);
200 auto nEntries = residHists->getBinCount(iDet * NBINSANGLEDIFF + iBin);
201 sumEntries += nEntries;
202 if (nEntries > 0) {
203 mFitFunctor.profiles[iDet]->Fill(2 * iBin - MAXIMPACTANGLE, angleDiffSum / nEntries, nEntries);
204 }
205 }
206 // Check if we have the minimum amount of entries
207 if (sumEntries < mMinEntriesChamber) {
208 LOGF(debug, "Chamber %d did not reach minimum amount of entries for refit: %d", iDet, sumEntries);
209 continue;
210 }
211 float laPreCorrTemp = mFitFunctor.laPreCorr[iDet];
212 float vdPreCorrTemp = mFitFunctor.vdPreCorr[iDet];
213 // Here we start from uncalibrated values, otherwise online calibration does not work properly
214 mFitFunctor.laPreCorr[iDet] = EXBDEFAULT;
215 mFitFunctor.vdPreCorr[iDet] = VDRIFTDEFAULT;
216
217 // Reset Start Parameter
218 mParamsStart[ParamIndex::LA] = 0.0;
219 mParamsStart[ParamIndex::VD] = 1.0;
220 mFitter.FitFCN();
221 auto fitResult = mFitter.Result();
222 if (fitResult.MinFcnValue() > 0.03) {
223 LOGF(debug, "Chamber %d fit did not converge properly, minimization value too high: %f", iDet, fitResult.MinFcnValue());
224 // The fit did not work properly, so we keep previous values
225 mFitFunctor.laPreCorr[iDet] = laPreCorrTemp;
226 mFitFunctor.vdPreCorr[iDet] = vdPreCorrTemp;
227 continue;
228 }
229 laFitResults[iDet] = fitResult.Parameter(ParamIndex::LA);
230 vdFitResults[iDet] = fitResult.Parameter(ParamIndex::VD);
231 LOGF(debug, "Fit result for chamber %i: vd=%f, la=%f, minimizer value=%f", iDet, vdFitResults[iDet], laFitResults[iDet] * TMath::RadToDeg(), fitResult.MinFcnValue());
232 // Update fit values for next fit
233 mFitFunctor.laPreCorr[iDet] = laFitResults[iDet];
234 mFitFunctor.vdPreCorr[iDet] = vdFitResults[iDet];
235 }
236 timer.Stop();
237 LOGF(info, "Done fitting angular residual histograms. CPU time: %f, real time: %f", timer.CpuTime(), timer.RealTime());
238
239 // Fill Tree and log to debug
240 if (mEnableOutput) {
241 mOutTree->Fill();
242 for (int iDet = 0; iDet < MAXCHAMBER; ++iDet) {
243 LOGF(debug, "Fit result for chamber %i: vd=%f, la=%f", iDet, vdFitResults[iDet], laFitResults[iDet] * TMath::RadToDeg());
244 }
245 }
246
247 // assemble CCDB object
248 CalVdriftExB calObject;
249 for (int iDet = 0; iDet < MAXCHAMBER; ++iDet) {
250 // For Chambers which did not have the minimum amount of entries in this slot e.g. missing, empty chambers.
251 // We will reuse the prevoius one. This would have been read either from the ccdb or come from a previous successful fit.
252 calObject.setVdrift(iDet, vdFitResults[iDet]);
253 calObject.setExB(iDet, laFitResults[iDet]);
254 }
255 auto clName = o2::utils::MemFileHelper::getClassName(calObject);
256 auto flName = o2::ccdb::CcdbApi::generateFileName(clName);
257 std::map<std::string, std::string> metadata; // TODO: do we want to store any meta data?
258 long startValidity = slot.getStartTimeMS() - 10 * o2::ccdb::CcdbObjectInfo::SECOND;
259 mInfoVector.emplace_back("TRD/Calib/CalVdriftExB", clName, flName, metadata, startValidity, startValidity + 1 * o2::ccdb::CcdbObjectInfo::HOUR);
260 mObjectVector.push_back(calObject);
261}
262
264{
265 auto& container = getSlots();
266 auto& slot = front ? container.emplace_front(tStart, tEnd) : container.emplace_back(tStart, tEnd);
267 slot.setContainer(std::make_unique<AngularResidHistos>());
268 return slot;
269}
270
272{
273 mEnableOutput = true;
274 mOutFile = std::make_unique<TFile>("trd_calibVdriftExB.root", "RECREATE");
275 if (mOutFile->IsZombie()) {
276 LOG(error) << "Failed to create output file!";
277 mEnableOutput = false;
278 return;
279 }
280 mOutTree = std::make_unique<TTree>("calib", "VDrift&ExB calibration");
281 LOG(info) << "Created output file trd_calibVdriftExB.root";
282}
283
285{
286 if (!mEnableOutput) {
287 return;
288 }
289
290 try {
291 mOutFile->cd();
292 mOutTree->Write();
293 mOutTree.reset();
294 mOutFile->Close();
295 mOutFile.reset();
296 } catch (std::exception const& e) {
297 LOG(error) << "Failed to write calibration data file, reason: " << e.what();
298 }
299 // after closing, we won't open a new file
300 mEnableOutput = false;
301}
302
303} // namespace o2::trd
TimeSlot-based calibration of vDrift and ExB.
Global TRD definitions and constants.
std::ostringstream debug
Definition of the Names Generator class.
uint16_t slope
Definition RawData.h:1
GPUd() value_type estimateLTFast(o2 static GPUd() float estimateLTIncrement(const o2 PropagatorImpl * Instance(bool uninitialized=false)
Definition Propagator.h:178
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:824
static constexpr long HOUR
static constexpr long SECOND
decltype(auto) get(R binding, int part=0) const
InputRecord & inputs()
The inputs associated with this processing context.
void setExB(int iDet, float exb)
void setVdrift(int iDet, float vd)
void retrievePrev(o2::framework::ProcessingContext &pc)
void finalizeSlot(Slot &slot) final
Slot & emplaceNewSlot(bool front, TFType tStart, TFType tEnd) final
float sum(float s, o2::dcs::DataPointValue v)
Definition dcs-ccdb.cxx:39
constexpr double VDRIFTMIN
min value for vDrift
Definition Constants.h:78
constexpr double EXBDEFAULT
default value for LorentzAngle
Definition Constants.h:80
constexpr double EXBMAX
max value for LorentzAngle
Definition Constants.h:82
constexpr int MAXCHAMBER
the maximum number of installed chambers
Definition Constants.h:30
constexpr double EXBMIN
min value for LorentzAngle
Definition Constants.h:81
constexpr double VDRIFTMAX
max value for vDrift
Definition Constants.h:79
constexpr double VDRIFTDEFAULT
default value for vDrift
Definition Constants.h:77
constexpr int NBINSANGLEDIFF
the number of bins for the track angle used for the vDrift and ExB calibration based on the tracking
Definition Constants.h:76
constexpr float MAXIMPACTANGLE
the maximum impact angle for tracks relative to the TRD detector plane to be considered for vDrift an...
Definition Constants.h:75
std::array< double, constants::MAXCHAMBER > laPreCorr
LorentzAngle from previous Run.
double mAnodePlane
distance of the TRD anode plane from the drift cathodes in m
int currDet
the current TRD chamber number
double operator()(const double *par) const
std::array< std::unique_ptr< TProfile >, constants::MAXCHAMBER > profiles
profile histograms for each TRD chamber
double calculateDeltaAlphaSim(double vdFit, double laFit, double impactAng) const
std::array< double, constants::MAXCHAMBER > vdPreCorr
vDrift from previous Run
float upperBoundAngleFit
upper bound for fit angle
float lowerBoundAngleFit
lower bound for fit angle
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