Project
Loading...
Searching...
No Matches
MisalignmentHits.cxx
Go to the documentation of this file.
1// Copyright 2020-2022 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
14#include "ITS3Base/ITS3Params.h"
17#include "Framework/Logger.h"
18
19#include "Math/Factory.h"
20#include "Math/UnaryOperators.h"
21#include "TGeoNode.h"
22#include "TGeoBBox.h"
23#include "TString.h"
24
25#include <memory>
26#include <string>
27#include <cstring>
28#include <algorithm>
29
30namespace o2::its3::align
31{
32
34{
35 if (o2::its3::ITS3Params::Instance().misalignmentHitsUseProp) {
36 mMethod = PropMethod::Propagator;
37 } else {
38 mMethod = PropMethod::Line;
39 }
40
42
43 mMinimizer.reset(ROOT::Math::Factory::CreateMinimizer("Minuit2", "Migrad"));
44 if (mMinimizer == nullptr) {
45 LOGP(fatal, "Cannot create minimizer");
46 }
47 mMinimizer->SetMaxFunctionCalls(1'000'000'000);
48 mMinimizer->SetStrategy(0);
49 mMinimizer->SetPrintLevel(0);
50
51 if (mMethod == PropMethod::Propagator) {
52 LOGP(info, "Using propagator to find intersection");
54 mMCReader = std::make_unique<o2::steer::MCKinematicsReader>(prefix, o2::steer::MCKinematicsReader::Mode::kMCKine);
55 mMinimizer->SetFunction(mPropagator);
56 } else {
57 LOGP(info, "Using local straight-line to find intersection");
58 mMinimizer->SetFunction(mLine);
59 }
60
61 resetStats();
62
63 if (auto file = o2::its3::ITS3Params::Instance().misalignmentHitsParams; file.empty()) {
64 LOGP(fatal, "No parameter file specified");
65 } else {
66 mDeformations.init(file);
67 }
68}
69
70std::optional<o2::itsmft::Hit> MisAlignmentHits::processHit(int iEvent, const o2::itsmft::Hit& hit)
71{
72 ++mStats[Stats::kHitTotal];
73
75 ++mStats[Stats::kHitIsOB];
76 return hit;
77 }
78 ++mStats[Stats::kHitIsIB];
79
80 // Set the working hits
81 mCurHit = hit;
82 mCurWorkingHits[WorkingHit::kEntering] = WorkingHit(iEvent, WorkingHit::kEntering, hit);
83 mCurWorkingHits[WorkingHit::kExiting] = WorkingHit(iEvent, WorkingHit::kExiting, hit);
84
85 // Do work
86 if (!deformHit(WorkingHit::kEntering) || !deformHit(WorkingHit::kExiting)) {
87 ++mStats[Stats::kHitDead];
88 return std::nullopt;
89 }
90 ++mStats[Stats::kHitAlive];
91
92 // Set the possibly new detectorIDs with mid point approximation
93 auto midPointOrig = mCurWorkingHits[WorkingHit::kEntering].mPoint + (mCurWorkingHits[WorkingHit::kExiting].mPoint - mCurWorkingHits[WorkingHit::kEntering].mPoint) * 0.5;
94 auto midPointDef = mCurWorkingHits[WorkingHit::kEntering].mPointDef + (mCurWorkingHits[WorkingHit::kExiting].mPointDef - mCurWorkingHits[WorkingHit::kEntering].mPointDef) * 0.5;
95 const short idDef = getDetID(midPointDef), idOrig = getDetID(midPointOrig);
96 if (idDef == -1) {
97 return std::nullopt;
98 }
99
100 if (idDef != idOrig) {
101 ++mStats[Stats::kHitMigrated];
102 } else {
103 ++mStats[Stats::kHitNotMigrated];
104 }
105
106 if constexpr (false) {
109 bool crossesBoundary{false};
110 TGeoNode *nEnt{nullptr}, *nExt{nullptr};
111 {
112 auto dirEnt = mCurWorkingHits[WorkingHit::kEntering].mPointDef - midPointDef;
113 auto stepEnt = std::min(static_cast<double>(dirEnt.R()), std::abs(dirEnt.R() - 5.e-4));
114 auto dirEntU = dirEnt.Unit();
115 gGeoManager->SetCurrentPoint(midPointDef.X(), midPointDef.Y(), midPointDef.Z());
116 gGeoManager->SetCurrentDirection(dirEntU.X(), dirEntU.Y(), dirEntU.Z());
117 nEnt = gGeoManager->FindNextBoundaryAndStep(stepEnt, false);
118 if (gGeoManager->IsOnBoundary()) {
119 ++mStats[Stats::kHitEntBoundary];
120 crossesBoundary = true;
121 }
122 }
123 {
124 auto dirExt = midPointDef - mCurWorkingHits[WorkingHit::kEntering].mPointDef;
125 auto stepExt = std::min(static_cast<double>(dirExt.R()), std::abs(dirExt.R() - 5.e-4));
126 auto dirExtU = dirExt.Unit();
127 gGeoManager->SetCurrentPoint(midPointDef.X(), midPointDef.Y(), midPointDef.Z());
128 gGeoManager->SetCurrentDirection(dirExtU.X(), dirExtU.Y(), dirExtU.Z());
129 nExt = gGeoManager->FindNextBoundaryAndStep(stepExt, false);
130 if (gGeoManager->IsOnBoundary()) {
131 ++mStats[Stats::kHitExtBoundary];
132 crossesBoundary = true;
133 }
134 }
135
136 if (crossesBoundary && nEnt != nullptr && nExt != nullptr) {
137 if (nEnt != nExt) {
138 return std::nullopt;
139 } else {
140 ++mStats[Stats::kHitSameBoundary]; // indicates that the step size is too large and we end up in the mother volume; just pretend that his fine for now
141 }
142 }
143 ++mStats[Stats::kHitNoBoundary];
144 }
145
146 // Get new postion
147 mCurHit.SetPosStart(mCurWorkingHits[WorkingHit::kEntering].mPointDef);
148 mCurHit.SetPos(mCurWorkingHits[WorkingHit::kExiting].mPointDef);
149 mCurHit.SetDetectorID(idDef);
150
151 ++mStats[Stats::kHitSuccess];
152 return mCurHit;
153}
154
155bool MisAlignmentHits::deformHit(WorkingHit::HitType t)
156{
157 auto& wHit = mCurWorkingHits[t];
158
159 mMinimizer->Clear(); // clear for next iteration
160 constexpr double minStep{1e-5};
161 constexpr double zMargin{4.0};
162 constexpr double phiMargin{0.4};
163 if (mMethod == PropMethod::Line) {
164 prepareLineMethod(t);
165 mMinimizer->SetVariable(0, "t", 0.0, minStep); // this is left as a free parameter on since t is very small since start and end of hit are close
166 } else {
167 if (!preparePropagatorMethod(t)) {
168 return false;
169 }
170 mMinimizer->SetVariable(0, "r", mPropagator.mTrack.getX(), minStep); // this is left as a free parameter on since t is very small since start and end of hit are close
171 }
172 mMinimizer->SetLimitedVariable(1, "phiStar", wHit.mPhi, minStep,
173 std::max(static_cast<double>(wHit.mPhiBorder1), static_cast<double>(wHit.mPhi) - phiMargin),
174 std::min(static_cast<double>(wHit.mPhiBorder2), static_cast<double>(wHit.mPhi) + phiMargin));
175 mMinimizer->SetLimitedVariable(2, "zStar", wHit.mPoint.Z(), minStep,
176 std::max(static_cast<double>(-constants::segment::lengthSensitive / 2.f), static_cast<double>(wHit.mPoint.Z()) - zMargin),
177 std::min(static_cast<double>(constants::segment::lengthSensitive / 2.f), static_cast<double>(wHit.mPoint.Z()) + zMargin));
178
179 mMinimizer->Minimize(); // perform the actual minimization
180
181 auto ss = mMinimizer->Status();
182 if (ss == 1) {
183 ++mStats[Stats::kMinimizerCovPos];
184 } else if (ss == 2) {
185 ++mStats[Stats::kMinimizerHesse];
186 } else if (ss == 3) {
187 ++mStats[Stats::kMinimizerEDM];
188 } else if (ss == 4) {
189 ++mStats[Stats::kMinimizerLimit];
190 } else if (ss == 5) {
191 ++mStats[Stats::kMinimizerOther];
192 } else {
193 ++mStats[Stats::kMinimizerConverged];
194 }
195
196 if (ss == 0 || ss == 1) { // for Minuit2 0=ok, 1=ok with pos. forced hesse
197 ++mStats[Stats::kMinimizerStatusOk];
198 if (mMinimizer->MinValue() < 2e-4) { // within 2 um considering the pixel pitch this good enough
199 ++mStats[Stats::kMinimizerValueOk];
200 } else {
201 ++mStats[Stats::kMinimizerValueBad];
202 return false;
203 }
204 } else {
205 ++mStats[Stats::kMinimizerStatusBad];
206 return false;
207 }
208
209 // Valid solution found; calculate new position on ideal geo
210 wHit.recalculateIdeal(static_cast<float>(mMinimizer->X()[1]), static_cast<float>(mMinimizer->X()[2]));
211
212 return true;
213}
214
215short MisAlignmentHits::getDetID(const o2::math_utils::Point3D<float>& point)
216{
217 // Do not modify the path, I do not know if this is needed but lets be safe
218 gGeoManager->PushPath();
219 auto id = getDetIDFromCords(point);
220 gGeoManager->PopPath();
221 return id;
222}
223
224short MisAlignmentHits::getDetIDFromCords(const o2::math_utils::Point3D<float>& point)
225{
226 // retrive if any the node which constains the point
227 const auto node = gGeoManager->FindNode(point.X(), point.Y(), point.Z());
228 if (node == nullptr) {
229 ++mStats[Stats::kFindNodeFailed];
230 return -1;
231 }
232 ++mStats[Stats::kFindNodeSuccess];
233
234 // check if this node is a sensitive volume
235 const std::string path = gGeoManager->GetPath();
236 if (path.find(o2::its::GeometryTGeo::getITS3SensorPattern()) == std::string::npos) {
237 ++mStats[Stats::kProjNonSensitive];
238 return -1;
239 }
240 ++mStats[Stats::kProjSensitive];
241
242 return getDetIDFromPath(path);
243}
244
245short MisAlignmentHits::getDetIDFromPath(const std::string& path) const
246{
247 static const std::regex pattern{R"(/cave_1/barrel_1/ITSV_2/ITSUWrapVol0_1/ITS3Layer(\d+)_(\d+)/ITS3CarbonForm(\d+)_(\d+)/ITS3Chip(\d+)_(\d+)/ITS3Segment(\d+)_(\d+)/ITS3RSU(\d+)_(\d+)/ITS3Tile(\d+)_(\d+)/ITS3PixelArray(\d+)_(\d+))"};
248 if (std::smatch matches; std::regex_search(path, matches, pattern)) {
249 if (matches.size() == 15) {
250 int iLayer = std::stoi(matches[1]);
251 int iCarbonForm = std::stoi(matches[4]);
252 int iSegment = std::stoi(matches[8]);
253 int iRSU = std::stoi(matches[10]);
254 int iTile = std::stoi(matches[12]);
255 return mGeo->getChipIndex(iLayer, iCarbonForm, 0, iSegment, iRSU, iTile);
256 } else {
257 LOGP(fatal, "Path did not contain expected number of matches ({})!", matches.size());
258 }
259 } else {
260 LOGP(fatal, "Path was not matched ({})!", path);
261 }
262 __builtin_unreachable();
263}
264
266{
267 auto makeFraction = [&](Stats n, Stats d) -> float { return static_cast<float>(mStats[n]) / static_cast<float>(mStats[d] + mStats[n]); };
268 LOGP(info, "Processed {} Hits (IB:{}; OB:{}) ({:.2f}%):", mStats[Stats::kHitTotal], mStats[Stats::kHitIsIB], mStats[Stats::kHitIsOB], makeFraction(Stats::kHitIsIB, Stats::kHitIsOB));
269 LOGP(info, " - Minimizer Status: {} ok {} bad ({:.2f}%)", mStats[Stats::kMinimizerStatusOk], mStats[Stats::kMinimizerStatusBad], makeFraction(Stats::kMinimizerStatusOk, Stats::kMinimizerStatusBad));
270 LOGP(info, " - Minimizer Value: {} ok {} bad ({:.2f}%)", mStats[Stats::kMinimizerValueOk], mStats[Stats::kMinimizerValueBad], makeFraction(Stats::kMinimizerValueOk, Stats::kMinimizerValueBad));
271 LOGP(info, " - Minimizer Detailed: {} Converged {} pos. forced Hesse ({:.2f}%)", mStats[Stats::kMinimizerConverged], mStats[Stats::kMinimizerHesse], makeFraction(Stats::kMinimizerConverged, Stats::kMinimizerHesse));
272 LOGP(info, " - Minimizer Detailed: {} EDM {} call limit {} other ({:.2f}%)", mStats[Stats::kMinimizerEDM], mStats[Stats::kMinimizerLimit], mStats[Stats::kMinimizerOther], makeFraction(Stats::kMinimizerEDM, Stats::kMinimizerLimit));
273 LOGP(info, " - FindNode: {} ok {} failed", mStats[Stats::kFindNodeSuccess], mStats[Stats::kFindNodeFailed]);
274 LOGP(info, " - IsSensitve: {} yes {} no ({:.2f}%)", mStats[Stats::kProjSensitive], mStats[Stats::kProjNonSensitive], makeFraction(Stats::kProjSensitive, Stats::kProjNonSensitive));
275 LOGP(info, " - IsAlive: {} yes {} no ({:.2f}%)", mStats[Stats::kHitAlive], mStats[Stats::kHitDead], makeFraction(Stats::kHitAlive, Stats::kHitDead));
276 LOGP(info, " - HasMigrated: {} yes {} no ({:.2f}%)", mStats[Stats::kHitMigrated], mStats[Stats::kHitNotMigrated], makeFraction(Stats::kHitMigrated, Stats::kHitNotMigrated));
277 // LOGP(info, " - Crosses Boundary: {} entering {} exiting {} same {} no", mStats[Stats::kHitEntBoundary], mStats[Stats::kHitExtBoundary], mStats[Stats::kHitSameBoundary], mStats[Stats::kHitNoBoundary]);
278 if (mMethod == PropMethod::Propagator) {
279 LOGP(info, " - Propagator: {} null track {} null pdg", mStats[Stats::kPropTrackNull], mStats[Stats::kPropPDGNull]);
280 }
281 LOGP(info, " --> Good Hits {} ({:.2f}%)", mStats[Stats::kHitSuccess], makeFraction(Stats::kHitSuccess, Stats::kHitIsIB));
282}
283
284void MisAlignmentHits::prepareLineMethod(WorkingHit::HitType from)
285{
286 // Set the starint point and radius
287 // always start from the entering hit that way t is always pos. defined
288 mLine.mStart = mCurWorkingHits[WorkingHit::kEntering].mPoint;
289 mLine.mRadius = mCurWorkingHits[from].mRadius;
290 mLine.mSensorID = mCurWorkingHits[from].mSensorID;
291 mLine.mPhiTot = mCurWorkingHits[from].mPhiBorder2 - mCurWorkingHits[from].mPhiBorder1;
292 mLine.mPhi1 = mCurWorkingHits[from].mPhiBorder1;
293 // Calculate the direction vector
294 mLine.mD[0] = mCurWorkingHits[WorkingHit::kExiting].mPoint.X() - mCurWorkingHits[WorkingHit::kEntering].mPoint.X();
295 mLine.mD[1] = mCurWorkingHits[WorkingHit::kExiting].mPoint.Y() - mCurWorkingHits[WorkingHit::kEntering].mPoint.Y();
296 mLine.mD[2] = mCurWorkingHits[WorkingHit::kExiting].mPoint.Z() - mCurWorkingHits[WorkingHit::kEntering].mPoint.Z();
297}
298
299double MisAlignmentHits::StraightLine::DoEval(const double* x) const
300{
301 const double t = x[0];
302 const double phi = x[1];
303 const double z = x[2];
304 const double nphi = std::clamp((phi - mPhi1) * 2.0 / mPhiTot - 1.0, -1.0, 1.0);
305 const double nz = std::clamp((z - (-constants::segment::lengthSensitive / 2.0)) * 2.0 / constants::segment::lengthSensitive - 1.0, -1.0, 1.0);
306
308 double xline = mStart.X() + t * mD[0],
309 yline = mStart.Y() + t * mD[1],
310 zline = mStart.Z() + t * mD[2];
311
312 // Find the point of the deformed geometry given a certain phi' and z'
313 double xideal = mRadius * std::cos(phi), yideal = mRadius * std::sin(phi),
314 zideal = z;
315 const auto [dx, dy, dz] = mMis->getDeformation(mSensorID, nphi, nz);
316 double xdef = xideal + dx, ydef = yideal + dy, zdef = zideal + dz;
317
318 // Minimize the euclidean distance of the line point and the deformed point
319 return std::hypot(xline - xdef, yline - ydef, zline - zdef);
320}
321
322bool MisAlignmentHits::preparePropagatorMethod(WorkingHit::HitType from)
323{
324 mPropagator.mRadius = mCurWorkingHits[from].mRadius;
325 mPropagator.mSensorID = mCurWorkingHits[from].mSensorID;
326 mPropagator.mPhiTot = mCurWorkingHits[from].mPhiBorder2 - mCurWorkingHits[from].mPhiBorder1;
327 mPropagator.mPhi1 = mCurWorkingHits[from].mPhiBorder1;
328 const auto mcTrack = mMCReader->getTrack(mCurWorkingHits[from].mEvent, mCurWorkingHits[from].mTrackID);
329 if (mcTrack == nullptr) {
330 ++mStats[Stats::kPropTrackNull];
331 return false;
332 }
333 const std::array<float, 3> xyz{(float)mcTrack->GetStartVertexCoordinatesX(), (float)mcTrack->GetStartVertexCoordinatesY(), (float)mcTrack->GetStartVertexCoordinatesZ()},
334 pxyz{(float)mcTrack->GetStartVertexMomentumX(), (float)mcTrack->GetStartVertexMomentumY(), (float)mcTrack->GetStartVertexMomentumZ()};
335 const TParticlePDG* pPDG = TDatabasePDG::Instance()->GetParticle(mcTrack->GetPdgCode());
336 if (pPDG == nullptr) {
337 ++mStats[Stats::kPropPDGNull];
338 return false;
339 }
340 mPropagator.mTrack = o2::track::TrackPar(xyz, pxyz, TMath::Nint(pPDG->Charge() / 3), false);
341 mPropagator.mBz = o2::base::Propagator::Instance()->getNominalBz();
342 return true;
343}
344
345double MisAlignmentHits::Propagator::DoEval(const double* x) const
346{
347 const double r = x[0];
348 const double phi = x[1];
349 const double z = x[2];
350 const double nphi = (phi - mPhi1) * 2.0 / mPhiTot - 1.0;
351 const double nz = (z - (-constants::segment::lengthSensitive / 2.0)) * 2.0 / constants::segment::lengthSensitive - 1.0;
352
353 auto trc = mTrack;
354 if (!trc.propagateTo(r, mBz)) {
355 return 999;
356 }
357 const auto glo = trc.getXYZGlo();
358
359 // Find the point of the deformed geometry given a certain phi' and z'
360 double xideal = mRadius * std::cos(phi), yideal = mRadius * std::sin(phi),
361 zideal = z;
362 const auto [dx, dy, dz] = mMis->getDeformation(mSensorID, nphi, nz);
363 double xdef = xideal + dx, ydef = yideal + dy, zdef = zideal + dz;
364
365 // Minimize the euclidean distance of the propagator point and the deformed point
366 return std::hypot(glo.X() - xdef, glo.Y() - ydef, glo.Z() - zdef);
367}
368
369} // namespace o2::its3::align
Definition of the SegmentationSuperAlpide class.
void SetDetectorID(short detID)
Definition BaseHits.h:78
short GetDetectorID() const
Definition BaseHits.h:73
void SetPos(math_utils::Point3D< T > const &p)
Definition BaseHits.h:88
GPUd() value_type estimateLTFast(o2 static GPUd() float estimateLTIncrement(const o2 PropagatorImpl * Instance(bool uninitialized=false)
Definition Propagator.h:143
void init(const std::filesystem::path &)
std::optional< o2::itsmft::Hit > processHit(int iEvent, const o2::itsmft::Hit &hit)
static const char * getITS3SensorPattern()
static GeometryTGeo * Instance()
int getChipIndex(int lay, int detInLay) const
void SetPosStart(const math_utils::Point3D< Float_t > &p)
Definition Hit.h:98
GLdouble n
Definition glcorearb.h:1982
GLint GLenum GLint x
Definition glcorearb.h:403
GLsizei const GLchar *const * path
Definition glcorearb.h:3591
GLboolean r
Definition glcorearb.h:1233
GLuint id
Definition glcorearb.h:650
GLdouble GLdouble GLdouble z
Definition glcorearb.h:843
bool isDetITS3(T detID)
Definition SpecsV2.h:161
TrackParF TrackPar
Definition Track.h:29
std::string digitizationgeometry_prefix
Definition DigiParams.h:30
std::array< uint16_t, 5 > pattern