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