Project
Loading...
Searching...
No Matches
Aligner.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
17#ifndef ALICEO2_MCH_ALIGNER
18#define ALICEO2_MCH_ALIGNER
19
20#include <string>
21#include <vector>
22
28#include "MCHTracking/Track.h"
29
30#include <TFile.h>
31#include <TGeoManager.h>
32#include <TGeoMatrix.h>
33#include <TObject.h>
34#include <TString.h>
35#include <TTree.h>
36
37namespace o2
38{
39
40namespace mch
41{
42
45{
46 public:
47 //* construction
48 LocalTrackParam() = default;
49 ~LocalTrackParam() = default;
50
51 // private:
52 //* y and z
53 double fTrackX = 0.0;
54 double fTrackY = 0.0;
55 double fTrackZ = 0.0;
56 double fTrackSlopeX = 0.0;
57 double fTrackSlopeY = 0.0;
58}; // class LocalTrackParam
59
62{
63 public:
64 //* construction
67
68 // private:
69 //* y and z
70 int fClDetElem = 0.0;
72
73 float fClusterX = 0.0;
74 float fClusterY = 0.0;
75 float fClusterZ = 0.0;
76 float fClusterXloc = 0.0;
77 float fClusterYloc = 0.0;
78
79 float fTrackX = 0.0;
80 float fTrackY = 0.0;
81 float fTrackZ = 0.0;
82 float fTrackXloc = 0.0;
83 float fTrackYloc = 0.0;
84
85 float fTrackSlopeX = 0.0;
86 float fTrackSlopeY = 0.0;
87
88 float fBendingMomentum = 0.0;
89
90 float fResiduXGlobal = 0.0;
91 float fResiduYGlobal = 0.0;
92
93 float fResiduXLocal = 0.0;
94 float fResiduYLocal = 0.0;
95
96 float fCharge = 0.0;
97
98 float fBx = 0.0;
99 float fBy = 0.0;
100 float fBz = 0.0;
101
102}; // class LocalTrackClusterResidual
103
104class Aligner : public TObject
105{
106
107 public:
108 Aligner();
109
110 ~Aligner();
111
112 // initialize
113 void init(TString DataRecFName = "millerecords.root", TString ConsRecFName = "milleconstraints.root");
114
115 // terminate
116 void terminate();
117
118 // array dimendions
119 enum {
121 fgNSt = 5,
122
124 fgNCh = 10,
125
128
131
134
138
140 fNLocal = 4, // t_x, t_y, x0, y0
141
143 fgNParCh = 4, // x,y,z,phi
144
147 };
148
150 static const int fgNDetElemCh[fgNCh];
151
153 static const int fgSNDetElemCh[fgNCh + 1];
154
156 static const int fgNDetElemHalfCh[fgNHalfCh];
157
160
173
186
187 void ProcessTrack(Track& track, const o2::mch::geo::TransformationCreator& transformation, Bool_t doAlignment, Double_t weight = 1);
188
189 //@name modifiers
191
193 void SetRunNumber(int id)
194 {
195 fRunNumber = id;
196 }
197
200 {
201 fBFieldOn = value;
202 }
203
206 {
207 fDoEvaluation = value;
208 }
209
212 {
213 fRefitStraightTracks = value;
214 }
215
216 void SetAllowedVariation(int iPar, double value);
217
218 void SetSigmaXY(double sigmaX, double sigmaY);
219
220 void FixAll(unsigned int parameterMask = ParAll);
221
222 void FixChamber(int iCh, unsigned int parameterMask = ParAll);
223
224 void FixDetElem(int iDetElemId, unsigned int parameterMask = ParAll);
225
226 void FixHalfSpectrometer(const bool* bChOnOff, unsigned int sidesMask = AllSides, unsigned int parameterMask = ParAll);
227
228 void FixParameter(int iPar);
229
230 void FixParameter(int iDetElem, int iPar)
231 {
232 FixParameter(iDetElem * fgNParCh + iPar);
233 }
234
236
237 //@name releasing detectors
239
240 void ReleaseChamber(int iCh, unsigned int parameterMask = ParAll);
241
242 void ReleaseDetElem(int iDetElemId, unsigned int parameterMask = ParAll);
243
244 void ReleaseParameter(int iPar);
245
246 void ReleaseParameter(int iDetElem, int iPar)
247 {
248 ReleaseParameter(iDetElem * fgNParCh + iPar);
249 }
250
252
253 //@name grouping detectors
255
256 void GroupChamber(int iCh, unsigned int parameterMask = ParAll);
257
258 void GroupHalfChamber(int iCh, int iHalf, unsigned int parameterMask = ParAll);
259
260 void GroupDetElems(int detElemMin, int detElemMax, unsigned int parameterMask = ParAll);
261
262 void GroupDetElems(const int* detElemList, int nDetElem, unsigned int parameterMask = ParAll);
263
265
266 //@name define non linearity
268
269 void SetChamberNonLinear(int iCh, unsigned int parameterMask);
270
271 void SetDetElemNonLinear(int iSt, unsigned int parameterMask);
272
273 void SetParameterNonLinear(int iPar);
274
275 void SetParameterNonLinear(int iDetElem, int iPar)
276 {
277 SetParameterNonLinear(iDetElem * fgNParCh + iPar);
278 }
279
281
282 //@name constraints
284
285 void AddConstraints(const bool* bChOnOff, unsigned int parameterMask);
286
287 void AddConstraints(const bool* bChOnOff, const bool* lVarXYT, unsigned int sidesMask = AllSides);
288
290
292 void InitGlobalParameters(double* par);
293
295 void GlobalFit(std::vector<double>& params, std::vector<double>& errors, std::vector<double>& pulls);
296
298 void PrintGlobalParameters(void) const;
299
301 double GetParError(int iPar) const;
302
303 o2::fwdalign::MillePedeRecord& GetRecord() { return fTrackRecord; }
304
305 void ReAlign(std::vector<o2::detectors::AlignParam>& params, std::vector<double>& misAlignments);
306
307 void SetAlignmentResolution(const TClonesArray* misAlignArray, int chId, double chResX, double chResY, double deResX, double deResY);
308
309 TTree* GetResTree()
310 {
311 return fTTree;
312 }
313
315 {
316 mRead = true;
317 }
318
320 {
321 fDisableRecordWriter = true;
322 }
323
324 private:
326 Aligner(const Aligner& right);
327
329 Aligner& operator=(const Aligner& right);
330
332 void SetLocalDerivative(int index, double value)
333 {
334 fLocalDerivatives[index] = value;
335 }
336
338 void SetGlobalDerivative(int index, double value)
339 {
340 fGlobalDerivatives[index] = value;
341 }
342
344 LocalTrackParam RefitStraightTrack(Track&, double) const;
345
346 void FillDetElemData(const Cluster*);
347
348 void FillRecPointData(const Cluster*);
349
350 void FillTrackParamData(const TrackParam*);
351
352 void LocalEquationX(const double* r);
353
354 void LocalEquationY(const double* r);
355
356 TGeoCombiTrans DeltaTransform(const double* detElemMisAlignment) const;
357
358 bool isMatrixConvertedToAngles(const double* rot, double& psi, double& theta, double& phi) const;
359
361
362
363 void AddConstraint(double* parameters, double value);
364
365 int GetChamberId(int iDetElemNumber) const;
366
367 bool DetElemIsValid(int iDetElemId) const;
368
369 int GetDetElemNumber(int iDetElemId) const;
370
371 TString GetParameterMaskString(unsigned int parameterMask) const;
372
373 TString GetSidesMaskString(unsigned int sidesMask) const;
374
376
378 bool fInitialized;
379
381 int fRunNumber;
382
384 bool fBFieldOn;
385
387 bool fRefitStraightTracks;
388
390 double fAllowVar[fgNParCh];
391
393 double fStartFac;
394
396 double fResCutInitial;
397
399 double fResCut;
400
402 o2::fwdalign::MillePede2* fMillepede; // AliMillePede2 implementation
403
405 int fNStdDev;
406
408 double fClustPos[3];
409
411 double fTrackSlope0[2];
412
414 double fTrackSlope[2];
415
417 double fTrackPos0[3];
418
420 double fTrackPos[3];
421
423 double fMeas[2];
424
426 double fSigma[2];
427
429 enum {
430 kFixedParId = -1,
431 kFreeParId = kFixedParId - 1,
432 kGroupBaseId = -10
433 };
434
437 std::vector<int> fGlobalParameterStatus;
438
440 std::vector<double> fGlobalDerivatives;
441
443 std::vector<double> fLocalDerivatives;
444
446 int fDetElemNumber;
447
450
452 o2::mch::geo::TransformationCreator fTransformCreator;
453
455 bool fDoEvaluation;
456
458 bool fDisableRecordWriter;
459
460 LocalTrackClusterResidual* fTrkClRes;
461
463 TFile* fTFile;
464
466 TTree* fTTree;
467
469 bool mRead;
470
471 long mNEntriesAutoSave = 10000;
472
473 o2::fwdalign::MilleRecordWriter* mRecordWriter;
474 bool mWithConstraintsRecWriter;
475 o2::fwdalign::MilleRecordWriter* mConstraintsRecWriter;
476
477 o2::fwdalign::MilleRecordReader* mRecordReader;
478 bool mWithConstraintsRecReader = false;
479 o2::fwdalign::MilleRecordReader* mConstraintsRecReader;
480
481}; // class Alignment
482
483} // namespace mch
484} // namespace o2
485#endif // ALICEO2_MCH_ALIGNER_H_
Definition of the base alignment parameters class.
Definition of the MCH cluster minimal structure.
Definition of the MCH track for internal use.
General class for alignment with large number of degrees of freedom, adapted from AliROOT.
Class to store the data of single track processing.
Store residuals and local/global deriavtives from a single track processing.
HMPID cluster implementation.
Definition Cluster.h:27
void ReleaseDetElem(int iDetElemId, unsigned int parameterMask=ParAll)
Definition Aligner.cxx:742
void SetRefitStraightTracks(bool value)
set to true to refit tracks
Definition Aligner.h:211
void SetBFieldOn(bool value)
Set flag for Magnetic field On/Off.
Definition Aligner.h:199
void SetAllowedVariation(int iPar, double value)
Definition Aligner.cxx:1267
void ReAlign(std::vector< o2::detectors::AlignParam > &params, std::vector< double > &misAlignments)
Definition Aligner.cxx:1322
void SetDetElemNonLinear(int iSt, unsigned int parameterMask)
Definition Aligner.cxx:878
void SetRunNumber(int id)
run number
Definition Aligner.h:193
static const int fgNDetElemCh[fgNCh]
Number of detection elements per chamber.
Definition Aligner.h:150
void FixParameter(int iPar)
Definition Aligner.cxx:698
o2::fwdalign::MillePedeRecord & GetRecord()
Definition Aligner.h:303
void SetAlignmentResolution(const TClonesArray *misAlignArray, int chId, double chResX, double chResY, double deResX, double deResY)
Definition Aligner.cxx:1404
void SetReadOnly()
Definition Aligner.h:314
void GlobalFit(std::vector< double > &params, std::vector< double > &errors, std::vector< double > &pulls)
perform global fit
Definition Aligner.cxx:1298
static const int fgDetElemHalfCh[fgNHalfCh][fgNDetHalfChMax]
list of detection elements per tracking module
Definition Aligner.h:159
@ fNGlobal
Number of global parameters.
Definition Aligner.h:146
@ fgNDetHalfChMax
max number of detector elements per half chamber
Definition Aligner.h:133
@ fgNSt
Number tracking stations.
Definition Aligner.h:121
@ fgNHalfCh
Number of half chambers.
Definition Aligner.h:130
@ fNLocal
Number of local parameters.
Definition Aligner.h:140
@ fgNCh
Number tracking chambers.
Definition Aligner.h:124
@ fgNParCh
Number of degrees of freedom per chamber.
Definition Aligner.h:143
@ fgNTrkMod
Number of tracking modules.
Definition Aligner.h:127
void AddConstraints(const bool *bChOnOff, unsigned int parameterMask)
Definition Aligner.cxx:909
void SetDoEvaluation(bool value)
set to true to do refit evaluation
Definition Aligner.h:205
void SetParameterNonLinear(int iPar)
Definition Aligner.cxx:897
void ReleaseParameter(int iDetElem, int iPar)
Definition Aligner.h:246
void FixHalfSpectrometer(const bool *bChOnOff, unsigned int sidesMask=AllSides, unsigned int parameterMask=ParAll)
Definition Aligner.cxx:628
void SetChamberNonLinear(int iCh, unsigned int parameterMask)
Definition Aligner.cxx:855
void init(TString DataRecFName="millerecords.root", TString ConsRecFName="milleconstraints.root")
Definition Aligner.cxx:195
void FixChamber(int iCh, unsigned int parameterMask=ParAll)
Definition Aligner.cxx:577
void PrintGlobalParameters(void) const
print global parameters
Definition Aligner.cxx:1310
SidesMask
detector sides bit set, used for selecting sides in constrains
Definition Aligner.h:175
static const int fgSNDetElemCh[fgNCh+1]
Sum of detection elements up to this chamber.
Definition Aligner.h:153
void ProcessTrack(Track &track, const o2::mch::geo::TransformationCreator &transformation, Bool_t doAlignment, Double_t weight=1)
Definition Aligner.cxx:376
double GetParError(int iPar) const
get error on a given parameter
Definition Aligner.cxx:1316
void ReleaseChamber(int iCh, unsigned int parameterMask=ParAll)
Definition Aligner.cxx:710
void GroupDetElems(int detElemMin, int detElemMax, unsigned int parameterMask=ParAll)
Definition Aligner.cxx:802
ParameterMask
global parameter bit set, used for masks
Definition Aligner.h:162
void InitGlobalParameters(double *par)
initialize global parameters to a give set of values
Definition Aligner.cxx:1256
void FixDetElem(int iDetElemId, unsigned int parameterMask=ParAll)
Definition Aligner.cxx:609
void GroupChamber(int iCh, unsigned int parameterMask=ParAll)
Definition Aligner.cxx:773
void FixAll(unsigned int parameterMask=ParAll)
Definition Aligner.cxx:554
void SetParameterNonLinear(int iDetElem, int iPar)
Definition Aligner.h:275
void FixParameter(int iDetElem, int iPar)
Definition Aligner.h:230
void SetSigmaXY(double sigmaX, double sigmaY)
Definition Aligner.cxx:1284
void ReleaseParameter(int iPar)
Definition Aligner.cxx:761
TTree * GetResTree()
Definition Aligner.h:309
void DisableRecordWriter()
Definition Aligner.h:319
void GroupHalfChamber(int iCh, int iHalf, unsigned int parameterMask=ParAll)
Definition Aligner.cxx:786
static const int fgNDetElemHalfCh[fgNHalfCh]
Number of detection element per tracking module.
Definition Aligner.h:156
local track residual, for tempoarary eval
Definition Aligner.h:62
local track parameters, for refit
Definition Aligner.h:45
track for internal use
Definition Track.h:33
GLuint index
Definition glcorearb.h:781
GLdouble GLdouble right
Definition glcorearb.h:4077
GLuint GLuint GLfloat weight
Definition glcorearb.h:5477
GLsizei const GLfloat * value
Definition glcorearb.h:819
GLenum const GLfloat * params
Definition glcorearb.h:272
GLboolean r
Definition glcorearb.h:1233
GLuint id
Definition glcorearb.h:650
std::function< o2::math_utils::Transform3D(int detElemId)> TransformationCreator
a couple of static helper functions to create timestamp values for CCDB queries or override obsolete ...
auto transformation