Project
Loading...
Searching...
No Matches
AlignableDetectorTRD.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
20#include "Align/Controller.h"
21#include "Align/AlignConfig.h"
22#include "TRDBase/Geometry.h"
30#include <TMath.h>
31#include <TGeoManager.h>
32
33namespace o2
34{
35namespace align
36{
37using namespace TMath;
38const char* AlignableDetectorTRD::CalibDOFName[AlignableDetectorTRD::NCalibParams] = {"DZdTglNRC", "DVDriftT"};
39
40//____________________________________________
45
46//____________________________________________
48{
49 // define TRD volumes
50 auto geo = o2::trd::Geometry::instance();
51 geo->createPadPlaneArray();
52 geo->createClusterMatrixArray(); // ideal T2L matrices
53
54 AlignableSensorTRD* chamb = nullptr;
56 AlignableVolume* volTRD = nullptr; // fictious envelope
57
58 addVolume(volTRD = new AlignableVolume("TRD_envelope", getDetLabel(), mController));
59 volTRD->setDummyEnvelope();
60
61 for (int ilr = 0; ilr < o2::trd::constants::NLAYER; ilr++) { // layer
62 for (int ich = 0; ich < o2::trd::constants::NSTACK * o2::trd::constants::NSECTOR; ich++) { // chamber
63 int isector = ich / o2::trd::constants::NSTACK;
64 int istack = ich % o2::trd::constants::NSTACK;
65 uint16_t sid = o2::trd::Geometry::getDetector(ilr, istack, isector);
66 const char* symname = Form("TRD/sm%02d/st%d/pl%d", isector, istack, ilr);
68 if (!gGeoManager->GetAlignableEntry(symname)) {
69 chamb->setDummy(true);
70 // continue;
71 }
72 if (!sect[isector]) {
73 sect[isector] = new AlignableVolume(Form("TRD/sm%02d", isector), getNonSensLabel(isector), mController);
74 sect[isector]->setParent(volTRD);
75 }
76 chamb->setParent(sect[isector]);
77 } // chamber
78 } // layer
79 //
80 for (int isc = 0; isc < o2::trd::constants::NSECTOR; isc++) {
81 if (sect[isc]) {
82 addVolume(sect[isc]);
83 }
84 }
85}
86
87//____________________________________________
88void AlignableDetectorTRD::Print(const Option_t* opt) const
89{
90 // print info
92 printf("Extra error for RC tracklets: Y:%e Z:%e\n", mExtraErrRC[0], mExtraErrRC[1]);
93}
94
95//____________________________________________
97{
98 // return calibration DOF name
99 return i < NCalibParams ? CalibDOFName[i] : nullptr;
100}
101
102//______________________________________________________
103void AlignableDetectorTRD::writePedeInfo(FILE* parOut, const Option_t* opt) const
104{
105 // contribute to params and constraints template files for PEDE
107 //
108 // write calibration parameters
109 enum { kOff,
110 kOn,
111 kOnOn };
112 const char* comment[3] = {" ", "! ", "!!"};
113 const char* kKeyParam = "parameter";
114 //
115 fprintf(parOut, "%s%s %s\t %d calibraction params for %s\n", comment[kOff], kKeyParam, comment[kOnOn],
116 getNCalibDOFs(), GetName());
117 //
118 for (int ip = 0; ip < getNCalibDOFs(); ip++) {
119 int cmt = isCondDOF(ip) ? kOff : kOn;
120 fprintf(parOut, "%s %9d %+e %+e\t%s %s p%d\n", comment[cmt], getParLab(ip),
121 -getParVal(ip), getParErr(ip), comment[kOnOn], isFreeDOF(ip) ? " " : "FX", ip);
122 }
123 //
124}
125
126//______________________________________________________
128{
129 //
131 //
132 for (int ip = 0; ip < getNCalibDOFs(); ip++) {
133 fprintf(parOut, "%9d %+e %+e\t! calib param %d of %s %s %s\n", getParLab(ip), -getParVal(ip), getParErr(ip), ip, GetName(),
134 isFreeDOF(ip) ? " " : "FXU", o2::align::utils::isZeroAbs(getParVal(ip)) ? "FXP" : " ");
135 }
136 //
137}
138
139//_______________________________________________________
141{
142 // return preset value of calibration dof
143 double val = 0;
144 switch (id) {
147 break;
148 case CalibDVT:
149 val = getCorrDVT();
150 break;
151 default:
152 break;
153 };
154 return val;
155}
156
157//_______________________________________________________
159{
160 // return preset value of calibration dof + mp correction
161 return getCalibDOFVal(id) + getParVal(id);
162}
163
164//____________________________________________
165int AlignableDetectorTRD::processPoints(GIndex gid, int npntCut, bool inv)
166{
167 // Extract the points corresponding to this detector, recalibrate/realign them to the
168 // level of the "starting point" for the alignment/calibration session.
169 // If inv==true, the track propagates in direction of decreasing tracking X
170 // (i.e. upper leg of cosmic track)
171 //
172 const auto& algConf = AlignConfig::Instance();
173 const auto recoData = mController->getRecoContainer();
174 const auto& trk = recoData->getTrack<o2::trd::TrackTRD>(gid);
175 if (trk.getNtracklets() < npntCut) {
176 return -1;
177 }
178 auto propagator = o2::base::Propagator::Instance(); // float version!
179 static bool firstCall = true;
180 if (firstCall) {
181 mRecoParam.init(propagator->getNominalBz());
182 firstCall = false;
183 }
184 const auto* transformer = mController->getTRDTransformer();
185 auto algTrack = mController->getAlgTrack();
186 const auto trackletsRaw = recoData->getTRDTracklets();
187 bool fail = false;
188 int nPntIni = algTrack->getNPoints(), npoints = 0;
189 o2::track::TrackPar trkParam = trk.getOuterParam(); // we refit outer param inward to get tracklet coordinates accounting for tilt
190 for (int il = o2::trd::constants::NLAYER; il--;) {
191 if (trk.getTrackletIndex(il) == -1) {
192 continue;
193 }
194 int trkltId = trk.getTrackletIndex(il);
195 const auto& trackletRaw = trackletsRaw[trkltId];
196 const auto trackletCalibLoc = transformer->transformTracklet(trackletRaw, false); // calibrated tracket in local frame !!!
197 int trkltDet = trackletRaw.getDetector();
198 auto* sensor = (AlignableSensorTRD*)getSensor(trkltDet);
199 if (sensor->isDummy()) {
200 LOGP(error, "Dummy sensor {} is referred by a track", trkltDet);
201 fail = true;
202 continue;
203 }
204 double locXYZ[3] = {trackletCalibLoc.getX(), trackletCalibLoc.getY(), trackletCalibLoc.getZ()}, locXYZC[3], traXYZ[3];
205 const auto& matAlg = sensor->getMatrixClAlg(); // local alignment matrix
206 matAlg.LocalToMaster(locXYZ, locXYZC); // aligned point in the local frame
207 const auto& mat = sensor->getMatrixT2L(); // RS FIXME check if correct
208 mat.MasterToLocal(locXYZC, traXYZ);
209 int trkltSec = sensor->getSector(); // trkltDet / (o2::trd::constants::NLAYER * o2::trd::constants::NSTACK);
210 float alpSens = sensor->getAlpTracking(); // o2::math_utils::sector2Angle(trkltSec);
211 if (trkltSec != o2::math_utils::angle2Sector(trkParam.getAlpha()) ||
212 !trkParam.rotateParam(alpSens) ||
213 !propagator->propagateTo(trkParam, traXYZ[0], propagator->getNominalBz(), o2::base::Propagator::MAX_SIN_PHI, 10., o2::base::Propagator::MatCorrType::USEMatCorrNONE)) { // we don't need high precision here
214 fail = true;
215 break;
216 }
217 const o2::trd::PadPlane* pad = o2::trd::Geometry::instance()->getPadPlane(trkltDet);
218 float tilt = std::tan(TMath::DegToRad() * pad->getTiltingAngle()); // tilt is signed! and returned in degrees
219 float tiltCorrUp = tilt * (traXYZ[2] - trkParam.getZ());
220 float zPosCorrUp = trkParam.getZ() + getNonRCCorrDzDtglWithCal() * trkParam.getTgl(); // + mRecoParam.getZCorrCoeffNRC() * trkParam.getTgl();
221 float padLength = pad->getRowSize(trackletRaw.getPadRow());
222 if (std::fabs(traXYZ[2] - trkParam.getZ()) < padLength) { // RS do we need this?
223 tiltCorrUp = 0.f;
224 }
225 std::array<float, 2> trkltPosUp{float(traXYZ[1] - tiltCorrUp), zPosCorrUp};
226 std::array<float, 3> trkltCovUp;
227 mRecoParam.recalcTrkltCov(tilt, trkParam.getSnp(), pad->getRowSize(trackletRaw.getPadRow()), trkltCovUp);
228 // Correction for DVT, equivalent to shift in X at which Y is evaluated: dY = tg_phi * dvt
229 {
230 auto dvt = getCorrDVTWithCal();
231 if (std::abs(dvt) > utils::AlmostZeroF) {
232 auto snp = trkParam.getSnp();
233 auto slpY = snp / std::sqrt((1.f - snp) * (1.f + snp));
234 trkltPosUp[0] += slpY * dvt;
235 }
236 }
237 auto& pnt = algTrack->addDetectorPoint();
238 const auto* sysE = sensor->getAddError(); // additional syst error
239 pnt.setYZErrTracking(trkltCovUp[0] + sysE[0] * sysE[0], trkltCovUp[1], trkltCovUp[2] + sysE[1] * sysE[1]);
240 if (getUseErrorParam()) { // errors will be calculated just before using the point in the fit, using track info
241 pnt.setNeedUpdateFromTrack();
242 }
243 pnt.setXYZTracking(traXYZ[0], trkltPosUp[0], trkltPosUp[1]);
244 pnt.setSensor(sensor);
245 pnt.setAlphaSens(alpSens);
246 pnt.setXSens(sensor->getXTracking());
247 pnt.setDetID(mDetID);
248 pnt.setSID(sensor->getSID());
249 pnt.setContainsMeasurement();
250 pnt.setInvDir(inv);
251 pnt.init();
252 npoints++;
253 }
254 if (fail) { // reset points to original start
255 algTrack->suppressLastPoints(npoints);
256 npoints = 0;
257 }
258 mNPoints += npoints;
259
260 return npoints;
261}
262
263} // namespace align
264} // namespace o2
Configuration file for global alignment.
TRD detector wrapper.
Base class of alignable volume.
Steering class for the global alignment.
Wrapper container for different reconstructed object types.
Global TRD definitions and constants.
int32_t i
void writePedeInfo(FILE *parOut, const Option_t *opt="") const final
void writeLabeledPedeResults(FILE *parOut) const final
const char * getCalibDOFName(int i) const final
double getCalibDOFVal(int id) const final
int processPoints(GIndex gid, int npntCut, bool inv) final
double getCalibDOFValWithCal(int id) const final
o2::gpu::GPUTRDRecoParam mRecoParam
void Print(const Option_t *opt="") const final
static const char * CalibDOFName[NCalibParams]
AlignableSensor * getSensor(int id) const
virtual void writeLabeledPedeResults(FILE *parOut) const
virtual void addVolume(AlignableVolume *vol)
virtual void writePedeInfo(FILE *parOut, const Option_t *opt="") const
void Print(const Option_t *opt="") const override
void setDummyEnvelope(bool v=true)
void setParent(AlignableVolume *par)
const o2::globaltracking::RecoContainer * getRecoContainer() const
Definition Controller.h:200
AlignmentTrack * getAlgTrack() const
Definition Controller.h:198
const o2::trd::TrackletTransformer * getTRDTransformer() const
Definition Controller.h:270
Controller * mController
Definition DOFSet.h:73
float getParVal(int par) const
Definition DOFSet.h:39
int getNCalibDOFs() const
Definition DOFSet.h:46
float getParErr(int par) const
Definition DOFSet.h:40
int getParLab(int par) const
Definition DOFSet.h:41
static int getSensID(o2::detectors::DetID detid, int sensid)
static constexpr float MAX_SIN_PHI
Definition Propagator.h:72
GPUd() value_type estimateLTFast(o2 static GPUd() float estimateLTIncrement(const o2 PropagatorImpl * Instance(bool uninitialized=false)
Definition Propagator.h:178
Static class with identifiers, bitmasks and names for ALICE detectors.
Definition DetID.h:58
void init(float bz, const GPUSettingsRec *rec=nullptr)
Load parameterization for given magnetic field.
static Geometry * instance()
Definition Geometry.h:33
GLuint GLfloat * val
Definition glcorearb.h:1582
constexpr bool isZeroAbs(double d) noexcept
Definition utils.h:63
constexpr float AlmostZeroF
Definition utils.h:39
int angle2Sector(float phi)
Definition Utils.h:183
void align(gsl::span< ElinkEncoder< BareFormat, CHARGESUM > > elinks)
constexpr int NLAYER
the number of layers
Definition Constants.h:27
constexpr int NSECTOR
the number of sectors
Definition Constants.h:25
constexpr int NSTACK
the number of stacks per sector
Definition Constants.h:26
a couple of static helper functions to create timestamp values for CCDB queries or override obsolete ...
const U & getTrack(int src, int id) const