]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/AliMUONAlignment.cxx
Fix for the dumps in the test script
[u/mrichter/AliRoot.git] / MUON / AliMUONAlignment.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /* $Id$ */
17
18 //-----------------------------------------------------------------------------
19 /// \class AliMUONAlignment
20 /// Alignment class for the ALICE DiMuon spectrometer
21 ///
22 /// MUON specific alignment class which interface to AliMillepede.
23 /// For each track ProcessTrack calculates the local and global derivatives
24 /// at each cluster and fill the corresponding local equations. Provide methods
25 /// for fixing or constraining detection elements for best results.
26 ///
27 /// \author Bruce Becker, Javier Castillo
28 //-----------------------------------------------------------------------------
29
30 #include "AliMUONAlignment.h"
31 #include "AliMUONTrack.h"
32 #include "AliMUONTrackParam.h"
33 #include "AliMUONVCluster.h"
34 #include "AliMUONGeometryTransformer.h"
35 #include "AliMUONGeometryModuleTransformer.h"
36 #include "AliMUONGeometryDetElement.h"
37 #include "AliMUONGeometryBuilder.h"
38 #include "AliMillepede.h"
39
40 #include "AliMpExMap.h"
41 #include "AliMpExMapIterator.h"
42
43 #include "AliAlignObjMatrix.h"
44 #include "AliLog.h"
45
46 #include "TMath.h"
47 #include "TMatrixDSym.h"
48 #include "TClonesArray.h"
49
50 /// \cond CLASSIMP
51 ClassImp(AliMUONAlignment)
52 /// \endcond
53
54 //_____________________________________________________________________
55 // static variables
56 Int_t AliMUONAlignment::fgNDetElem = 4*2+4*2+18*2+26*2+26*2;
57 Int_t AliMUONAlignment::fgNDetElemCh[10] = {4,4,4,4,18,18,26,26,26,26};
58 Int_t AliMUONAlignment::fgSNDetElemCh[10] = {4,8,12,16,34,52,78,104,130,156};
59 Int_t AliMUONAlignment::fgNParCh = 4;
60 Int_t AliMUONAlignment::fgNTrkMod = 16;
61 Int_t AliMUONAlignment::fgNCh = 10;
62 Int_t AliMUONAlignment::fgNSt = 5;
63
64 //_____________________________________________________________________
65 AliMUONAlignment::AliMUONAlignment()
66   : TObject(),
67     fBFieldOn(kTRUE),
68     fStartFac(256.),
69     fResCutInitial(100.),
70     fResCut(100.),
71     fMillepede(0),
72     fTrack(0),
73     fCluster(0),
74     fTrackParam(0),
75     fNGlobal(fgNDetElem*fgNParCh),
76     fNLocal(4),
77     fNStdDev(3),
78     fDetElemId(0),
79     fDetElemNumber(0),
80     fPhi(0.),
81     fCosPhi(1.),
82     fSinPhi(0.),
83     fTrackRecord(),
84     fTransform(0)
85 {
86
87   fSigma[0] = 1.5e-1;
88   fSigma[1] = 1.0e-2;
89
90   AliInfo(Form("fSigma[0]: %f\t fSigma[1]: %f",fSigma[0],fSigma[1]));
91
92   fDoF[0] = kTRUE;  fDoF[1] = kTRUE;  fDoF[2] = kTRUE;  fDoF[3] = kTRUE;
93   fAllowVar[0] = 0.05;  fAllowVar[1] = 0.05;  fAllowVar[2] = 0.001;  fAllowVar[3] = 0.5;
94
95   AliInfo(Form("fAllowVar[0]: %f\t fAllowVar[1]: %f\t fPhi: %f\t fgNDetElem: %i\t fNGlobal: %i\t fNLocal: %i",fAllowVar[0],fAllowVar[1],fPhi,fgNDetElem,fNGlobal,fNLocal));
96
97   fMillepede = new AliMillepede();
98
99   Init(fNGlobal, fNLocal, fNStdDev);
100
101   ResetLocalEquation();
102   AliInfo("Parameters initialized to zero");
103
104 }
105
106 //_____________________________________________________________________
107 AliMUONAlignment::~AliMUONAlignment()
108 {}
109
110 //_____________________________________________________________________
111 void AliMUONAlignment::Init(
112   Int_t nGlobal,  /* number of global paramers */
113   Int_t nLocal,   /* number of local parameters */
114   Int_t nStdDev   /* std dev cut */ )
115 {
116   /// Initialization of AliMillepede. Fix parameters, define constraints ...
117   fMillepede->InitMille(nGlobal,nLocal,nStdDev,fResCut,fResCutInitial);
118
119   //  Bool_t bStOnOff[5] = {kTRUE,kTRUE,kTRUE,kTRUE,kTRUE};
120   //  Bool_t bChOnOff[10] = {kTRUE,kTRUE,kTRUE,kTRUE,kTRUE,kTRUE,kTRUE,kTRUE,kTRUE,kTRUE};
121   //  Bool_t bSpecLROnOff[2] = {kTRUE,kTRUE};
122
123   //  AllowVariations(bChOnOff);
124
125   // Fix parameters or add constraints here
126   //   for (Int_t iSt=0; iSt<5; iSt++)
127   //   { if (!bStOnOff[iSt]) FixStation(iSt+1); }
128
129   //   for (Int_t iCh=0; iCh<10; iCh++)
130   //   { if (!bChOnOff[iCh]) FixChamber(iCh+1); }
131
132 //   FixHalfSpectrometer(bChOnOff,bSpecLROnOff);
133
134   ResetConstraints();
135
136   // Define global constrains to be applied
137   // X, Y, P, XvsZ, YvsZ, PvsZ, XvsY, YvsY, PvsY
138   Bool_t bVarXYT[9] = {kTRUE,kTRUE,kTRUE,kTRUE,kTRUE,kTRUE,kTRUE,kTRUE,kTRUE};
139   Bool_t bDetTLBR[4] = {kFALSE,kTRUE,kFALSE,kTRUE};
140   //  AddConstraints(bStOnOff,bVarXYT,bDetTLBR,bSpecLROnOff);
141
142   // Other possible way to add constrains
143   bVarXYT[0] = kFALSE; bVarXYT[1] = kFALSE; bVarXYT[2] = kTRUE;
144   bDetTLBR[0] = kFALSE; bDetTLBR[1] = kTRUE; bDetTLBR[2] = kFALSE; bDetTLBR[3] = kFALSE;
145   //  AddConstraints(bStOnOff,bVarXYT,bDetTLBR);
146
147   bVarXYT[0] = kTRUE; bVarXYT[1] = kTRUE; bVarXYT[2] = kFALSE;
148   //  AddConstraints(bStOnOff,bVarXYT);
149
150   // Set iterations
151   if (fStartFac>1) fMillepede->SetIterations(fStartFac);
152 }
153
154 //_____________________________________________________________________
155 void AliMUONAlignment::FixStation(Int_t iSt)
156 {
157   /// Fix all detection elements of station iSt
158   Int_t iDetElemFirst = (iSt>1) ? fgSNDetElemCh[2*(iSt-1)-1] : 0;
159   Int_t iDetElemLast = fgSNDetElemCh[2*(iSt)-1];
160   for (Int_t i = iDetElemFirst; i < iDetElemLast; i++)
161   {
162     FixParameter(i*fgNParCh+0, 0.0);
163     FixParameter(i*fgNParCh+1, 0.0);
164     FixParameter(i*fgNParCh+2, 0.0);
165     FixParameter(i*fgNParCh+3, 0.0);
166   }
167
168 }
169
170 //_____________________________________________________________________
171 void AliMUONAlignment::FixChamber(Int_t iCh)
172 {
173   /// Fix all detection elements of chamber iCh
174   Int_t iDetElemFirst = (iCh>1) ? fgSNDetElemCh[iCh-2] : 0;
175   Int_t iDetElemLast = fgSNDetElemCh[iCh-1];
176   for (Int_t i = iDetElemFirst; i < iDetElemLast; i++)
177   {
178     FixParameter(i*fgNParCh+0, 0.0);
179     FixParameter(i*fgNParCh+1, 0.0);
180     FixParameter(i*fgNParCh+2, 0.0);
181     FixParameter(i*fgNParCh+3, 0.0);
182   }
183 }
184
185 //_____________________________________________________________________
186 void AliMUONAlignment::FixDetElem(Int_t iDetElemId, TString sVarXYT)
187 {
188
189   /// Fix a given detection element
190   Int_t iDetElemNumber = iDetElemId%100;
191   for (int iCh=0; iCh<iDetElemId/100-1; iCh++)
192   {
193     iDetElemNumber += fgNDetElemCh[iCh];
194   }
195
196   if (sVarXYT.Contains("X"))
197   {
198     // X constraint
199     FixParameter(iDetElemNumber*fgNParCh+0, 0.0);
200   }
201
202   if (sVarXYT.Contains("Y"))
203   {
204     // Y constraint
205     FixParameter(iDetElemNumber*fgNParCh+1, 0.0);
206   }
207
208   if (sVarXYT.Contains("T"))
209   {
210     // T constraint
211     FixParameter(iDetElemNumber*fgNParCh+2, 0.0);
212   }
213
214   if (sVarXYT.Contains("Z"))
215   {
216     // T constraint
217     FixParameter(iDetElemNumber*fgNParCh+3, 0.0);
218   }
219
220 }
221
222 //_____________________________________________________________________
223 void AliMUONAlignment::FixHalfSpectrometer(const Bool_t *lChOnOff, const Bool_t *lSpecLROnOff)
224 {
225
226   /// Fix left or right detector
227   for (Int_t i = 0; i < fgNDetElem; i++)
228   {
229
230     Int_t iCh=0;
231     for (iCh=1; iCh<=fgNCh; iCh++)
232     { if (i<fgSNDetElemCh[iCh-1]) break; }
233
234     if (lChOnOff[iCh-1])
235     {
236
237       Int_t lDetElemNumber = (iCh==1) ? i : i-fgSNDetElemCh[iCh-2];
238       if (iCh>=1 && iCh<=4)
239       {
240
241         if ((lDetElemNumber==1 || lDetElemNumber==2) && !lSpecLROnOff[0])
242         {
243           // From track crossings
244           FixParameter(i*fgNParCh+0, 0.0);
245           FixParameter(i*fgNParCh+1, 0.0);
246           FixParameter(i*fgNParCh+2, 0.0);
247           FixParameter(i*fgNParCh+3, 0.0);
248         }
249
250         if ((lDetElemNumber==0 || lDetElemNumber==3) && !lSpecLROnOff[1])
251         {
252           // From track crossings
253           FixParameter(i*fgNParCh+0, 0.0);
254           FixParameter(i*fgNParCh+1, 0.0);
255           FixParameter(i*fgNParCh+2, 0.0);
256           FixParameter(i*fgNParCh+3, 0.0);
257         }
258
259       }
260
261       if (iCh>=5 && iCh<=6)
262       {
263
264         if ((lDetElemNumber>=5&&lDetElemNumber<=13) && !lSpecLROnOff[0])
265         {
266           FixParameter(i*fgNParCh+0, 0.0);
267           FixParameter(i*fgNParCh+1, 0.0);
268           FixParameter(i*fgNParCh+2, 0.0);
269           FixParameter(i*fgNParCh+3, 0.0);
270         }
271
272         if (((lDetElemNumber>=0&&lDetElemNumber<=4) ||
273           (lDetElemNumber>=14&&lDetElemNumber<=17)) && !lSpecLROnOff[1])
274         {
275
276           FixParameter(i*fgNParCh+0, 0.0);
277           FixParameter(i*fgNParCh+1, 0.0);
278           FixParameter(i*fgNParCh+2, 0.0);
279           FixParameter(i*fgNParCh+3, 0.0);
280         }
281
282       }
283
284       if (iCh>=7 && iCh<=10)
285       {
286
287         if ((lDetElemNumber>=7&&lDetElemNumber<=19) && !lSpecLROnOff[0])
288         {
289           FixParameter(i*fgNParCh+0, 0.0);
290           FixParameter(i*fgNParCh+1, 0.0);
291           FixParameter(i*fgNParCh+2, 0.0);
292           FixParameter(i*fgNParCh+3, 0.0);
293         }
294
295         if (((lDetElemNumber>=0&&lDetElemNumber<=6) ||
296           (lDetElemNumber>=20&&lDetElemNumber<=25)) && !lSpecLROnOff[1])
297         {
298           FixParameter(i*fgNParCh+0, 0.0);
299           FixParameter(i*fgNParCh+1, 0.0);
300           FixParameter(i*fgNParCh+2, 0.0);
301           FixParameter(i*fgNParCh+3, 0.0);
302         }
303
304       }
305
306     }
307
308   }
309
310 }
311
312 //______________________________________________________________________
313 void AliMUONAlignment::SetNonLinear(const Bool_t *lChOnOff, const Bool_t *lVarXYT)
314 {
315
316   /// Set non linear parameter flag selected chambers and degrees of freedom
317   for (Int_t i = 0; i < fgNDetElem; i++)
318   {
319
320     Int_t iCh=0;
321     for (iCh=1; iCh<=fgNCh; iCh++)
322     { if (i<fgSNDetElemCh[iCh-1]) break; }
323
324     if (lChOnOff[iCh-1])
325     {
326
327       if (lVarXYT[0])
328       {
329         // X constraint
330         SetNonLinear(i*fgNParCh+0);
331       }
332
333       if (lVarXYT[1])
334       {
335         // Y constraint
336         SetNonLinear(i*fgNParCh+1);
337       }
338
339       if (lVarXYT[2])
340       {
341         // T constraint
342         SetNonLinear(i*fgNParCh+2);
343       }
344
345       if (lVarXYT[3])
346       {
347         // Z constraint
348         SetNonLinear(i*fgNParCh+3);
349       }
350
351     }
352
353   }
354
355 }
356
357 //______________________________________________________________________
358 void AliMUONAlignment::AddConstraints(const Bool_t *lChOnOff, const Bool_t *lVarXYT)
359 {
360
361   /// Add constraint equations for selected chambers and degrees of freedom
362   for (Int_t i = 0; i < fgNDetElem; i++)
363   {
364
365     Int_t iCh=0;
366     for (iCh=1; iCh<=fgNCh; iCh++)
367     {
368       if (i<fgSNDetElemCh[iCh-1]) break;
369     }
370
371     if (lChOnOff[iCh-1])
372     {
373       if (lVarXYT[0])
374       {
375         // X constraint
376         fConstraintX[i*fgNParCh+0]=1.0;
377       }
378
379       if (lVarXYT[1])
380       {
381         // Y constraint
382         fConstraintY[i*fgNParCh+1]=1.0;
383       }
384
385       if (lVarXYT[2])
386       {
387         // T constraint
388         fConstraintP[i*fgNParCh+2]=1.0;
389       }
390       //       if (lVarXYT[3]) { // Z constraint
391       //        fConstraintP[i*fgNParCh+3]=1.0;
392       //       }
393
394     }
395   }
396
397   if (lVarXYT[0])
398   {
399     // X constraint
400     AddConstraint(fConstraintX,0.0);
401   }
402
403   if (lVarXYT[1])
404   {
405     // Y constraint
406     AddConstraint(fConstraintY,0.0);
407   }
408
409   if (lVarXYT[2])
410   {
411     // T constraint
412     AddConstraint(fConstraintP,0.0);
413   }
414
415 //   if (lVarXYT[3]) { // Z constraint
416 //     AddConstraint(fConstraintP,0.0);
417 //   }
418 }
419
420 //______________________________________________________________________
421 void AliMUONAlignment::AddConstraints(const Bool_t *lChOnOff, const Bool_t *lVarXYT, const Bool_t *lDetTLBR, const Bool_t *lSpecLROnOff)
422 {
423   /// Add constraint equations for selected chambers, degrees of freedom and detector half
424   Double_t lDetElemLocX = 0.;
425   Double_t lDetElemLocY = 0.;
426   Double_t lDetElemLocZ = 0.;
427   Double_t lDetElemGloX = 0.;
428   Double_t lDetElemGloY = 0.;
429   Double_t lDetElemGloZ = 0.;
430   Double_t lMeanY = 0.;
431   Double_t lSigmaY = 0.;
432   Double_t lMeanZ = 0.;
433   Double_t lSigmaZ = 0.;
434   Int_t lNDetElem = 0;
435   for (Int_t i = 0; i < fgNDetElem; i++)
436   {
437
438     Int_t iCh=0;
439     for (iCh=1; iCh<=fgNCh; iCh++){
440       if (i<fgSNDetElemCh[iCh-1]) break;
441     }
442     if (lChOnOff[iCh-1]){
443       Int_t lDetElemNumber = (iCh==1) ? i : i-fgSNDetElemCh[iCh-2];
444       Int_t lDetElemId = iCh*100+lDetElemNumber;
445       fTransform->Local2Global(lDetElemId,lDetElemLocX,lDetElemLocY,lDetElemLocZ,
446              lDetElemGloX,lDetElemGloY,lDetElemGloZ);
447       if (iCh>=1 && iCh<=4){
448   if ((lDetElemNumber==1 || lDetElemNumber==2) && lSpecLROnOff[0]){ // From track crossings
449     lMeanY += lDetElemGloY;
450     lSigmaY += lDetElemGloY*lDetElemGloY;
451     lMeanZ += lDetElemGloZ;
452     lSigmaZ += lDetElemGloZ*lDetElemGloZ;
453     lNDetElem++;
454   }
455   if ((lDetElemNumber==0 || lDetElemNumber==3) && lSpecLROnOff[1]){ // From track crossings
456     lMeanY += lDetElemGloY;
457     lSigmaY += lDetElemGloY*lDetElemGloY;
458     lMeanZ += lDetElemGloZ;
459     lSigmaZ += lDetElemGloZ*lDetElemGloZ;
460     lNDetElem++;
461   }
462       }
463       if (iCh>=5 && iCh<=6){
464   if ((lDetElemNumber>=5&&lDetElemNumber<=13) && lSpecLROnOff[0]){
465     lMeanY += lDetElemGloY;
466     lSigmaY += lDetElemGloY*lDetElemGloY;
467     lMeanZ += lDetElemGloZ;
468     lSigmaZ += lDetElemGloZ*lDetElemGloZ;
469     lNDetElem++;
470   }
471   if (((lDetElemNumber>=0&&lDetElemNumber<=4) ||
472        (lDetElemNumber>=14&&lDetElemNumber<=17)) && lSpecLROnOff[1]){
473     lMeanY += lDetElemGloY;
474     lSigmaY += lDetElemGloY*lDetElemGloY;
475     lMeanZ += lDetElemGloZ;
476     lSigmaZ += lDetElemGloZ*lDetElemGloZ;
477     lNDetElem++;
478   }
479       }
480       if (iCh>=7 && iCh<=10){
481   if ((lDetElemNumber>=7&&lDetElemNumber<=19) && lSpecLROnOff[0]){
482     lMeanY += lDetElemGloY;
483     lSigmaY += lDetElemGloY*lDetElemGloY;
484     lMeanZ += lDetElemGloZ;
485     lSigmaZ += lDetElemGloZ*lDetElemGloZ;
486     lNDetElem++;
487   }
488   if (((lDetElemNumber>=0&&lDetElemNumber<=6) ||
489        (lDetElemNumber>=20&&lDetElemNumber<=25)) && lSpecLROnOff[1]){
490     lMeanY += lDetElemGloY;
491     lSigmaY += lDetElemGloY*lDetElemGloY;
492     lMeanZ += lDetElemGloZ;
493     lSigmaZ += lDetElemGloZ*lDetElemGloZ;
494     lNDetElem++;
495   }
496       }
497     }
498   }
499   lMeanY /= lNDetElem;
500   lSigmaY /= lNDetElem;
501   lSigmaY = TMath::Sqrt(lSigmaY-lMeanY*lMeanY);
502   lMeanZ /= lNDetElem;
503   lSigmaZ /= lNDetElem;
504   lSigmaZ = TMath::Sqrt(lSigmaZ-lMeanZ*lMeanZ);
505   AliInfo(Form("Used %i DetElem, MeanZ= %f , SigmaZ= %f", lNDetElem,lMeanZ,lSigmaZ));
506
507   for (Int_t i = 0; i < fgNDetElem; i++){
508     Int_t iCh=0;
509     for (iCh=1; iCh<=fgNCh; iCh++){
510       if (i<fgSNDetElemCh[iCh-1]) break;
511     }
512     if (lChOnOff[iCh-1]){
513       Int_t lDetElemNumber = (iCh==1) ? i : i-fgSNDetElemCh[iCh-2];
514       Int_t lDetElemId = iCh*100+lDetElemNumber;
515       fTransform->Local2Global(lDetElemId,lDetElemLocX,lDetElemLocY,lDetElemLocZ,
516              lDetElemGloX,lDetElemGloY,lDetElemGloZ);
517       if (lVarXYT[0]) { // X constraint
518   if (lDetTLBR[0]) ConstrainT(i,iCh,fConstraintXT,0); // Top half
519   if (lDetTLBR[1]) ConstrainL(i,iCh,fConstraintXL,0); // Left half
520   if (lDetTLBR[2]) ConstrainB(i,iCh,fConstraintXB,0); // Bottom half
521   if (lDetTLBR[3]) ConstrainR(i,iCh,fConstraintXR,0); // Right half
522       }
523       if (lVarXYT[1]) { // Y constraint
524   if (lDetTLBR[0]) ConstrainT(i,iCh,fConstraintYT,1); // Top half
525   if (lDetTLBR[1]) ConstrainL(i,iCh,fConstraintYL,1); // Left half
526   if (lDetTLBR[2]) ConstrainB(i,iCh,fConstraintYB,1); // Bottom half
527   if (lDetTLBR[3]) ConstrainR(i,iCh,fConstraintYR,1); // Right half
528       }
529       if (lVarXYT[2]) { // P constraint
530   if (lDetTLBR[0]) ConstrainT(i,iCh,fConstraintPT,2); // Top half
531   if (lDetTLBR[1]) ConstrainL(i,iCh,fConstraintPL,2); // Left half
532   if (lDetTLBR[2]) ConstrainB(i,iCh,fConstraintPB,2); // Bottom half
533   if (lDetTLBR[3]) ConstrainR(i,iCh,fConstraintPR,2); // Right half
534       }
535       if (lVarXYT[3]) { // X-Z shearing
536   if (lDetTLBR[0]) ConstrainT(i,iCh,fConstraintXZT,0,(lDetElemGloZ-lMeanZ)/lSigmaZ); // Top half
537   if (lDetTLBR[1]) ConstrainL(i,iCh,fConstraintXZL,0,(lDetElemGloZ-lMeanZ)/lSigmaZ); // Left half
538   if (lDetTLBR[2]) ConstrainB(i,iCh,fConstraintXZB,0,(lDetElemGloZ-lMeanZ)/lSigmaZ); // Bottom half
539   if (lDetTLBR[3]) ConstrainR(i,iCh,fConstraintXZR,0,(lDetElemGloZ-lMeanZ)/lSigmaZ); // Right half
540       }
541       if (lVarXYT[4]) { // Y-Z shearing
542   if (lDetTLBR[0]) ConstrainT(i,iCh,fConstraintYZT,1,(lDetElemGloZ-lMeanZ)/lSigmaZ); // Top half
543   if (lDetTLBR[1]) ConstrainL(i,iCh,fConstraintYZL,1,(lDetElemGloZ-lMeanZ)/lSigmaZ); // Left half
544   if (lDetTLBR[2]) ConstrainB(i,iCh,fConstraintYZB,1,(lDetElemGloZ-lMeanZ)/lSigmaZ); // Bottom half
545   if (lDetTLBR[3]) ConstrainR(i,iCh,fConstraintYZR,1,(lDetElemGloZ-lMeanZ)/lSigmaZ); // Right half
546       }
547       if (lVarXYT[5]) { // P-Z rotation
548   if (lDetTLBR[0]) ConstrainT(i,iCh,fConstraintPZT,2,(lDetElemGloZ-lMeanZ)/lSigmaZ); // Top half
549   if (lDetTLBR[1]) ConstrainL(i,iCh,fConstraintPZL,2,(lDetElemGloZ-lMeanZ)/lSigmaZ); // Left half
550   if (lDetTLBR[2]) ConstrainB(i,iCh,fConstraintPZB,2,(lDetElemGloZ-lMeanZ)/lSigmaZ); // Bottom half
551   if (lDetTLBR[3]) ConstrainR(i,iCh,fConstraintPZR,2,(lDetElemGloZ-lMeanZ)/lSigmaZ); // Right half
552       }
553       if (lVarXYT[6]) { // X-Y shearing
554   if (lDetTLBR[0]) ConstrainT(i,iCh,fConstraintXYT,0,(lDetElemGloY-lMeanY)/lSigmaY); // Top half
555   if (lDetTLBR[1]) ConstrainL(i,iCh,fConstraintXYL,0,(lDetElemGloY-lMeanY)/lSigmaY); // Left half
556   if (lDetTLBR[2]) ConstrainB(i,iCh,fConstraintXYB,0,(lDetElemGloY-lMeanY)/lSigmaY); // Bottom half
557   if (lDetTLBR[3]) ConstrainR(i,iCh,fConstraintXYR,0,(lDetElemGloY-lMeanY)/lSigmaY); // Right half
558       }
559       if (lVarXYT[7]) { // Y-Y scaling
560   if (lDetTLBR[0]) ConstrainT(i,iCh,fConstraintYYT,1,(lDetElemGloY-lMeanY)/lSigmaY); // Top half
561   if (lDetTLBR[1]) ConstrainL(i,iCh,fConstraintYYL,1,(lDetElemGloY-lMeanY)/lSigmaY); // Left half
562   if (lDetTLBR[2]) ConstrainB(i,iCh,fConstraintYYB,1,(lDetElemGloY-lMeanY)/lSigmaY); // Bottom half
563   if (lDetTLBR[3]) ConstrainR(i,iCh,fConstraintYYR,1,(lDetElemGloY-lMeanY)/lSigmaY); // Right half
564       }
565       if (lVarXYT[8]) { // P-Y rotation
566   if (lDetTLBR[0]) ConstrainT(i,iCh,fConstraintPYT,2,(lDetElemGloY-lMeanY)/lSigmaY); // Top half
567   if (lDetTLBR[1]) ConstrainL(i,iCh,fConstraintPYL,2,(lDetElemGloY-lMeanY)/lSigmaY); // Left half
568   if (lDetTLBR[2]) ConstrainB(i,iCh,fConstraintPYB,2,(lDetElemGloY-lMeanY)/lSigmaY); // Bottom half
569   if (lDetTLBR[3]) ConstrainR(i,iCh,fConstraintPYR,2,(lDetElemGloY-lMeanY)/lSigmaY); // Right half
570       }
571     }
572   }
573   if (lVarXYT[0]) { // X constraint
574     if (lDetTLBR[0]) AddConstraint(fConstraintXT,0.0); // Top half
575     if (lDetTLBR[1]) AddConstraint(fConstraintXL,0.0); // Left half
576     if (lDetTLBR[2]) AddConstraint(fConstraintXB,0.0); // Bottom half
577     if (lDetTLBR[3]) AddConstraint(fConstraintXR,0.0); // Right half
578   }
579   if (lVarXYT[1]) { // Y constraint
580     if (lDetTLBR[0]) AddConstraint(fConstraintYT,0.0); // Top half
581     if (lDetTLBR[1]) AddConstraint(fConstraintYL,0.0); // Left half
582     if (lDetTLBR[2]) AddConstraint(fConstraintYB,0.0); // Bottom half
583     if (lDetTLBR[3]) AddConstraint(fConstraintYR,0.0); // Right half
584   }
585   if (lVarXYT[2]) { // T constraint
586     if (lDetTLBR[0]) AddConstraint(fConstraintPT,0.0); // Top half
587     if (lDetTLBR[1]) AddConstraint(fConstraintPL,0.0); // Left half
588     if (lDetTLBR[2]) AddConstraint(fConstraintPB,0.0); // Bottom half
589     if (lDetTLBR[3]) AddConstraint(fConstraintPR,0.0); // Right half
590   }
591   if (lVarXYT[3]) { // X-Z constraint
592     if (lDetTLBR[0]) AddConstraint(fConstraintXZT,0.0); // Top half
593     if (lDetTLBR[1]) AddConstraint(fConstraintXZL,0.0); // Left half
594     if (lDetTLBR[2]) AddConstraint(fConstraintXZB,0.0); // Bottom half
595     if (lDetTLBR[3]) AddConstraint(fConstraintXZR,0.0); // Right half
596   }
597   if (lVarXYT[4]) { // Y-Z constraint
598     if (lDetTLBR[0]) AddConstraint(fConstraintYZT,0.0); // Top half
599     if (lDetTLBR[1]) AddConstraint(fConstraintYZL,0.0); // Left half
600     if (lDetTLBR[2]) AddConstraint(fConstraintYZB,0.0); // Bottom half
601     if (lDetTLBR[3]) AddConstraint(fConstraintYZR,0.0); // Right half
602   }
603   if (lVarXYT[5]) { // P-Z constraint
604     if (lDetTLBR[0]) AddConstraint(fConstraintPZT,0.0); // Top half
605     if (lDetTLBR[1]) AddConstraint(fConstraintPZL,0.0); // Left half
606     if (lDetTLBR[2]) AddConstraint(fConstraintPZB,0.0); // Bottom half
607     if (lDetTLBR[3]) AddConstraint(fConstraintPZR,0.0); // Right half
608   }
609   if (lVarXYT[6]) { // X-Y constraint
610     if (lDetTLBR[0]) AddConstraint(fConstraintXYT,0.0); // Top half
611     if (lDetTLBR[1]) AddConstraint(fConstraintXYL,0.0); // Left half
612     if (lDetTLBR[2]) AddConstraint(fConstraintXYB,0.0); // Bottom half
613     if (lDetTLBR[3]) AddConstraint(fConstraintXYR,0.0); // Right half
614   }
615   if (lVarXYT[7]) { // Y-Y constraint
616     if (lDetTLBR[0]) AddConstraint(fConstraintYYT,0.0); // Top half
617     if (lDetTLBR[1]) AddConstraint(fConstraintYYL,0.0); // Left half
618     if (lDetTLBR[2]) AddConstraint(fConstraintYYB,0.0); // Bottom half
619     if (lDetTLBR[3]) AddConstraint(fConstraintYYR,0.0); // Right half
620   }
621   if (lVarXYT[8]) { // P-Y constraint
622     if (lDetTLBR[0]) AddConstraint(fConstraintPYT,0.0); // Top half
623     if (lDetTLBR[1]) AddConstraint(fConstraintPYL,0.0); // Left half
624     if (lDetTLBR[2]) AddConstraint(fConstraintPYB,0.0); // Bottom half
625     if (lDetTLBR[3]) AddConstraint(fConstraintPYR,0.0); // Right half
626   }
627 }
628
629 void AliMUONAlignment::ConstrainT(Int_t lDetElem, Int_t lCh, Double_t *lConstraintT, Int_t iVar, Double_t /*lWeight*/) const{
630   /// Set constrain equation for top half of spectrometer
631   Int_t lDetElemNumber = (lCh==1) ? lDetElem : lDetElem-fgSNDetElemCh[lCh-2];
632   if (lCh>=1 && lCh<=4){
633     if (lDetElemNumber==0 || lDetElemNumber==1){ // From track crossings
634       lConstraintT[lDetElem*fgNParCh+iVar]=1.0;
635     }
636   }
637   if (lCh>=5 && lCh<=6){
638     if (lDetElemNumber>=0&&lDetElemNumber<=9){
639       lConstraintT[lDetElem*fgNParCh+iVar]=1.0;
640     }
641   }
642   if (lCh>=7 && lCh<=10){
643     if (lDetElemNumber>=0&&lDetElemNumber<=13){
644       lConstraintT[lDetElem*fgNParCh+iVar]=1.0;
645     }
646   }
647 }
648
649 //______________________________________________________________________
650 void AliMUONAlignment::ConstrainL(Int_t lDetElem, Int_t lCh, Double_t *lConstraintL, Int_t iVar, Double_t lWeight) const{
651   /// Set constrain equation for left half of spectrometer
652   Int_t lDetElemNumber = (lCh==1) ? lDetElem : lDetElem-fgSNDetElemCh[lCh-2];
653   if (lCh>=1 && lCh<=4){
654     if (lDetElemNumber==1 || lDetElemNumber==2){ // From track crossings
655       lConstraintL[lDetElem*fgNParCh+iVar]=lWeight;
656     }
657   }
658   if (lCh>=5 && lCh<=6){
659     if (lDetElemNumber>=5&&lDetElemNumber<=13){
660       lConstraintL[lDetElem*fgNParCh+iVar]=lWeight;
661     }
662   }
663   if (lCh>=7 && lCh<=10){
664     if (lDetElemNumber>=7&&lDetElemNumber<=19){
665       lConstraintL[lDetElem*fgNParCh+iVar]=lWeight;
666     }
667   }
668 }
669
670 //______________________________________________________________________
671 void AliMUONAlignment::ConstrainB(Int_t lDetElem, Int_t lCh, Double_t *lConstraintB, Int_t iVar, Double_t /*lWeight*/) const{
672   /// Set constrain equation for bottom half of spectrometer
673   Int_t lDetElemNumber = (lCh==1) ? lDetElem : lDetElem-fgSNDetElemCh[lCh-2];
674   if (lCh>=1 && lCh<=4){
675     if (lDetElemNumber==2 || lDetElemNumber==3){ // From track crossings
676       lConstraintB[lDetElem*fgNParCh+iVar]=1.0;
677     }
678   }
679   if (lCh>=5 && lCh<=6){
680     if ((lDetElemNumber>=9&&lDetElemNumber<=17) ||
681   (lDetElemNumber==0)){
682       lConstraintB[lDetElem*fgNParCh+iVar]=1.0;
683     }
684   }
685   if (lCh>=7 && lCh<=10){
686     if ((lDetElemNumber>=13&&lDetElemNumber<=25) ||
687   (lDetElemNumber==0)){
688       lConstraintB[lDetElem*fgNParCh+iVar]=1.0;
689     }
690   }
691 }
692
693 //______________________________________________________________________
694 void AliMUONAlignment::ConstrainR(Int_t lDetElem, Int_t lCh, Double_t *lConstraintR, Int_t iVar, Double_t lWeight) const{
695   /// Set constrain equation for right half of spectrometer
696   Int_t lDetElemNumber = (lCh==1) ? lDetElem : lDetElem-fgSNDetElemCh[lCh-2];
697   if (lCh>=1 && lCh<=4){
698     if (lDetElemNumber==0 || lDetElemNumber==3){ // From track crossings
699       lConstraintR[lDetElem*fgNParCh+iVar]=lWeight;
700     }
701   }
702   if (lCh>=5 && lCh<=6){
703     if ((lDetElemNumber>=0&&lDetElemNumber<=4) ||
704   (lDetElemNumber>=14&&lDetElemNumber<=17)){
705       lConstraintR[lDetElem*fgNParCh+iVar]=lWeight;
706     }
707   }
708   if (lCh>=7 && lCh<=10){
709     if ((lDetElemNumber>=0&&lDetElemNumber<=6) ||
710   (lDetElemNumber>=20&&lDetElemNumber<=25)){
711       lConstraintR[lDetElem*fgNParCh+iVar]=lWeight;
712     }
713   }
714 }
715
716 //______________________________________________________________________
717 void AliMUONAlignment::ResetConstraints(){
718   /// Reset all constraint equations
719   for (Int_t i = 0; i < fgNDetElem; i++){
720     fConstraintX[i*fgNParCh+0]=0.0;
721     fConstraintX[i*fgNParCh+1]=0.0;
722     fConstraintX[i*fgNParCh+2]=0.0;
723     fConstraintY[i*fgNParCh+0]=0.0;
724     fConstraintY[i*fgNParCh+1]=0.0;
725     fConstraintY[i*fgNParCh+2]=0.0;
726     fConstraintP[i*fgNParCh+0]=0.0;
727     fConstraintP[i*fgNParCh+1]=0.0;
728     fConstraintP[i*fgNParCh+2]=0.0;
729     fConstraintXT[i*fgNParCh+0]=0.0;
730     fConstraintXT[i*fgNParCh+1]=0.0;
731     fConstraintXT[i*fgNParCh+2]=0.0;
732     fConstraintYT[i*fgNParCh+0]=0.0;
733     fConstraintYT[i*fgNParCh+1]=0.0;
734     fConstraintYT[i*fgNParCh+2]=0.0;
735     fConstraintPT[i*fgNParCh+0]=0.0;
736     fConstraintPT[i*fgNParCh+1]=0.0;
737     fConstraintPT[i*fgNParCh+2]=0.0;
738     fConstraintXZT[i*fgNParCh+0]=0.0;
739     fConstraintXZT[i*fgNParCh+1]=0.0;
740     fConstraintXZT[i*fgNParCh+2]=0.0;
741     fConstraintYZT[i*fgNParCh+0]=0.0;
742     fConstraintYZT[i*fgNParCh+1]=0.0;
743     fConstraintYZT[i*fgNParCh+2]=0.0;
744     fConstraintPZT[i*fgNParCh+0]=0.0;
745     fConstraintPZT[i*fgNParCh+1]=0.0;
746     fConstraintPZT[i*fgNParCh+2]=0.0;
747     fConstraintXYT[i*fgNParCh+0]=0.0;
748     fConstraintXYT[i*fgNParCh+1]=0.0;
749     fConstraintXYT[i*fgNParCh+2]=0.0;
750     fConstraintYYT[i*fgNParCh+0]=0.0;
751     fConstraintYYT[i*fgNParCh+1]=0.0;
752     fConstraintYYT[i*fgNParCh+2]=0.0;
753     fConstraintPYT[i*fgNParCh+0]=0.0;
754     fConstraintPYT[i*fgNParCh+1]=0.0;
755     fConstraintPYT[i*fgNParCh+2]=0.0;
756     fConstraintXL[i*fgNParCh+0]=0.0;
757     fConstraintXL[i*fgNParCh+1]=0.0;
758     fConstraintXL[i*fgNParCh+2]=0.0;
759     fConstraintYL[i*fgNParCh+0]=0.0;
760     fConstraintYL[i*fgNParCh+1]=0.0;
761     fConstraintYL[i*fgNParCh+2]=0.0;
762     fConstraintPL[i*fgNParCh+0]=0.0;
763     fConstraintPL[i*fgNParCh+1]=0.0;
764     fConstraintPL[i*fgNParCh+2]=0.0;
765     fConstraintXZL[i*fgNParCh+0]=0.0;
766     fConstraintXZL[i*fgNParCh+1]=0.0;
767     fConstraintXZL[i*fgNParCh+2]=0.0;
768     fConstraintYZL[i*fgNParCh+0]=0.0;
769     fConstraintYZL[i*fgNParCh+1]=0.0;
770     fConstraintYZL[i*fgNParCh+2]=0.0;
771     fConstraintPZL[i*fgNParCh+0]=0.0;
772     fConstraintPZL[i*fgNParCh+1]=0.0;
773     fConstraintPZL[i*fgNParCh+2]=0.0;
774     fConstraintXYL[i*fgNParCh+0]=0.0;
775     fConstraintXYL[i*fgNParCh+1]=0.0;
776     fConstraintXYL[i*fgNParCh+2]=0.0;
777     fConstraintYYL[i*fgNParCh+0]=0.0;
778     fConstraintYYL[i*fgNParCh+1]=0.0;
779     fConstraintYYL[i*fgNParCh+2]=0.0;
780     fConstraintPYL[i*fgNParCh+0]=0.0;
781     fConstraintPYL[i*fgNParCh+1]=0.0;
782     fConstraintPYL[i*fgNParCh+2]=0.0;
783     fConstraintXB[i*fgNParCh+0]=0.0;
784     fConstraintXB[i*fgNParCh+1]=0.0;
785     fConstraintXB[i*fgNParCh+2]=0.0;
786     fConstraintYB[i*fgNParCh+0]=0.0;
787     fConstraintYB[i*fgNParCh+1]=0.0;
788     fConstraintYB[i*fgNParCh+2]=0.0;
789     fConstraintPB[i*fgNParCh+0]=0.0;
790     fConstraintPB[i*fgNParCh+1]=0.0;
791     fConstraintPB[i*fgNParCh+2]=0.0;
792     fConstraintXZB[i*fgNParCh+0]=0.0;
793     fConstraintXZB[i*fgNParCh+1]=0.0;
794     fConstraintXZB[i*fgNParCh+2]=0.0;
795     fConstraintYZB[i*fgNParCh+0]=0.0;
796     fConstraintYZB[i*fgNParCh+1]=0.0;
797     fConstraintYZB[i*fgNParCh+2]=0.0;
798     fConstraintPZB[i*fgNParCh+0]=0.0;
799     fConstraintPZB[i*fgNParCh+1]=0.0;
800     fConstraintPZB[i*fgNParCh+2]=0.0;
801     fConstraintXYB[i*fgNParCh+0]=0.0;
802     fConstraintXYB[i*fgNParCh+1]=0.0;
803     fConstraintXYB[i*fgNParCh+2]=0.0;
804     fConstraintYYB[i*fgNParCh+0]=0.0;
805     fConstraintYYB[i*fgNParCh+1]=0.0;
806     fConstraintYYB[i*fgNParCh+2]=0.0;
807     fConstraintPYB[i*fgNParCh+0]=0.0;
808     fConstraintPYB[i*fgNParCh+1]=0.0;
809     fConstraintPYB[i*fgNParCh+2]=0.0;
810     fConstraintXR[i*fgNParCh+0]=0.0;
811     fConstraintXR[i*fgNParCh+1]=0.0;
812     fConstraintXR[i*fgNParCh+2]=0.0;
813     fConstraintYR[i*fgNParCh+0]=0.0;
814     fConstraintYR[i*fgNParCh+1]=0.0;
815     fConstraintYR[i*fgNParCh+2]=0.0;
816     fConstraintPR[i*fgNParCh+0]=0.0;
817     fConstraintPR[i*fgNParCh+1]=0.0;
818     fConstraintPR[i*fgNParCh+2]=0.0;
819     fConstraintXZR[i*fgNParCh+0]=0.0;
820     fConstraintXZR[i*fgNParCh+1]=0.0;
821     fConstraintXZR[i*fgNParCh+2]=0.0;
822     fConstraintYZR[i*fgNParCh+0]=0.0;
823     fConstraintYZR[i*fgNParCh+1]=0.0;
824     fConstraintYZR[i*fgNParCh+2]=0.0;
825     fConstraintPZR[i*fgNParCh+0]=0.0;
826     fConstraintPZR[i*fgNParCh+1]=0.0;
827     fConstraintPZR[i*fgNParCh+2]=0.0;
828     fConstraintPZR[i*fgNParCh+0]=0.0;
829     fConstraintPZR[i*fgNParCh+1]=0.0;
830     fConstraintPZR[i*fgNParCh+2]=0.0;
831     fConstraintXYR[i*fgNParCh+0]=0.0;
832     fConstraintXYR[i*fgNParCh+1]=0.0;
833     fConstraintXYR[i*fgNParCh+2]=0.0;
834     fConstraintYYR[i*fgNParCh+0]=0.0;
835     fConstraintYYR[i*fgNParCh+1]=0.0;
836     fConstraintYYR[i*fgNParCh+2]=0.0;
837     fConstraintPYR[i*fgNParCh+0]=0.0;
838     fConstraintPYR[i*fgNParCh+1]=0.0;
839     fConstraintPYR[i*fgNParCh+2]=0.0;
840   }
841 }
842
843 //______________________________________________________________________
844 void AliMUONAlignment::AddConstraint(Double_t *par, Double_t value) {
845   /// Constrain equation defined by par to value
846   fMillepede->SetGlobalConstraint(par, value);
847   AliInfo("Adding constraint");
848 }
849
850 //______________________________________________________________________
851 void AliMUONAlignment::InitGlobalParameters(Double_t *par) {
852   /// Initialize global parameters with par array
853   fMillepede->SetGlobalParameters(par);
854   AliInfo("Init Global Parameters");
855 }
856
857 //______________________________________________________________________
858 void AliMUONAlignment::FixParameter(Int_t iPar, Double_t value) {
859   /// Parameter iPar is encourage to vary in [-value;value].
860   /// If value == 0, parameter is fixed
861   fMillepede->SetParSigma(iPar, value);
862   if (TMath::Abs(value)<1e-4) AliInfo(Form("Parameter %i Fixed", iPar));
863 }
864
865 //______________________________________________________________________
866 void AliMUONAlignment::ResetLocalEquation()
867 {
868   /// Reset the derivative vectors
869   for(int i=0; i<fNLocal; i++) {
870     fLocalDerivatives[i] = 0.0;
871   }
872   for(int i=0; i<fNGlobal; i++) {
873     fGlobalDerivatives[i] = 0.0;
874   }
875 }
876
877 //______________________________________________________________________
878 void AliMUONAlignment::AllowVariations(const Bool_t *bChOnOff)
879 {
880
881   /// Set allowed variation for selected chambers based on fDoF and fAllowVar
882   for (Int_t iCh=1; iCh<=10; iCh++)
883   {
884     if (bChOnOff[iCh-1])
885     {
886
887       Int_t iDetElemFirst = (iCh>1) ? fgSNDetElemCh[iCh-2] : 0;
888       Int_t iDetElemLast = fgSNDetElemCh[iCh-1];
889       for (int i=0; i<fgNParCh; i++)
890       {
891         AliDebug(1,Form("fDoF[%d]= %d",i,fDoF[i]));
892         if (fDoF[i])
893         {
894           for (Int_t j=iDetElemFirst; j<iDetElemLast; j++){
895             FixParameter(j*fgNParCh+i, fAllowVar[i]);
896           }
897         }
898
899       }
900
901     }
902
903   }
904
905 }
906
907 //______________________________________________________________________
908 void AliMUONAlignment::SetNonLinear(Int_t iPar  /* set non linear flag */ )
909 {
910   /// Set nonlinear flag for parameter iPar
911   fMillepede->SetNonLinear(iPar);
912   AliInfo(Form("Parameter %i set to non linear", iPar));
913 }
914
915 //______________________________________________________________________
916 void AliMUONAlignment::SetSigmaXY(Double_t sigmaX, Double_t sigmaY)
917 {
918
919   /// Set expected measurement resolution
920   fSigma[0] = sigmaX;   fSigma[1] = sigmaY;
921   AliInfo(Form("Using fSigma[0]=%f and fSigma[1]=%f",fSigma[0],fSigma[1]));
922 }
923
924
925 //______________________________________________________________________
926 void AliMUONAlignment::LocalEquationX( Bool_t doAlignment )
927 {
928
929
930   // local cluster record
931   AliMUONAlignmentClusterRecord clusterRecord;
932
933   // store detector and measurement
934   clusterRecord.SetDetElemId( fDetElemId );
935   clusterRecord.SetDetElemNumber( fDetElemNumber );
936   clusterRecord.SetMeas( fMeas[0] );
937   clusterRecord.SetSigma( fSigma[0] );
938
939   // store local derivatives
940   clusterRecord.SetLocalDerivative( 0, fCosPhi );
941   clusterRecord.SetLocalDerivative( 1, fCosPhi*(fTrackPos[2] - fTrackPos0[2]) );
942   clusterRecord.SetLocalDerivative( 2, fSinPhi );
943   clusterRecord.SetLocalDerivative( 3, fSinPhi*(fTrackPos[2] - fTrackPos0[2]) );
944
945   // store global derivatives
946   clusterRecord.SetGlobalDerivative( 0, -fCosPhi );
947   clusterRecord.SetGlobalDerivative( 1, -fSinPhi );
948
949   if (fBFieldOn)
950   {
951
952     clusterRecord.SetGlobalDerivative(
953       2,
954       -fSinPhi*(fTrackPos[0]-fDetElemPos[0])
955       +fCosPhi*(fTrackPos[1]-fDetElemPos[1]) );
956
957   } else {
958
959     clusterRecord.SetGlobalDerivative(
960       2,
961       -fSinPhi*(fTrackPos0[0]+fTrackSlope0[0]*(fTrackPos[2]-fTrackPos0[2])-fDetElemPos[0])
962       +fCosPhi*(fTrackPos0[1]+fTrackSlope0[1]*(fTrackPos[2]-fTrackPos0[2])-fDetElemPos[1]) );
963
964   }
965
966   clusterRecord.SetGlobalDerivative( 3, fCosPhi*fTrackSlope0[0]+fSinPhi*fTrackSlope0[1] );
967
968   // append to trackRecord
969   fTrackRecord.AddClusterRecord( clusterRecord );
970
971   // store local equation
972   if( doAlignment ) LocalEquation( clusterRecord );
973
974 }
975
976 //______________________________________________________________________
977 void AliMUONAlignment::LocalEquationY(Bool_t doAlignment )
978 {
979
980   // local cluster record
981   AliMUONAlignmentClusterRecord clusterRecord;
982
983   // store detector and measurement
984   clusterRecord.SetDetElemId( fDetElemId );
985   clusterRecord.SetDetElemNumber( fDetElemNumber );
986   clusterRecord.SetMeas( fMeas[1] );
987   clusterRecord.SetSigma( fSigma[1] );
988
989   // store local derivatives
990   clusterRecord.SetLocalDerivative( 0, -fSinPhi );
991   clusterRecord.SetLocalDerivative( 1, -fSinPhi*(fTrackPos[2] - fTrackPos0[2] ) );
992   clusterRecord.SetLocalDerivative( 2, fCosPhi );
993   clusterRecord.SetLocalDerivative( 3, fCosPhi*(fTrackPos[2] - fTrackPos0[2] ) );
994
995   // set global derivatives
996   clusterRecord.SetGlobalDerivative( 0,  fSinPhi);
997   clusterRecord.SetGlobalDerivative( 1, -fCosPhi);
998
999   if (fBFieldOn)
1000   {
1001
1002     clusterRecord.SetGlobalDerivative(
1003        2,
1004       -fCosPhi*(fTrackPos[0]-fDetElemPos[0])
1005       -fSinPhi*(fTrackPos[1]-fDetElemPos[1]));
1006
1007   } else {
1008
1009     clusterRecord.SetGlobalDerivative(
1010        2,
1011       -fCosPhi*(fTrackPos0[0]+fTrackSlope0[0]*(fTrackPos[2]-fTrackPos0[2])-fDetElemPos[0])
1012       -fSinPhi*(fTrackPos0[1]+fTrackSlope0[1]*(fTrackPos[2]-fTrackPos0[2])-fDetElemPos[1]));
1013   }
1014
1015   clusterRecord.SetGlobalDerivative( 3, -fSinPhi*fTrackSlope0[0]+fCosPhi*fTrackSlope0[1]);
1016
1017   // append to trackRecord
1018   fTrackRecord.AddClusterRecord( clusterRecord );
1019
1020   // store local equation
1021   if( doAlignment ) LocalEquation( clusterRecord );
1022
1023 }
1024
1025 //_____________________________________________________
1026 void AliMUONAlignment::LocalEquation( const AliMUONAlignmentClusterRecord& clusterRecord )
1027 {
1028
1029   // copy to local variables
1030   for( Int_t index = 0; index < 4; ++index )
1031   {
1032     SetLocalDerivative( index, clusterRecord.GetLocalDerivative( index ) );
1033     SetGlobalDerivative( clusterRecord.GetDetElemNumber()*fgNParCh + index, clusterRecord.GetGlobalDerivative( index ) );
1034   }
1035
1036   // pass equation parameters to millepede
1037   fMillepede->SetLocalEquation( fGlobalDerivatives, fLocalDerivatives, clusterRecord.GetMeas(), clusterRecord.GetSigma() );
1038
1039 }
1040
1041 //______________________________________________________________________
1042 void AliMUONAlignment::FillRecPointData()
1043 {
1044
1045   /// Get information of current cluster
1046   fClustPos[0] = fCluster->GetX();
1047   fClustPos[1] = fCluster->GetY();
1048   fClustPos[2] = fCluster->GetZ();
1049   fTransform->Global2Local(
1050     fDetElemId,fClustPos[0],fClustPos[1],fClustPos[2],
1051     fClustPosLoc[0],fClustPosLoc[1],fClustPosLoc[2]);
1052
1053 }
1054
1055 //______________________________________________________________________
1056 void AliMUONAlignment::FillTrackParamData()
1057 {
1058
1059   /// Get information of current track at current cluster
1060   fTrackPos[0] = fTrackParam->GetNonBendingCoor();
1061   fTrackPos[1] = fTrackParam->GetBendingCoor();
1062   fTrackPos[2] = fTrackParam->GetZ();
1063   fTrackSlope[0] = fTrackParam->GetNonBendingSlope();
1064   fTrackSlope[1] = fTrackParam->GetBendingSlope();
1065   fTransform->Global2Local(
1066     fDetElemId,fTrackPos[0],fTrackPos[1],fTrackPos[2],
1067     fTrackPosLoc[0],fTrackPosLoc[1],fTrackPosLoc[2]);
1068
1069 }
1070
1071 //_____________________________________________________
1072 void AliMUONAlignment::FillDetElemData()
1073 {
1074
1075   /// Get information of current detection element
1076   Double_t lDetElemLocX = 0.;
1077   Double_t lDetElemLocY = 0.;
1078   Double_t lDetElemLocZ = 0.;
1079   fDetElemId = fCluster->GetDetElemId();
1080   fDetElemNumber = fDetElemId%100;
1081   for (int iCh=0; iCh<fDetElemId/100-1; iCh++)
1082   { fDetElemNumber += fgNDetElemCh[iCh]; }
1083
1084   fTransform->Local2Global(
1085     fDetElemId,lDetElemLocX,lDetElemLocY,lDetElemLocZ,
1086     fDetElemPos[0],fDetElemPos[1],fDetElemPos[2]);
1087
1088 }
1089
1090 //_____________________________________________________
1091 AliMUONAlignmentTrackRecord* AliMUONAlignment::ProcessTrack( AliMUONTrack* track, Bool_t doAlignment )
1092 {
1093
1094   // store current track in running member.
1095   fTrack = track;
1096
1097   // clear track record
1098   fTrackRecord.Clear();
1099
1100   // get number of tracks
1101   Int_t nTrackParam = fTrack->GetTrackParamAtCluster()->GetEntries();
1102   AliDebug(1,Form("Number of track param entries : %i ", nTrackParam));
1103
1104   Bool_t first( kTRUE );
1105   for(Int_t iCluster=0; iCluster<nTrackParam; iCluster++)
1106   {
1107
1108     // and get new pointers
1109     fTrackParam = (AliMUONTrackParam *) fTrack->GetTrackParamAtCluster()->At(iCluster);
1110     if ( ! fTrackParam ) continue;
1111     fCluster = fTrackParam->GetClusterPtr();
1112     if ( ! fCluster ) continue;
1113     // if (fDetElemId<500) continue;
1114
1115     // fill local variables for this position --> one measurement
1116     FillDetElemData();
1117     FillRecPointData();
1118     FillTrackParamData();
1119
1120     if( first )
1121     {
1122
1123       // for first valid cluster, save track position as "starting" values
1124       first = kFALSE;
1125
1126       fTrackPos0[0] = fTrackPos[0];
1127       fTrackPos0[1] = fTrackPos[1];
1128       fTrackPos0[2] = fTrackPos[2];
1129       fTrackSlope0[0] = fTrackSlope[0];
1130       fTrackSlope0[1] = fTrackSlope[1];
1131
1132     }
1133
1134     // calculate measurements
1135     fCosPhi = TMath::Cos(fPhi);
1136     fSinPhi = TMath::Sin(fPhi);
1137     if( fBFieldOn )
1138     {
1139
1140       fMeas[0] = fTrackPos[0] - fClustPos[0];
1141       fMeas[1] = fTrackPos[1] - fClustPos[1];
1142
1143     } else {
1144
1145       fMeas[0] = - fClustPos[0];
1146       fMeas[1] = - fClustPos[1];
1147
1148     }
1149
1150     // soùe debugging
1151     AliDebug(1,Form("cluster: %i", iCluster));
1152     AliDebug(1,Form("x: %f\t y: %f\t z: %f\t DetElemID: %i\t ", fClustPos[0], fClustPos[1], fClustPos[2], fDetElemId));
1153     AliDebug(1,Form("fDetElemPos[0]: %f\t fDetElemPos[1]: %f\t fDetElemPos[2]: %f\t DetElemID: %i\t ", fDetElemPos[0],fDetElemPos[1],fDetElemPos[2], fDetElemId));
1154
1155     AliDebug(1,Form("Track Parameter: %i", iCluster));
1156     AliDebug(1,Form("x: %f\t y: %f\t z: %f\t slopex: %f\t slopey: %f", fTrackPos[0], fTrackPos[1], fTrackPos[2], fTrackSlope[0], fTrackSlope[1]));
1157     AliDebug(1,Form("x0: %f\t y0: %f\t z0: %f\t slopex0: %f\t slopey0: %f", fTrackPos0[0], fTrackPos0[1], fTrackPos0[2], fTrackSlope0[0], fTrackSlope0[1]));
1158
1159     AliDebug(1,Form("fMeas[0]: %f\t fMeas[1]: %f\t fSigma[0]: %f\t fSigma[1]: %f", fMeas[0], fMeas[1], fSigma[0], fSigma[1]));
1160
1161     // Set local equations
1162     LocalEquationX( doAlignment );
1163     LocalEquationY( doAlignment );
1164
1165   }
1166
1167   return &fTrackRecord;
1168
1169 }
1170
1171 //______________________________________________________________________________
1172 void AliMUONAlignment::ProcessTrack( AliMUONAlignmentTrackRecord* track, Bool_t doAlignment )
1173 {
1174
1175   if( !( track && doAlignment ) ) return;
1176
1177   // loop over clusters
1178   for( Int_t index = 0; index < track->GetNRecords(); ++index )
1179   { if( AliMUONAlignmentClusterRecord* record = track->GetRecord( index ) ) LocalEquation( *record ); }
1180
1181   return;
1182
1183 }
1184
1185 //_____________________________________________________
1186 void AliMUONAlignment::LocalFit(Int_t iTrack, Double_t *lTrackParam, Int_t lSingleFit)
1187 {
1188
1189   /// Call local fit for this tracks
1190   Int_t iRes = fMillepede->LocalFit(iTrack,lTrackParam,lSingleFit);
1191   if (iRes && !lSingleFit)
1192   { fMillepede->SetNLocalEquations(fMillepede->GetNLocalEquations()+1); }
1193
1194 }
1195
1196 //_____________________________________________________
1197 void AliMUONAlignment::GlobalFit(Double_t *parameters,Double_t *errors,Double_t *pulls)
1198 {
1199
1200   /// Call global fit; Global parameters are stored in parameters
1201   fMillepede->GlobalFit(parameters,errors,pulls);
1202
1203   AliInfo("Done fitting global parameters!");
1204   for (int iGlob=0; iGlob<fgNDetElem; iGlob++)
1205   { printf("%d\t %f\t %f\t %f\t %f \n",iGlob,parameters[iGlob*fgNParCh+0],parameters[iGlob*fgNParCh+1],parameters[iGlob*fgNParCh+3],parameters[iGlob*fgNParCh+2]); }
1206
1207
1208 }
1209
1210 //_____________________________________________________
1211 Double_t AliMUONAlignment::GetParError(Int_t iPar)
1212 {
1213   /// Get error of parameter iPar
1214   Double_t lErr = fMillepede->GetParError(iPar);
1215   return lErr;
1216 }
1217
1218 //_____________________________________________________
1219 void AliMUONAlignment::PrintGlobalParameters()
1220 {
1221   /// Print global parameters
1222   fMillepede->PrintGlobalParameters();
1223 }
1224
1225 //_________________________________________________________________________
1226 TGeoCombiTrans AliMUONAlignment::ReAlign(const TGeoCombiTrans & transform, const double *lMisAlignment) const
1227 {
1228   /// Realign given transformation by given misalignment and return the misaligned transformation
1229
1230   Double_t cartMisAlig[3] = {0,0,0};
1231   Double_t angMisAlig[3] = {0,0,0};
1232 //   const Double_t *trans = transform.GetTranslation();
1233 //   TGeoRotation *rot;
1234 //   // check if the rotation we obtain is not NULL
1235 //   if (transform.GetRotation()) {
1236 //     rot = transform.GetRotation();
1237 //   }
1238 //   else {
1239 //     rot = new TGeoRotation("rot");
1240 //   }                  // default constructor.
1241
1242   cartMisAlig[0] = -TMath::Sign(1.0,transform.GetRotationMatrix()[0])*lMisAlignment[0];
1243   cartMisAlig[1] = -TMath::Sign(1.0,transform.GetRotationMatrix()[4])*lMisAlignment[1];
1244   cartMisAlig[2] = -TMath::Sign(1.0,transform.GetRotationMatrix()[8])*lMisAlignment[3];
1245   angMisAlig[2] = -TMath::Sign(1.0,transform.GetRotationMatrix()[0]*transform.GetRotationMatrix()[4])*lMisAlignment[2]*180./TMath::Pi();
1246
1247   TGeoTranslation deltaTrans(cartMisAlig[0], cartMisAlig[1], cartMisAlig[2]);
1248   TGeoRotation deltaRot;
1249   deltaRot.RotateX(angMisAlig[0]);
1250   deltaRot.RotateY(angMisAlig[1]);
1251   deltaRot.RotateZ(angMisAlig[2]);
1252
1253   TGeoCombiTrans deltaTransf(deltaTrans,deltaRot);
1254   TGeoHMatrix newTransfMat = transform * deltaTransf;
1255
1256   return TGeoCombiTrans(newTransfMat);
1257 }
1258
1259 //______________________________________________________________________
1260 AliMUONGeometryTransformer *
1261 AliMUONAlignment::ReAlign(const AliMUONGeometryTransformer * transformer,
1262           const double *misAlignments, Bool_t verbose)
1263
1264 {
1265   /// Returns a new AliMUONGeometryTransformer with the found misalignments
1266   /// applied.
1267
1268   // Takes the internal geometry module transformers, copies them
1269   // and gets the Detection Elements from them.
1270   // Takes misalignment parameters and applies these
1271   // to the local transform of the Detection Element
1272   // Obtains the global transform by multiplying the module transformer
1273   // transformation with the local transformation
1274   // Applies the global transform to a new detection element
1275   // Adds the new detection element to a new module transformer
1276   // Adds the new module transformer to a new geometry transformer
1277   // Returns the new geometry transformer
1278
1279   Double_t lModuleMisAlignment[4] = {0.,0.,0.,0.};
1280   Double_t lDetElemMisAlignment[4] = {0.,0.,0.,0.};
1281   Int_t iDetElemId = 0;
1282   Int_t iDetElemNumber = 0;
1283
1284   AliMUONGeometryTransformer *newGeometryTransformer =
1285     new AliMUONGeometryTransformer();
1286   for (Int_t iMt = 0; iMt < transformer->GetNofModuleTransformers(); iMt++) {
1287     // module transformers
1288     const AliMUONGeometryModuleTransformer *kModuleTransformer =
1289       transformer->GetModuleTransformer(iMt, true);
1290
1291     AliMUONGeometryModuleTransformer *newModuleTransformer =
1292       new AliMUONGeometryModuleTransformer(iMt);
1293     newGeometryTransformer->AddModuleTransformer(newModuleTransformer);
1294
1295     TGeoCombiTrans moduleTransform =
1296       TGeoCombiTrans(*kModuleTransformer->GetTransformation());
1297     // New module transformation
1298     TGeoCombiTrans newModuleTransform = ReAlign(moduleTransform,lModuleMisAlignment);
1299     newModuleTransformer->SetTransformation(newModuleTransform);
1300
1301     // Get delta transformation:
1302     // Tdelta = Tnew * Told.inverse
1303     TGeoHMatrix deltaModuleTransform =
1304       AliMUONGeometryBuilder::Multiply(newModuleTransform,
1305                kModuleTransformer->GetTransformation()->Inverse());
1306     // Create module mis alignment matrix
1307     newGeometryTransformer
1308       ->AddMisAlignModule(kModuleTransformer->GetModuleId(), deltaModuleTransform);
1309
1310     AliMpExMap *detElements = kModuleTransformer->GetDetElementStore();
1311
1312     if (verbose)
1313       AliInfo(Form("%i DEs in old GeometryStore  %i",detElements->GetSize(), iMt));
1314
1315     TIter next(detElements->CreateIterator());
1316     AliMUONGeometryDetElement* detElement;
1317     Int_t iDe(-1);
1318     while ( ( detElement = static_cast<AliMUONGeometryDetElement*>(next()) ) )
1319     {
1320       ++iDe;
1321       // make a new detection element
1322       AliMUONGeometryDetElement *newDetElement =
1323   new AliMUONGeometryDetElement(detElement->GetId(),
1324               detElement->GetVolumePath());
1325       TString lDetElemName(detElement->GetDEName());
1326       lDetElemName.ReplaceAll("DE","");
1327       iDetElemId = lDetElemName.Atoi();
1328       iDetElemNumber = iDetElemId%100;
1329       for (int iCh=0; iCh<iDetElemId/100-1; iCh++){
1330   iDetElemNumber += fgNDetElemCh[iCh];
1331       }
1332       for (int i=0; i<fgNParCh; i++) {
1333   lDetElemMisAlignment[i] = 0.0;
1334   if (iMt<fgNTrkMod) {
1335     AliInfo(Form("iMt %i, iCh %i, iDe %i, iDeId %i, iDeNb %i, iPar %i",iMt, iDetElemId/100, iDe, iDetElemId, iDetElemNumber, iDetElemNumber*fgNParCh+i));
1336     lDetElemMisAlignment[i] =  misAlignments[iDetElemNumber*fgNParCh+i];
1337   }
1338       }
1339       // local transformation of this detection element.
1340       TGeoCombiTrans localTransform
1341   = TGeoCombiTrans(*detElement->GetLocalTransformation());
1342       TGeoCombiTrans newLocalTransform = ReAlign(localTransform,lDetElemMisAlignment);
1343       newDetElement->SetLocalTransformation(newLocalTransform);
1344
1345       // global transformation
1346       TGeoHMatrix newGlobalTransform =
1347   AliMUONGeometryBuilder::Multiply(newModuleTransform,
1348            newLocalTransform);
1349       newDetElement->SetGlobalTransformation(newGlobalTransform);
1350
1351       // add this det element to module
1352       newModuleTransformer->GetDetElementStore()->Add(newDetElement->GetId(),
1353                   newDetElement);
1354
1355       // In the Alice Alignment Framework misalignment objects store
1356       // global delta transformation
1357       // Get detection "intermediate" global transformation
1358       TGeoHMatrix newOldGlobalTransform = newModuleTransform * localTransform;
1359       // Get detection element global delta transformation:
1360       // Tdelta = Tnew * Told.inverse
1361       TGeoHMatrix  deltaGlobalTransform
1362   = AliMUONGeometryBuilder::Multiply(newGlobalTransform,
1363              newOldGlobalTransform.Inverse());
1364
1365       // Create mis alignment matrix
1366       newGeometryTransformer
1367   ->AddMisAlignDetElement(detElement->GetId(), deltaGlobalTransform);
1368     }
1369
1370     if (verbose)
1371       AliInfo(Form("Added module transformer %i to the transformer", iMt));
1372     newGeometryTransformer->AddModuleTransformer(newModuleTransformer);
1373   }
1374   return newGeometryTransformer;
1375 }
1376
1377 //______________________________________________________________________
1378 void AliMUONAlignment::SetAlignmentResolution(const TClonesArray* misAlignArray, Int_t rChId, Double_t rChResX, Double_t rChResY, Double_t rDeResX, Double_t rDeResY)
1379 {
1380
1381   /// Set alignment resolution to misalign objects to be stored in CDB
1382   Int_t chIdMin = (rChId<0)? 0 : rChId;
1383   Int_t chIdMax = (rChId<0)? 9 : rChId;
1384   Double_t chResX = rChResX;
1385   Double_t chResY = rChResY;
1386   Double_t deResX = rDeResX;
1387   Double_t deResY = rDeResY;
1388
1389   TMatrixDSym mChCorrMatrix(6);
1390   mChCorrMatrix[0][0]=chResX*chResX;
1391   mChCorrMatrix[1][1]=chResY*chResY;
1392   //  mChCorrMatrix.Print();
1393
1394   TMatrixDSym mDECorrMatrix(6);
1395   mDECorrMatrix[0][0]=deResX*deResX;
1396   mDECorrMatrix[1][1]=deResY*deResY;
1397   //  mDECorrMatrix.Print();
1398
1399   AliAlignObjMatrix *alignMat = 0x0;
1400
1401   for(Int_t chId=chIdMin; chId<=chIdMax; chId++) {
1402     TString chName1;
1403     TString chName2;
1404     if (chId<4){
1405       chName1 = Form("GM%d",chId);
1406       chName2 = Form("GM%d",chId);
1407     } else {
1408       chName1 = Form("GM%d",4+(chId-4)*2);
1409       chName2 = Form("GM%d",4+(chId-4)*2+1);
1410     }
1411
1412     for (int i=0; i<misAlignArray->GetEntries(); i++) {
1413       alignMat = (AliAlignObjMatrix*)misAlignArray->At(i);
1414       TString volName(alignMat->GetSymName());
1415       if((volName.Contains(chName1)&&
1416     ((volName.Last('/')==volName.Index(chName1)+chName1.Length())||
1417      (volName.Length()==volName.Index(chName1)+chName1.Length())))||
1418    (volName.Contains(chName2)&&
1419     ((volName.Last('/')==volName.Index(chName2)+chName2.Length())||
1420      (volName.Length()==volName.Index(chName2)+chName2.Length())))){
1421   volName.Remove(0,volName.Last('/')+1);
1422   if (volName.Contains("GM")) {
1423     //  alignMat->Print("NULL");
1424     alignMat->SetCorrMatrix(mChCorrMatrix);
1425   } else if (volName.Contains("DE")) {
1426     //  alignMat->Print("NULL");
1427     alignMat->SetCorrMatrix(mDECorrMatrix);
1428   }
1429       }
1430     }
1431   }
1432 }