]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/AliMUONAlignment.h
Added documentation about how to operate the Alignment task
[u/mrichter/AliRoot.git] / MUON / AliMUONAlignment.h
1 #ifndef ALIMUONALIGNMENT_H
2 #define ALIMUONALIGNMENT_H
3 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * See cxx source for full Copyright notice                               */
5
6 /* $Id: AliMUONAlignment.h 50858 2011-07-29 10:58:20Z ivana $ */
7
8 /// \ingroup rec
9 /// \class AliMUONAlignment
10 /// \brief Class for alignment of muon spectrometer
11 //
12 // Authors: Bruce Becker, Javier Castillo, Hugo Pereira Da Costa
13
14 #include "AliMillePedeRecord.h"
15
16 #include <TObject.h>
17 #include <TString.h>
18 #include <TGeoMatrix.h>
19
20 class TClonesArray;
21 class AliMillePede2;
22 class AliMUONGeometryTransformer;
23 class AliMUONTrack;
24 class AliMUONTrackParam;
25 class AliMUONVCluster;
26
27 class AliMUONAlignment:public TObject
28 {
29
30   public:
31
32   AliMUONAlignment();
33
34   virtual ~AliMUONAlignment();
35
36   void Init( void );
37
38   // array dimendions
39   enum
40   {
41     /// Number tracking stations
42     fgNSt = 5,
43
44     /// Number tracking chambers
45     fgNCh = 10,
46
47     /// Number of tracking modules (4 ch + 6*2 half-ch)
48     fgNTrkMod = 16,
49
50     /// Total number of detection elements
51     /// (4*2 + 4*2 + 18*2 + 26*2 + 26*2)
52     fgNDetElem = 156,
53
54     /// Number of local parameters
55     fNLocal = 4,
56
57     /// Number of degrees of freedom per chamber
58     fgNParCh = 4,
59
60     /// Number of global parameters
61     fNGlobal = fgNParCh*fgNDetElem
62   };
63
64   /// Number of detection elements per chamber
65   static const Int_t fgNDetElemCh[fgNCh];
66
67   /// Sum of detection elements up to this chamber
68   static const Int_t fgSNDetElemCh[fgNCh+1];
69
70   /// global parameter bit set, used for masks
71   enum ParameterMask
72   {
73     ParX = 1<<0,
74     ParY = 1<<1,
75     ParZ = 1<<2,
76     ParTZ = 1<<3,
77
78     ParAllTranslations = ParX|ParY|ParZ,
79     ParAllRotations = ParTZ,
80     ParAll = ParAllTranslations|ParAllRotations
81
82   };
83
84   /// detector sides bit set, used for selecting sides in constrains
85   enum SidesMask
86   {
87     SideTop = 1<<0,
88     SideLeft = 1<<1,
89     SideBottom = 1<<2,
90     SideRight = 1<<3,
91     SideTopLeft = SideTop|SideLeft,
92     SideTopRight = SideTop|SideRight,
93     SideBottomLeft = SideBottom|SideLeft,
94     SideBottomRight = SideBottom|SideRight,
95     AllSides = SideTop|SideBottom|SideLeft|SideRight
96   };
97
98   AliMillePedeRecord* ProcessTrack( AliMUONTrack* track, Bool_t doAlignment, Double_t weight = 1 );
99
100   void ProcessTrack( AliMillePedeRecord* );
101
102   //@name modifiers
103   //@{
104
105   /// run number
106   void SetRunNumber( Int_t id )
107   {fRunNumber = id;}
108
109   /// Set flag for Magnetic field On/Off
110   void SetBFieldOn( Bool_t value )
111   { fBFieldOn = value; }
112
113   /// use unbiased residuals
114   void SetUnbias( Bool_t value )
115   { fUnbias = value; }
116
117   void SetAllowedVariation( Int_t iPar, Double_t value );
118
119   void SetSigmaXY(Double_t sigmaX, Double_t sigmaY);
120
121   /// Set geometry transformer
122   void SetGeometryTransformer( AliMUONGeometryTransformer * transformer )
123   { fTransform = transformer; }
124
125   //@}
126
127   //@name fixing detectors
128   //@{
129
130   void FixAll( UInt_t parameterMask = ParAll );
131
132   void FixChamber( Int_t iCh, UInt_t parameterMask = ParAll );
133
134   void FixDetElem( Int_t iDetElemId, UInt_t parameterMask = ParAll );
135
136   void FixHalfSpectrometer( const Bool_t* bChOnOff, UInt_t sidesMask = AllSides, UInt_t parameterMask = ParAll );
137
138   void FixParameter( Int_t iPar );
139
140   void FixParameter( Int_t iDetElem, Int_t iPar )
141   { FixParameter( iDetElem*fgNParCh + iPar ); }
142
143   //@}
144
145   //@name releasing detectors
146   //@{
147
148   void ReleaseChamber( Int_t iCh, UInt_t parameterMask = ParAll );
149
150   void ReleaseDetElem( Int_t iDetElemId, UInt_t parameterMask = ParAll );
151
152   void ReleaseParameter( Int_t iPar );
153
154   void ReleaseParameter( Int_t iDetElem, Int_t iPar )
155   { ReleaseParameter( iDetElem*fgNParCh + iPar ); }
156
157   //@}
158
159   //@name grouping detectors
160   //@{
161
162   void GroupChamber( Int_t iCh, UInt_t parameterMask = ParAll );
163
164   void GroupDetElems( Int_t detElemMin, Int_t detElemMax, UInt_t parameterMask = ParAll );
165
166   void GroupDetElems( Int_t *detElemList, Int_t nDetElem, UInt_t parameterMask = ParAll );
167
168   //@}
169
170   //@name define non linearity
171   //@{
172
173   void SetChamberNonLinear( Int_t iCh, UInt_t parameterMask);
174
175   void SetDetElemNonLinear( Int_t iSt, UInt_t parameterMask);
176
177   void SetParameterNonLinear( Int_t iPar );
178
179   void SetParameterNonLinear( Int_t iDetElem, Int_t iPar )
180   { SetParameterNonLinear( iDetElem*fgNParCh + iPar ); }
181
182   //@}
183
184   //@name constraints
185   //@{
186
187   void AddConstraints( const Bool_t *bChOnOff, UInt_t parameterMask );
188
189   void AddConstraints( const Bool_t *bChOnOff, const Bool_t *lVarXYT, UInt_t sidesMask = AllSides );
190
191   //@}
192
193   /// initialize global parameters to a give set of values
194   void InitGlobalParameters( Double_t *par );
195
196   /// perform global fit
197   void GlobalFit( Double_t *parameters,Double_t *errors,Double_t *pulls );
198
199   /// print global parameters
200   void PrintGlobalParameters( void ) const;
201
202   /// get error on a given parameter
203   Double_t GetParError( Int_t iPar ) const;
204
205   AliMUONGeometryTransformer* ReAlign( const AliMUONGeometryTransformer * transformer, const double *misAlignments, Bool_t verbose );
206
207   void SetAlignmentResolution( const TClonesArray* misAlignArray, Int_t chId, Double_t chResX, Double_t chResY, Double_t deResX, Double_t deResY );
208
209   private:
210
211   /// Not implemented
212   AliMUONAlignment(const AliMUONAlignment& right);
213
214   /// Not implemented
215   AliMUONAlignment&  operator = (const AliMUONAlignment& right);
216
217   /// Set array of local derivatives
218   void SetLocalDerivative(Int_t index, Double_t value)
219   { fLocalDerivatives[index] = value; }
220
221   /// Set array of global derivatives
222   void SetGlobalDerivative(Int_t index, Double_t value)
223   { fGlobalDerivatives[index] = value; }
224
225   void FillDetElemData( AliMUONVCluster* );
226
227   void FillRecPointData( AliMUONVCluster* );
228
229   void FillTrackParamData( AliMUONTrackParam* );
230
231   Bool_t UnbiasTrackParamData( AliMUONTrackParam* ) const;
232
233   void LocalEquationX( void );
234
235   void LocalEquationY( void );
236
237   TGeoCombiTrans DeltaTransform(const double *detElemMisAlignment) const;
238
239   ///@name utilities
240   //@{
241
242   void AddConstraint(Double_t* parameters, Double_t value );
243
244   Int_t GetChamberId( Int_t iDetElemNumber ) const;
245
246   Bool_t DetElemIsValid( Int_t iDetElemId ) const;
247
248   Int_t GetDetElemNumber( Int_t iDetElemId ) const;
249
250   TString GetParameterMaskString( UInt_t parameterMask ) const;
251
252   TString GetSidesMaskString( UInt_t sidesMask ) const;
253
254   //@}
255
256   /// true when initialized
257   Bool_t fInitialized;
258
259   /// current run id
260   Int_t fRunNumber;
261
262   /// Flag for Magnetic filed On/Off
263   Bool_t fBFieldOn;
264
265   /// "Encouraged" variation for degrees of freedom
266   Double_t fAllowVar[fgNParCh];
267
268   /// Initial value for chi2 cut
269   /** if > 1 Iterations in AliMillepede are turned on */
270   Double_t fStartFac;
271
272   /// Cut on residual for first iteration
273   Double_t fResCutInitial;
274
275   /// Cut on residual for other iterations
276   Double_t fResCut;
277
278   /// Detector independent alignment class
279   AliMillePede2* fMillepede;
280
281   /// running AliMUONVCluster
282   AliMUONVCluster* fCluster;
283
284   /// Number of standard deviations for chi2 cut
285   Int_t fNStdDev;
286
287   /// Cluster (global) position
288   Double_t fClustPos[3];
289
290   /// Track slope at reference point
291   Double_t fTrackSlope0[2];
292
293   /// Track slope at current point
294   Double_t fTrackSlope[2];
295
296   /// Track intersection at reference point
297   Double_t fTrackPos0[3];
298
299   /// Track intersection at current point
300   Double_t fTrackPos[3];
301
302   /// Current measurement (depend on B field On/Off)
303   Double_t fMeas[2];
304
305   /// Estimated resolution on measurement
306   Double_t fSigma[2];
307
308   /// degrees of freedom
309   enum
310   {
311     kFixedParId = -1,
312     kFreeParId = kFixedParId-1,
313     kGroupBaseId = -10
314   };
315
316   /// Array of effective degrees of freedom
317   /// it is used to fix detectors, parameters, etc.
318   Int_t fGlobalParameterStatus[fNGlobal];
319
320   /// Array of global derivatives
321   Double_t fGlobalDerivatives[fNGlobal];
322
323   /// Array of local derivatives
324   Double_t fLocalDerivatives[fNLocal];
325
326   /// current detection element number
327   Int_t fDetElemNumber;
328
329   /// unbias
330   Bool_t fUnbias;
331
332   /// running Track record
333   AliMillePedeRecord fTrackRecord;
334
335   /// Geometry transformation
336   AliMUONGeometryTransformer *fTransform;
337   TGeoCombiTrans fGeoCombiTransInverse;
338
339   ClassDef(AliMUONAlignment, 2)
340
341 };
342
343 #endif