Project
Loading...
Searching...
No Matches
AlignableVolume.h
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
24#ifndef ALIGNABLEVOLUME_H
25#define ALIGNABLEVOLUME_H
26
27#include <TNamed.h>
28#include <TObjArray.h>
29#include <TGeoMatrix.h>
30#include <cstdio>
32#include "Align/DOFSet.h"
33#include <vector>
34
35class TObjArray;
36class TH1;
37
38namespace o2
39{
40namespace align
41{
42
43class Controller;
44
45class AlignableVolume : public DOFSet
46{
47 public:
56 enum { kDOFBitTX = BIT(kDOFTX),
61 kDOFBitPH = BIT(kDOFPH) };
62 enum { kNDOFMax = 32 };
63
64 enum Frame_t { kLOC,
66 kNVarFrames }; // variation frames defined
67 enum { kInitDOFsDoneBit = BIT(14),
68 kSkipBit = BIT(15),
70 enum { kDefChildConstr = 0xff };
71 //
72 AlignableVolume() = default;
73 AlignableVolume(const char* symname, int iid, Controller* ctr);
74 ~AlignableVolume() override;
75 //
76 const char* getSymName() const { return GetName(); }
77 //
78 int getVolID() const { return (int)GetUniqueID(); }
79 void setVolID(int v) { SetUniqueID(v); }
80 int getInternalID() const { return mIntID; }
81 void setInternalID(int v) { mIntID = v; }
82 //
83 //
84 void assignDOFs();
85 void initDOFs();
86 //
87 void getParValGeom(double* delta) const;
88
89 Frame_t getVarFrame() const { return mVarFrame; }
91 bool isFrameTRA() const { return mVarFrame == kTRA; }
92 bool isFrameLOC() const { return mVarFrame == kLOC; }
93 //
94 void setFreeDOF(int dof)
95 {
96 mDOF |= 0x1 << dof;
97 calcFree();
98 }
99 void fixDOF(int dof)
100 {
101 mDOF &= ~(0x1 << dof);
102 calcFree();
103 }
104 void setFreeDOFPattern(uint32_t pat)
105 {
106 mDOF = pat;
107 calcFree();
108 }
109
110 void setMeasuredDOFPattern(uint32_t pat)
111 {
112 mDOFAsMeas = pat;
113 }
114 bool isNameMatching(const std::string& regexStr) const;
115 bool isFreeDOF(int dof) const { return (mDOF & (0x1 << dof)) != 0; }
116 bool isMeasuredDOF(int dof) const { return isFreeDOF(dof) && ((mDOFAsMeas & (0x1 << dof)) != 0); }
117 bool isCondDOF(int dof) const;
118 uint32_t getFreeDOFPattern() const { return mDOF; }
119 uint32_t getFreeDOFGeomPattern() const { return mDOF & kAllGeomDOF; }
120 //
121 void addAutoConstraints();
122 bool isChildrenDOFConstrained(int dof) const { return mConstrChild & 0x1 << dof; }
123 uint8_t getChildrenConstraintPattern() const { return mConstrChild; }
124 void constrainChildrenDOF(int dof) { mConstrChild |= 0x1 << dof; }
125 void unConstrainChildrenDOF(int dof) { mConstrChild &= ~(0x1 << dof); }
126 void setChildrenConstrainPattern(uint32_t pat) { mConstrChild = pat; }
127 bool hasChildrenConstraint() const { return mConstrChild; }
128 //
129 AlignableVolume* getParent() const { return mParent; }
131 {
132 mParent = par;
133 if (par) {
134 par->addChild(this);
135 }
136 }
137 int countParents() const;
138 //
139 int getNChildren() const { return mChildren ? mChildren->GetEntriesFast() : 0; }
140 AlignableVolume* getChild(int i) const { return mChildren ? (AlignableVolume*)mChildren->UncheckedAt(i) : nullptr; }
141 virtual void addChild(AlignableVolume* ch);
142 //
143 double getXTracking() const { return mX; }
144 double getAlpTracking() const { return mAlp; }
145 //
146 int getNProcessedPoints() const { return mNProcPoints; }
147 virtual int finalizeStat();
148 //
149 int getNDOFGeomFree() const { return mNDOFGeomFree; }
150 //
151 virtual void prepareMatrixT2L();
152 //
153 const TGeoHMatrix& getMatrixL2G() const { return mMatL2G; }
154 const TGeoHMatrix& getMatrixL2GIdeal() const { return mMatL2GIdeal; }
155 const TGeoHMatrix& getMatrixL2GReco() const { return mMatL2GReco; }
156 const TGeoHMatrix& getGlobalDeltaRef() const { return mMatDeltaRefGlo; }
157 const TGeoHMatrix& getMatrixT2L() const { return mMatT2L; }
158
159 void setMatrixL2G(const TGeoHMatrix& m) { mMatL2G = m; }
160 void setMatrixL2GIdeal(const TGeoHMatrix& m) { mMatL2GIdeal = m; }
161 void setMatrixL2GReco(const TGeoHMatrix& m) { mMatL2GReco = m; }
162 void setGlobalDeltaRef(const TGeoHMatrix& mat) { mMatDeltaRefGlo = mat; }
163 void setMatrixT2L(const TGeoHMatrix& m) { mMatT2L = m; }
164
165 //
166 virtual void prepareMatrixL2G(bool reco = false);
167 virtual void prepareMatrixL2GIdeal();
168 virtual void updateL2GRecoMatrices(const std::vector<o2::detectors::AlignParam>& algArr, const TGeoHMatrix* cumulDelta);
169 //
170 void getMatrixT2G(TGeoHMatrix& m) const;
171 //
172 void delta2Matrix(TGeoHMatrix& deltaM, const double* delta) const;
173 //
174 // preparation of variation matrices
175 void getDeltaT2LmodLOC(TGeoHMatrix& matMod, const double* delta) const;
176 void getDeltaT2LmodTRA(TGeoHMatrix& matMod, const double* delta) const;
177 void getDeltaT2LmodLOC(TGeoHMatrix& matMod, const double* delta, const TGeoHMatrix& relMat) const;
178 void getDeltaT2LmodTRA(TGeoHMatrix& matMod, const double* delta, const TGeoHMatrix& relMat) const;
179 //
180 // creation of global matrices for storage
181 bool createGloDeltaMatrix(TGeoHMatrix& deltaM) const;
182 bool createLocDeltaMatrix(TGeoHMatrix& deltaM) const; // return true if the matrix is not unity
183 void createPreGloDeltaMatrix(TGeoHMatrix& deltaM) const;
184 void createPreLocDeltaMatrix(TGeoHMatrix& deltaM) const;
185 void createAlignmenMatrix(TGeoHMatrix& alg, const TGeoHMatrix* envelopeDelta = nullptr) const;
186 void createAlignmentObjects(std::vector<o2::detectors::AlignParam>& arr, const TGeoHMatrix* envelopeDelta = nullptr) const;
187 //
188 void setSkip(bool v = true) { SetBit(kSkipBit, v); }
189 bool getSkip() const { return TestBit(kSkipBit); }
190 //
193 //
195 bool getInitDOFsDone() const { return TestBit(kInitDOFsDoneBit); }
196 //
197 bool ownsDOFID(int id) const;
198 AlignableVolume* getVolOfDOFID(int id) const;
199 //
200 bool isDummyEnvelope() const { return mIsDummyEnvelope; }
201 void setDummyEnvelope(bool v = true) { mIsDummyEnvelope = v; }
202 //
203 bool isDummy() const { return mIsDummy; }
204 void setDummy(bool v) { mIsDummy = v; }
205 //
206 virtual bool isSensor() const { return false; }
207 //
208 virtual const char* getDOFName(int i) const;
209 void Print(const Option_t* opt = "") const override;
210 virtual void writePedeInfo(FILE* parOut, const Option_t* opt = "") const;
211 virtual void writeLabeledPedeResults(FILE* parOut) const;
212 //
213 static const char* getGeomDOFName(int i) { return i < kNDOFGeom ? sDOFName[i] : nullptr; }
214 static void setDefGeomFree(uint8_t patt) { sDefGeomFree = patt; }
215 static uint8_t getDefGeomFree() { return sDefGeomFree; }
216 //
217 protected:
218 void calcFree(bool condFree = true);
219 //
220 // ------- dummies -------
223 //
224 protected:
225 //
226 Frame_t mVarFrame = kTRA; // Variation frame for this volume
227 int mIntID = -1; // internal id within the detector
228 double mX = 0.; // tracking frame X offset
229 double mAlp = 0.; // tracking frame alpa
230 //
231 uint32_t mDOF = 0; // pattern of DOFs
232 uint32_t mDOFAsMeas = 0; // consider DOF as measured with presigma error
233 bool mIsDummy = false; // placeholder (e.g. inactive), used to have the numbering corresponding to position in the container
234 bool mIsDummyEnvelope = false; // some volumes are dummy envelopes for their children
235 //
236 char mNDOFGeomFree = 0; // number of free geom degrees of freedom
237 uint8_t mConstrChild = 0; // bitpattern for constraints on children corrections
238 //
239 AlignableVolume* mParent = nullptr; // parent volume
240 TObjArray* mChildren = nullptr; // array of childrens
241 //
242 int mNProcPoints = 0; // n of processed points
243
244 TGeoHMatrix mMatL2GReco{}; // local to global matrix used for reco of data being processed
245 TGeoHMatrix mMatL2G{}; // local to global matrix, including current alignment
246 TGeoHMatrix mMatL2GIdeal{}; // local to global matrix, ideal
247 TGeoHMatrix mMatT2L{}; // tracking to local matrix (ideal)
248 TGeoHMatrix mMatDeltaRefGlo{}; // global reference delta from Align/Data
249 //
250 static const char* sDOFName[kNDOFGeom];
251 static const char* sFrameName[kNVarFrames];
252 static uint32_t sDefGeomFree;
253 //
255};
256
257//___________________________________________________________
258inline void AlignableVolume::getMatrixT2G(TGeoHMatrix& m) const
259{
260 // compute tracking to global matrix, i.e. glo = T2G*tra = L2G*loc = L2G*T2L*tra
262 m *= getMatrixT2L();
263}
264} // namespace align
265} // namespace o2
266#endif
Definition of the base alignment parameters class.
Interface to contiguous set of DOFs in the controller class.
int32_t i
void calcFree(bool condFree=true)
void setDummyEnvelope(bool v=true)
virtual void prepareMatrixL2G(bool reco=false)
uint32_t getFreeDOFGeomPattern() const
virtual void updateL2GRecoMatrices(const std::vector< o2::detectors::AlignParam > &algArr, const TGeoHMatrix *cumulDelta)
void createPreGloDeltaMatrix(TGeoHMatrix &deltaM) const
bool isNameMatching(const std::string &regexStr) const
static const char * sFrameName[kNVarFrames]
AlignableVolume * getChild(int i) const
const TGeoHMatrix & getMatrixT2L() const
virtual void writePedeInfo(FILE *parOut, const Option_t *opt="") const
bool isChildrenDOFConstrained(int dof) const
virtual bool isSensor() const
static void setDefGeomFree(uint8_t patt)
void Print(const Option_t *opt="") const override
static const char * sDOFName[kNDOFGeom]
void getDeltaT2LmodLOC(TGeoHMatrix &matMod, const double *delta) const
void setMeasuredDOFPattern(uint32_t pat)
bool isMeasuredDOF(int dof) const
void delta2Matrix(TGeoHMatrix &deltaM, const double *delta) const
void setMatrixL2GIdeal(const TGeoHMatrix &m)
void getParValGeom(double *delta) const
void getDeltaT2LmodTRA(TGeoHMatrix &matMod, const double *delta) const
bool getExcludeFromParentConstraint() const
ClassDefOverride(AlignableVolume, 2)
bool isCondDOF(int dof) const
const TGeoHMatrix & getMatrixL2G() const
void setGlobalDeltaRef(const TGeoHMatrix &mat)
virtual void writeLabeledPedeResults(FILE *parOut) const
const TGeoHMatrix & getMatrixL2GReco() const
bool isFreeDOF(int dof) const
bool createGloDeltaMatrix(TGeoHMatrix &deltaM) const
uint32_t getFreeDOFPattern() const
virtual void addChild(AlignableVolume *ch)
uint8_t getChildrenConstraintPattern() const
bool createLocDeltaMatrix(TGeoHMatrix &deltaM) const
void createAlignmenMatrix(TGeoHMatrix &alg, const TGeoHMatrix *envelopeDelta=nullptr) const
void setParent(AlignableVolume *par)
void setMatrixL2GReco(const TGeoHMatrix &m)
void getMatrixT2G(TGeoHMatrix &m) const
static const char * getGeomDOFName(int i)
AlignableVolume * getVolOfDOFID(int id) const
AlignableVolume & operator=(const AlignableVolume &)
void setChildrenConstrainPattern(uint32_t pat)
AlignableVolume(const AlignableVolume &)
void setMatrixL2G(const TGeoHMatrix &m)
void createAlignmentObjects(std::vector< o2::detectors::AlignParam > &arr, const TGeoHMatrix *envelopeDelta=nullptr) const
const TGeoHMatrix & getMatrixL2GIdeal() const
AlignableVolume * getParent() const
void createPreLocDeltaMatrix(TGeoHMatrix &deltaM) const
const TGeoHMatrix & getGlobalDeltaRef() const
void setMatrixT2L(const TGeoHMatrix &m)
const char * getSymName() const
static uint8_t getDefGeomFree()
virtual const char * getDOFName(int i) const
void excludeFromParentConstraint(bool v=true)
void setFreeDOFPattern(uint32_t pat)
const GLfloat * m
Definition glcorearb.h:4066
const GLdouble * v
Definition glcorearb.h:832
GLdouble f
Definition glcorearb.h:310
void align(gsl::span< ElinkEncoder< BareFormat, CHARGESUM > > elinks)
a couple of static helper functions to create timestamp values for CCDB queries or override obsolete ...