]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/AliMUONAlignment.cxx
PreReading of MC information on demand.
[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 "AliMUONConstants.h"
39 #include "AliMillepede.h"
40
41 #include "AliMpExMap.h"
42 #include "AliMpExMapIterator.h"
43
44 #include "AliAlignObjMatrix.h"
45 #include "AliLog.h"
46
47 #include "TMath.h"
48 #include "TMatrixDSym.h"
49 #include "TSystem.h"
50
51 /// \cond CLASSIMP
52 ClassImp(AliMUONAlignment)
53 /// \endcond
54
55   Int_t AliMUONAlignment::fgNDetElem = 4*2+4*2+18*2+26*2+26*2;
56   Int_t AliMUONAlignment::fgNDetElemCh[10] = {4,4,4,4,18,18,26,26,26,26};
57   Int_t AliMUONAlignment::fgSNDetElemCh[10] = {4,8,12,16,34,52,78,104,130,156};
58   Int_t AliMUONAlignment::fgNParCh = 4;
59   Int_t AliMUONAlignment::fgNTrkMod = 16;
60   Int_t AliMUONAlignment::fgNCh = 10;
61   Int_t AliMUONAlignment::fgNSt = 5;
62
63 AliMUONAlignment::AliMUONAlignment() 
64   : TObject(),
65     fBFieldOn(kTRUE),
66     fStartFac(256.), 
67     fResCutInitial(100.), 
68     fResCut(100.),
69     fMillepede(0),
70     fTrackParamAtCluster(0),
71     fTrack(0),
72     fCluster(0),
73     fTrackParam(0),
74     fNGlobal(fgNDetElem*fgNParCh),
75     fNLocal(4),
76     fNStdDev(3),
77     fDetElemId(0),
78     fDetElemNumber(0),
79     fPhi(0.),
80     fCosPhi(1.),
81     fSinPhi(0.),
82     fTransform(0)
83 {
84   /// Default constructor, setting define alignment parameters
85   fSigma[0] = 1.0e-1;
86   fSigma[1] = 1.0e-2;
87
88   fDoF[0] = kTRUE;  fDoF[1] = kTRUE;  fDoF[2] = kTRUE;  fDoF[3] = kTRUE;
89   fAllowVar[0] = 0.05;  fAllowVar[1] = 0.05;  fAllowVar[2] = 0.001;  fAllowVar[3] = 0.5;
90   
91   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));
92
93   fMillepede = new AliMillepede();
94
95   Init(fNGlobal, fNLocal, fNStdDev);
96
97   ResetLocalEquation();
98   AliInfo("Parameters initialized to zero");
99
100 }
101
102 AliMUONAlignment::~AliMUONAlignment() {
103   /// Destructor
104 }
105
106 void AliMUONAlignment::Init(Int_t nGlobal,  /* number of global paramers */
107                             Int_t nLocal,   /* number of local parameters */
108                             Int_t nStdDev   /* std dev cut */ )
109 {
110   /// Initialization of AliMillepede. Fix parameters, define constraints ...
111   fMillepede->InitMille(nGlobal,nLocal,nStdDev,fResCut,fResCutInitial);
112
113 //  Bool_t bStOnOff[5] = {kTRUE,kTRUE,kTRUE,kTRUE,kTRUE};
114 //  Bool_t bChOnOff[10] = {kTRUE,kTRUE,kTRUE,kTRUE,kTRUE,kTRUE,kTRUE,kTRUE,kTRUE,kTRUE};
115 //  Bool_t bSpecLROnOff[2] = {kTRUE,kTRUE};
116
117 //   AllowVariations(bChOnOff);
118
119   // Fix parameters or add constraints here
120 //   for (Int_t iSt=0; iSt<5; iSt++)
121 //     if (!bStOnOff[iSt]) FixStation(iSt+1);
122 //   for (Int_t iCh=0; iCh<10; iCh++)
123 //     if (!bChOnOff[iCh]) FixChamber(iCh+1);
124
125 //   FixHalfSpectrometer(bChOnOff,bSpecLROnOff);
126
127   ResetConstraints();
128   
129   // Define global constrains to be applied
130   // X, Y, P, XvsZ, YvsZ, PvsZ, XvsY, YvsY, PvsY
131   Bool_t bVarXYT[9] = {kTRUE,kTRUE,kTRUE,kTRUE,kTRUE,kTRUE,kTRUE,kTRUE,kTRUE};
132   Bool_t bDetTLBR[4] = {kFALSE,kTRUE,kFALSE,kTRUE};
133   //  AddConstraints(bStOnOff,bVarXYT,bDetTLBR,bSpecLROnOff);
134
135   // Other possible way to add constrains
136   bVarXYT[0] = kFALSE; bVarXYT[1] = kFALSE; bVarXYT[2] = kTRUE;
137   bDetTLBR[0] = kFALSE; bDetTLBR[1] = kTRUE; bDetTLBR[2] = kFALSE; bDetTLBR[3] = kFALSE;
138 //   AddConstraints(bStOnOff,bVarXYT,bDetTLBR);
139
140   bVarXYT[0] = kTRUE; bVarXYT[1] = kTRUE; bVarXYT[2] = kFALSE;
141   //  AddConstraints(bStOnOff,bVarXYT);
142   
143   // Set iterations
144   if (fStartFac>1) fMillepede->SetIterations(fStartFac);          
145 }
146
147 void AliMUONAlignment::FixStation(Int_t iSt){
148   /// Fix all detection elements of station iSt
149   Int_t iDetElemFirst = (iSt>1) ? fgSNDetElemCh[2*(iSt-1)-1] : 0; 
150   Int_t iDetElemLast = fgSNDetElemCh[2*(iSt)-1]; 
151   for (Int_t i = iDetElemFirst; i < iDetElemLast; i++){    
152     FixParameter(i*fgNParCh+0, 0.0);
153     FixParameter(i*fgNParCh+1, 0.0);
154     FixParameter(i*fgNParCh+2, 0.0);
155     FixParameter(i*fgNParCh+3, 0.0);
156   }
157 }
158
159 void AliMUONAlignment::FixChamber(Int_t iCh){
160   /// Fix all detection elements of chamber iCh
161   Int_t iDetElemFirst = (iCh>1) ? fgSNDetElemCh[iCh-2] : 0; 
162   Int_t iDetElemLast = fgSNDetElemCh[iCh-1]; 
163   for (Int_t i = iDetElemFirst; i < iDetElemLast; i++){    
164     FixParameter(i*fgNParCh+0, 0.0);
165     FixParameter(i*fgNParCh+1, 0.0);
166     FixParameter(i*fgNParCh+2, 0.0);
167     FixParameter(i*fgNParCh+3, 0.0);
168   }
169 }
170
171 void AliMUONAlignment::FixDetElem(Int_t iDetElemId, TString sVarXYT){
172   /// Fix a given detection element
173   Int_t iDetElemNumber = iDetElemId%100;
174   for (int iCh=0; iCh<iDetElemId/100-1; iCh++){
175     iDetElemNumber += fgNDetElemCh[iCh];
176   }
177   if (sVarXYT.Contains("X")) { // X constraint
178     FixParameter(iDetElemNumber*fgNParCh+0, 0.0);
179   }
180   if (sVarXYT.Contains("Y")) { // Y constraint
181     FixParameter(iDetElemNumber*fgNParCh+1, 0.0);
182   }
183   if (sVarXYT.Contains("T")) { // T constraint
184     FixParameter(iDetElemNumber*fgNParCh+2, 0.0);
185   }
186   if (sVarXYT.Contains("Z")) { // T constraint
187     FixParameter(iDetElemNumber*fgNParCh+3, 0.0);
188   }
189 }
190
191 void AliMUONAlignment::FixHalfSpectrometer(Bool_t *lChOnOff, Bool_t *lSpecLROnOff){
192   /// Fix left or right detector
193   for (Int_t i = 0; i < fgNDetElem; i++){    
194     Int_t iCh=0;
195     for (iCh=1; iCh<=fgNCh; iCh++){
196       if (i<fgSNDetElemCh[iCh-1]) break;
197     }
198     if (lChOnOff[iCh-1]){
199       Int_t lDetElemNumber = (iCh==1) ? i : i-fgSNDetElemCh[iCh-2];
200       if (iCh>=1 && iCh<=4){
201         if ((lDetElemNumber==1 || lDetElemNumber==2) && !lSpecLROnOff[0]){ // From track crossings
202           FixParameter(i*fgNParCh+0, 0.0);
203           FixParameter(i*fgNParCh+1, 0.0);
204           FixParameter(i*fgNParCh+2, 0.0);
205           FixParameter(i*fgNParCh+3, 0.0);
206         }
207         if ((lDetElemNumber==0 || lDetElemNumber==3) && !lSpecLROnOff[1]){ // From track crossings
208           FixParameter(i*fgNParCh+0, 0.0);
209           FixParameter(i*fgNParCh+1, 0.0);
210           FixParameter(i*fgNParCh+2, 0.0);
211           FixParameter(i*fgNParCh+3, 0.0);
212         }
213       }
214       if (iCh>=5 && iCh<=6){
215         if ((lDetElemNumber>=5&&lDetElemNumber<=13) && !lSpecLROnOff[0]){
216           FixParameter(i*fgNParCh+0, 0.0);
217           FixParameter(i*fgNParCh+1, 0.0);
218           FixParameter(i*fgNParCh+2, 0.0);
219           FixParameter(i*fgNParCh+3, 0.0);
220         }
221         if (((lDetElemNumber>=0&&lDetElemNumber<=4) || 
222              (lDetElemNumber>=14&&lDetElemNumber<=17)) && !lSpecLROnOff[1]){
223           FixParameter(i*fgNParCh+0, 0.0);
224           FixParameter(i*fgNParCh+1, 0.0);
225           FixParameter(i*fgNParCh+2, 0.0);
226           FixParameter(i*fgNParCh+3, 0.0);
227         }
228       }
229       if (iCh>=7 && iCh<=10){
230         if ((lDetElemNumber>=7&&lDetElemNumber<=19) && !lSpecLROnOff[0]){
231           FixParameter(i*fgNParCh+0, 0.0);
232           FixParameter(i*fgNParCh+1, 0.0);
233           FixParameter(i*fgNParCh+2, 0.0);
234           FixParameter(i*fgNParCh+3, 0.0);
235         }
236         if (((lDetElemNumber>=0&&lDetElemNumber<=6) || 
237              (lDetElemNumber>=20&&lDetElemNumber<=25)) && !lSpecLROnOff[1]){
238           FixParameter(i*fgNParCh+0, 0.0);
239           FixParameter(i*fgNParCh+1, 0.0);
240           FixParameter(i*fgNParCh+2, 0.0);
241           FixParameter(i*fgNParCh+3, 0.0);
242         }
243       }
244     }
245   }
246 }
247
248 void AliMUONAlignment::SetNonLinear(Bool_t *lChOnOff,Bool_t *lVarXYT){
249   /// Set non linear parameter flag selected chambers and degrees of freedom
250   for (Int_t i = 0; i < fgNDetElem; i++){    
251     Int_t iCh=0;
252     for (iCh=1; iCh<=fgNCh; iCh++){
253       if (i<fgSNDetElemCh[iCh-1]) break;
254     }
255     if (lChOnOff[iCh-1]){
256       if (lVarXYT[0]) { // X constraint
257         SetNonLinear(i*fgNParCh+0);
258       }
259       if (lVarXYT[1]) { // Y constraint
260         SetNonLinear(i*fgNParCh+1);
261       }
262       if (lVarXYT[2]) { // T constraint
263         SetNonLinear(i*fgNParCh+2);
264       }
265       if (lVarXYT[3]) { // Z constraint
266         SetNonLinear(i*fgNParCh+3);
267       }
268     }
269   }
270 }
271
272 void AliMUONAlignment::AddConstraints(Bool_t *lChOnOff,Bool_t *lVarXYT){
273   /// Add constraint equations for selected chambers and degrees of freedom 
274   for (Int_t i = 0; i < fgNDetElem; i++){    
275     Int_t iCh=0;
276     for (iCh=1; iCh<=fgNCh; iCh++){
277       if (i<fgSNDetElemCh[iCh-1]) break;
278     }
279     if (lChOnOff[iCh-1]){
280       if (lVarXYT[0]) { // X constraint
281         fConstraintX[i*fgNParCh+0]=1.0;
282       }
283       if (lVarXYT[1]) { // Y constraint
284         fConstraintY[i*fgNParCh+1]=1.0;
285       }
286       if (lVarXYT[2]) { // T constraint
287         fConstraintP[i*fgNParCh+2]=1.0;
288       }
289 //       if (lVarXYT[3]) { // Z constraint
290 //      fConstraintP[i*fgNParCh+3]=1.0;
291 //       }
292     }
293   }
294   if (lVarXYT[0]) { // X constraint
295     AddConstraint(fConstraintX,0.0);
296   }
297   if (lVarXYT[1]) { // Y constraint
298     AddConstraint(fConstraintY,0.0);
299   }
300   if (lVarXYT[2]) { // T constraint
301     AddConstraint(fConstraintP,0.0);
302   }
303 //   if (lVarXYT[3]) { // Z constraint
304 //     AddConstraint(fConstraintP,0.0);
305 //   }
306 }
307
308 void AliMUONAlignment::AddConstraints(Bool_t *lChOnOff,Bool_t *lVarXYT, Bool_t *lDetTLBR, Bool_t *lSpecLROnOff){
309   /// Add constraint equations for selected chambers, degrees of freedom and detector half 
310   Double_t lDetElemLocX = 0.;
311   Double_t lDetElemLocY = 0.;
312   Double_t lDetElemLocZ = 0.;
313   Double_t lDetElemGloX = 0.;
314   Double_t lDetElemGloY = 0.;
315   Double_t lDetElemGloZ = 0.;
316   Double_t lMeanY = 0.;
317   Double_t lSigmaY = 0.;
318   Double_t lMeanZ = 0.;
319   Double_t lSigmaZ = 0.;
320   Int_t lNDetElem = 0;
321   for (Int_t i = 0; i < fgNDetElem; i++){    
322     Int_t iCh=0;
323     for (iCh=1; iCh<=fgNCh; iCh++){
324       if (i<fgSNDetElemCh[iCh-1]) break;
325     }
326     if (lChOnOff[iCh-1]){ 
327       Int_t lDetElemNumber = (iCh==1) ? i : i-fgSNDetElemCh[iCh-2];
328       Int_t lDetElemId = iCh*100+lDetElemNumber;
329       fTransform->Local2Global(lDetElemId,lDetElemLocX,lDetElemLocY,lDetElemLocZ,
330                                lDetElemGloX,lDetElemGloY,lDetElemGloZ);
331       if (iCh>=1 && iCh<=4){
332         if ((lDetElemNumber==1 || lDetElemNumber==2) && lSpecLROnOff[0]){ // From track crossings
333           lMeanY += lDetElemGloY;
334           lSigmaY += lDetElemGloY*lDetElemGloY;
335           lMeanZ += lDetElemGloZ;
336           lSigmaZ += lDetElemGloZ*lDetElemGloZ;
337           lNDetElem++;
338         }
339         if ((lDetElemNumber==0 || lDetElemNumber==3) && lSpecLROnOff[1]){ // From track crossings
340           lMeanY += lDetElemGloY;
341           lSigmaY += lDetElemGloY*lDetElemGloY;
342           lMeanZ += lDetElemGloZ;
343           lSigmaZ += lDetElemGloZ*lDetElemGloZ;
344           lNDetElem++;
345         }
346       }
347       if (iCh>=5 && iCh<=6){
348         if ((lDetElemNumber>=5&&lDetElemNumber<=13) && lSpecLROnOff[0]){
349           lMeanY += lDetElemGloY;
350           lSigmaY += lDetElemGloY*lDetElemGloY;
351           lMeanZ += lDetElemGloZ;
352           lSigmaZ += lDetElemGloZ*lDetElemGloZ;
353           lNDetElem++;
354         }
355         if (((lDetElemNumber>=0&&lDetElemNumber<=4) || 
356              (lDetElemNumber>=14&&lDetElemNumber<=17)) && lSpecLROnOff[1]){
357           lMeanY += lDetElemGloY;
358           lSigmaY += lDetElemGloY*lDetElemGloY;
359           lMeanZ += lDetElemGloZ;
360           lSigmaZ += lDetElemGloZ*lDetElemGloZ;
361           lNDetElem++;
362         }
363       }
364       if (iCh>=7 && iCh<=10){
365         if ((lDetElemNumber>=7&&lDetElemNumber<=19) && lSpecLROnOff[0]){
366           lMeanY += lDetElemGloY;
367           lSigmaY += lDetElemGloY*lDetElemGloY;
368           lMeanZ += lDetElemGloZ;
369           lSigmaZ += lDetElemGloZ*lDetElemGloZ;
370           lNDetElem++;
371         }
372         if (((lDetElemNumber>=0&&lDetElemNumber<=6) || 
373              (lDetElemNumber>=20&&lDetElemNumber<=25)) && lSpecLROnOff[1]){
374           lMeanY += lDetElemGloY;
375           lSigmaY += lDetElemGloY*lDetElemGloY;
376           lMeanZ += lDetElemGloZ;
377           lSigmaZ += lDetElemGloZ*lDetElemGloZ;
378           lNDetElem++;
379         }
380       }
381     }
382   }
383   lMeanY /= lNDetElem;
384   lSigmaY /= lNDetElem;
385   lSigmaY = TMath::Sqrt(lSigmaY-lMeanY*lMeanY);
386   lMeanZ /= lNDetElem;
387   lSigmaZ /= lNDetElem;
388   lSigmaZ = TMath::Sqrt(lSigmaZ-lMeanZ*lMeanZ);
389   AliInfo(Form("Used %i DetElem, MeanZ= %f , SigmaZ= %f", lNDetElem,lMeanZ,lSigmaZ));  
390
391   for (Int_t i = 0; i < fgNDetElem; i++){    
392     Int_t iCh=0;
393     for (iCh=1; iCh<=fgNCh; iCh++){
394       if (i<fgSNDetElemCh[iCh-1]) break;
395     }
396     if (lChOnOff[iCh-1]){
397       Int_t lDetElemNumber = (iCh==1) ? i : i-fgSNDetElemCh[iCh-2];
398       Int_t lDetElemId = iCh*100+lDetElemNumber;
399       fTransform->Local2Global(lDetElemId,lDetElemLocX,lDetElemLocY,lDetElemLocZ,
400                                lDetElemGloX,lDetElemGloY,lDetElemGloZ);
401       if (lVarXYT[0]) { // X constraint
402         if (lDetTLBR[0]) ConstrainT(i,iCh,fConstraintXT,0); // Top half
403         if (lDetTLBR[1]) ConstrainL(i,iCh,fConstraintXL,0); // Left half
404         if (lDetTLBR[2]) ConstrainB(i,iCh,fConstraintXB,0); // Bottom half
405         if (lDetTLBR[3]) ConstrainR(i,iCh,fConstraintXR,0); // Right half
406       }
407       if (lVarXYT[1]) { // Y constraint
408         if (lDetTLBR[0]) ConstrainT(i,iCh,fConstraintYT,1); // Top half
409         if (lDetTLBR[1]) ConstrainL(i,iCh,fConstraintYL,1); // Left half
410         if (lDetTLBR[2]) ConstrainB(i,iCh,fConstraintYB,1); // Bottom half
411         if (lDetTLBR[3]) ConstrainR(i,iCh,fConstraintYR,1); // Right half
412       }
413       if (lVarXYT[2]) { // P constraint
414         if (lDetTLBR[0]) ConstrainT(i,iCh,fConstraintPT,2); // Top half
415         if (lDetTLBR[1]) ConstrainL(i,iCh,fConstraintPL,2); // Left half
416         if (lDetTLBR[2]) ConstrainB(i,iCh,fConstraintPB,2); // Bottom half
417         if (lDetTLBR[3]) ConstrainR(i,iCh,fConstraintPR,2); // Right half
418       }
419       if (lVarXYT[3]) { // X-Z shearing
420         if (lDetTLBR[0]) ConstrainT(i,iCh,fConstraintXZT,0,(lDetElemGloZ-lMeanZ)/lSigmaZ); // Top half
421         if (lDetTLBR[1]) ConstrainL(i,iCh,fConstraintXZL,0,(lDetElemGloZ-lMeanZ)/lSigmaZ); // Left half
422         if (lDetTLBR[2]) ConstrainB(i,iCh,fConstraintXZB,0,(lDetElemGloZ-lMeanZ)/lSigmaZ); // Bottom half
423         if (lDetTLBR[3]) ConstrainR(i,iCh,fConstraintXZR,0,(lDetElemGloZ-lMeanZ)/lSigmaZ); // Right half
424       }
425       if (lVarXYT[4]) { // Y-Z shearing
426         if (lDetTLBR[0]) ConstrainT(i,iCh,fConstraintYZT,1,(lDetElemGloZ-lMeanZ)/lSigmaZ); // Top half
427         if (lDetTLBR[1]) ConstrainL(i,iCh,fConstraintYZL,1,(lDetElemGloZ-lMeanZ)/lSigmaZ); // Left half
428         if (lDetTLBR[2]) ConstrainB(i,iCh,fConstraintYZB,1,(lDetElemGloZ-lMeanZ)/lSigmaZ); // Bottom half
429         if (lDetTLBR[3]) ConstrainR(i,iCh,fConstraintYZR,1,(lDetElemGloZ-lMeanZ)/lSigmaZ); // Right half
430       }
431       if (lVarXYT[5]) { // P-Z rotation
432         if (lDetTLBR[0]) ConstrainT(i,iCh,fConstraintPZT,2,(lDetElemGloZ-lMeanZ)/lSigmaZ); // Top half
433         if (lDetTLBR[1]) ConstrainL(i,iCh,fConstraintPZL,2,(lDetElemGloZ-lMeanZ)/lSigmaZ); // Left half
434         if (lDetTLBR[2]) ConstrainB(i,iCh,fConstraintPZB,2,(lDetElemGloZ-lMeanZ)/lSigmaZ); // Bottom half
435         if (lDetTLBR[3]) ConstrainR(i,iCh,fConstraintPZR,2,(lDetElemGloZ-lMeanZ)/lSigmaZ); // Right half
436       }
437       if (lVarXYT[6]) { // X-Y shearing
438         if (lDetTLBR[0]) ConstrainT(i,iCh,fConstraintXYT,0,(lDetElemGloY-lMeanY)/lSigmaY); // Top half
439         if (lDetTLBR[1]) ConstrainL(i,iCh,fConstraintXYL,0,(lDetElemGloY-lMeanY)/lSigmaY); // Left half
440         if (lDetTLBR[2]) ConstrainB(i,iCh,fConstraintXYB,0,(lDetElemGloY-lMeanY)/lSigmaY); // Bottom half
441         if (lDetTLBR[3]) ConstrainR(i,iCh,fConstraintXYR,0,(lDetElemGloY-lMeanY)/lSigmaY); // Right half
442       }
443       if (lVarXYT[7]) { // Y-Y scaling
444         if (lDetTLBR[0]) ConstrainT(i,iCh,fConstraintYYT,1,(lDetElemGloY-lMeanY)/lSigmaY); // Top half
445         if (lDetTLBR[1]) ConstrainL(i,iCh,fConstraintYYL,1,(lDetElemGloY-lMeanY)/lSigmaY); // Left half
446         if (lDetTLBR[2]) ConstrainB(i,iCh,fConstraintYYB,1,(lDetElemGloY-lMeanY)/lSigmaY); // Bottom half
447         if (lDetTLBR[3]) ConstrainR(i,iCh,fConstraintYYR,1,(lDetElemGloY-lMeanY)/lSigmaY); // Right half
448       }
449       if (lVarXYT[8]) { // P-Y rotation
450         if (lDetTLBR[0]) ConstrainT(i,iCh,fConstraintPYT,2,(lDetElemGloY-lMeanY)/lSigmaY); // Top half
451         if (lDetTLBR[1]) ConstrainL(i,iCh,fConstraintPYL,2,(lDetElemGloY-lMeanY)/lSigmaY); // Left half
452         if (lDetTLBR[2]) ConstrainB(i,iCh,fConstraintPYB,2,(lDetElemGloY-lMeanY)/lSigmaY); // Bottom half
453         if (lDetTLBR[3]) ConstrainR(i,iCh,fConstraintPYR,2,(lDetElemGloY-lMeanY)/lSigmaY); // Right half
454       }
455     }
456   }
457   if (lVarXYT[0]) { // X constraint
458     if (lDetTLBR[0]) AddConstraint(fConstraintXT,0.0); // Top half
459     if (lDetTLBR[1]) AddConstraint(fConstraintXL,0.0); // Left half
460     if (lDetTLBR[2]) AddConstraint(fConstraintXB,0.0); // Bottom half
461     if (lDetTLBR[3]) AddConstraint(fConstraintXR,0.0); // Right half
462   }
463   if (lVarXYT[1]) { // Y constraint
464     if (lDetTLBR[0]) AddConstraint(fConstraintYT,0.0); // Top half
465     if (lDetTLBR[1]) AddConstraint(fConstraintYL,0.0); // Left half
466     if (lDetTLBR[2]) AddConstraint(fConstraintYB,0.0); // Bottom half
467     if (lDetTLBR[3]) AddConstraint(fConstraintYR,0.0); // Right half
468   }
469   if (lVarXYT[2]) { // T constraint
470     if (lDetTLBR[0]) AddConstraint(fConstraintPT,0.0); // Top half
471     if (lDetTLBR[1]) AddConstraint(fConstraintPL,0.0); // Left half
472     if (lDetTLBR[2]) AddConstraint(fConstraintPB,0.0); // Bottom half
473     if (lDetTLBR[3]) AddConstraint(fConstraintPR,0.0); // Right half
474   }
475   if (lVarXYT[3]) { // X-Z constraint
476     if (lDetTLBR[0]) AddConstraint(fConstraintXZT,0.0); // Top half
477     if (lDetTLBR[1]) AddConstraint(fConstraintXZL,0.0); // Left half
478     if (lDetTLBR[2]) AddConstraint(fConstraintXZB,0.0); // Bottom half
479     if (lDetTLBR[3]) AddConstraint(fConstraintXZR,0.0); // Right half
480   }
481   if (lVarXYT[4]) { // Y-Z constraint
482     if (lDetTLBR[0]) AddConstraint(fConstraintYZT,0.0); // Top half
483     if (lDetTLBR[1]) AddConstraint(fConstraintYZL,0.0); // Left half
484     if (lDetTLBR[2]) AddConstraint(fConstraintYZB,0.0); // Bottom half
485     if (lDetTLBR[3]) AddConstraint(fConstraintYZR,0.0); // Right half
486   }
487   if (lVarXYT[5]) { // P-Z constraint
488     if (lDetTLBR[0]) AddConstraint(fConstraintPZT,0.0); // Top half
489     if (lDetTLBR[1]) AddConstraint(fConstraintPZL,0.0); // Left half
490     if (lDetTLBR[2]) AddConstraint(fConstraintPZB,0.0); // Bottom half
491     if (lDetTLBR[3]) AddConstraint(fConstraintPZR,0.0); // Right half
492   }
493   if (lVarXYT[6]) { // X-Y constraint
494     if (lDetTLBR[0]) AddConstraint(fConstraintXYT,0.0); // Top half
495     if (lDetTLBR[1]) AddConstraint(fConstraintXYL,0.0); // Left half
496     if (lDetTLBR[2]) AddConstraint(fConstraintXYB,0.0); // Bottom half
497     if (lDetTLBR[3]) AddConstraint(fConstraintXYR,0.0); // Right half
498   }
499   if (lVarXYT[7]) { // Y-Y constraint
500     if (lDetTLBR[0]) AddConstraint(fConstraintYYT,0.0); // Top half
501     if (lDetTLBR[1]) AddConstraint(fConstraintYYL,0.0); // Left half
502     if (lDetTLBR[2]) AddConstraint(fConstraintYYB,0.0); // Bottom half
503     if (lDetTLBR[3]) AddConstraint(fConstraintYYR,0.0); // Right half
504   }
505   if (lVarXYT[8]) { // P-Y constraint
506     if (lDetTLBR[0]) AddConstraint(fConstraintPYT,0.0); // Top half
507     if (lDetTLBR[1]) AddConstraint(fConstraintPYL,0.0); // Left half
508     if (lDetTLBR[2]) AddConstraint(fConstraintPYB,0.0); // Bottom half
509     if (lDetTLBR[3]) AddConstraint(fConstraintPYR,0.0); // Right half
510   }
511 }
512
513 void AliMUONAlignment::ConstrainT(Int_t lDetElem, Int_t lCh, Double_t *lConstraintT, Int_t iVar, Double_t /*lWeight*/){
514   /// Set constrain equation for top half of spectrometer
515   Int_t lDetElemNumber = (lCh==1) ? lDetElem : lDetElem-fgSNDetElemCh[lCh-2];
516   if (lCh>=1 && lCh<=4){
517     if (lDetElemNumber==0 || lDetElemNumber==1){ // From track crossings
518       lConstraintT[lDetElem*fgNParCh+iVar]=1.0;
519     }
520   }
521   if (lCh>=5 && lCh<=6){
522     if (lDetElemNumber>=0&&lDetElemNumber<=9){
523       lConstraintT[lDetElem*fgNParCh+iVar]=1.0;
524     }
525   }
526   if (lCh>=7 && lCh<=10){
527     if (lDetElemNumber>=0&&lDetElemNumber<=13){
528       lConstraintT[lDetElem*fgNParCh+iVar]=1.0;
529     }
530   }
531 }
532
533 void AliMUONAlignment::ConstrainL(Int_t lDetElem, Int_t lCh, Double_t *lConstraintL, Int_t iVar, Double_t lWeight){
534   /// Set constrain equation for left half of spectrometer
535   Int_t lDetElemNumber = (lCh==1) ? lDetElem : lDetElem-fgSNDetElemCh[lCh-2];
536   if (lCh>=1 && lCh<=4){
537     if (lDetElemNumber==1 || lDetElemNumber==2){ // From track crossings
538       lConstraintL[lDetElem*fgNParCh+iVar]=lWeight;
539     }
540   }
541   if (lCh>=5 && lCh<=6){
542     if (lDetElemNumber>=5&&lDetElemNumber<=13){
543       lConstraintL[lDetElem*fgNParCh+iVar]=lWeight;
544     }
545   }
546   if (lCh>=7 && lCh<=10){
547     if (lDetElemNumber>=7&&lDetElemNumber<=19){
548       lConstraintL[lDetElem*fgNParCh+iVar]=lWeight;
549     }
550   }
551 }
552
553 void AliMUONAlignment::ConstrainB(Int_t lDetElem, Int_t lCh, Double_t *lConstraintB, Int_t iVar, Double_t /*lWeight*/){
554   /// Set constrain equation for bottom half of spectrometer
555   Int_t lDetElemNumber = (lCh==1) ? lDetElem : lDetElem-fgSNDetElemCh[lCh-2];
556   if (lCh>=1 && lCh<=4){
557     if (lDetElemNumber==2 || lDetElemNumber==3){ // From track crossings
558       lConstraintB[lDetElem*fgNParCh+iVar]=1.0;
559     }
560   }
561   if (lCh>=5 && lCh<=6){
562     if ((lDetElemNumber>=9&&lDetElemNumber<=17) || 
563         (lDetElemNumber==0)){
564       lConstraintB[lDetElem*fgNParCh+iVar]=1.0;
565     }
566   }
567   if (lCh>=7 && lCh<=10){
568     if ((lDetElemNumber>=13&&lDetElemNumber<=25) || 
569         (lDetElemNumber==0)){
570       lConstraintB[lDetElem*fgNParCh+iVar]=1.0;
571     }
572   }
573 }
574
575 void AliMUONAlignment::ConstrainR(Int_t lDetElem, Int_t lCh, Double_t *lConstraintR, Int_t iVar, Double_t lWeight){
576   /// Set constrain equation for right half of spectrometer
577   Int_t lDetElemNumber = (lCh==1) ? lDetElem : lDetElem-fgSNDetElemCh[lCh-2];
578   if (lCh>=1 && lCh<=4){
579     if (lDetElemNumber==0 || lDetElemNumber==3){ // From track crossings
580       lConstraintR[lDetElem*fgNParCh+iVar]=lWeight;
581     }
582   }
583   if (lCh>=5 && lCh<=6){
584     if ((lDetElemNumber>=0&&lDetElemNumber<=4) || 
585         (lDetElemNumber>=14&&lDetElemNumber<=17)){
586       lConstraintR[lDetElem*fgNParCh+iVar]=lWeight;
587     }
588   }
589   if (lCh>=7 && lCh<=10){
590     if ((lDetElemNumber>=0&&lDetElemNumber<=6) || 
591         (lDetElemNumber>=20&&lDetElemNumber<=25)){
592       lConstraintR[lDetElem*fgNParCh+iVar]=lWeight;
593     }
594   }
595 }
596
597 void AliMUONAlignment::ResetConstraints(){
598   /// Reset all constraint equations
599   for (Int_t i = 0; i < fgNDetElem; i++){    
600     fConstraintX[i*fgNParCh+0]=0.0;
601     fConstraintX[i*fgNParCh+1]=0.0;
602     fConstraintX[i*fgNParCh+2]=0.0;
603     fConstraintY[i*fgNParCh+0]=0.0;
604     fConstraintY[i*fgNParCh+1]=0.0;
605     fConstraintY[i*fgNParCh+2]=0.0;
606     fConstraintP[i*fgNParCh+0]=0.0;
607     fConstraintP[i*fgNParCh+1]=0.0;
608     fConstraintP[i*fgNParCh+2]=0.0;
609     fConstraintXT[i*fgNParCh+0]=0.0;
610     fConstraintXT[i*fgNParCh+1]=0.0;
611     fConstraintXT[i*fgNParCh+2]=0.0;
612     fConstraintYT[i*fgNParCh+0]=0.0;
613     fConstraintYT[i*fgNParCh+1]=0.0;
614     fConstraintYT[i*fgNParCh+2]=0.0;
615     fConstraintPT[i*fgNParCh+0]=0.0;
616     fConstraintPT[i*fgNParCh+1]=0.0;
617     fConstraintPT[i*fgNParCh+2]=0.0;
618     fConstraintXZT[i*fgNParCh+0]=0.0;
619     fConstraintXZT[i*fgNParCh+1]=0.0;
620     fConstraintXZT[i*fgNParCh+2]=0.0;
621     fConstraintYZT[i*fgNParCh+0]=0.0;
622     fConstraintYZT[i*fgNParCh+1]=0.0;
623     fConstraintYZT[i*fgNParCh+2]=0.0;
624     fConstraintPZT[i*fgNParCh+0]=0.0;
625     fConstraintPZT[i*fgNParCh+1]=0.0;
626     fConstraintPZT[i*fgNParCh+2]=0.0;
627     fConstraintXYT[i*fgNParCh+0]=0.0;
628     fConstraintXYT[i*fgNParCh+1]=0.0;
629     fConstraintXYT[i*fgNParCh+2]=0.0;
630     fConstraintYYT[i*fgNParCh+0]=0.0;
631     fConstraintYYT[i*fgNParCh+1]=0.0;
632     fConstraintYYT[i*fgNParCh+2]=0.0;
633     fConstraintPYT[i*fgNParCh+0]=0.0;
634     fConstraintPYT[i*fgNParCh+1]=0.0;
635     fConstraintPYT[i*fgNParCh+2]=0.0;
636     fConstraintXL[i*fgNParCh+0]=0.0;
637     fConstraintXL[i*fgNParCh+1]=0.0;
638     fConstraintXL[i*fgNParCh+2]=0.0;
639     fConstraintYL[i*fgNParCh+0]=0.0;
640     fConstraintYL[i*fgNParCh+1]=0.0;
641     fConstraintYL[i*fgNParCh+2]=0.0;
642     fConstraintPL[i*fgNParCh+0]=0.0;
643     fConstraintPL[i*fgNParCh+1]=0.0;
644     fConstraintPL[i*fgNParCh+2]=0.0;
645     fConstraintXZL[i*fgNParCh+0]=0.0;
646     fConstraintXZL[i*fgNParCh+1]=0.0;
647     fConstraintXZL[i*fgNParCh+2]=0.0;
648     fConstraintYZL[i*fgNParCh+0]=0.0;
649     fConstraintYZL[i*fgNParCh+1]=0.0;
650     fConstraintYZL[i*fgNParCh+2]=0.0;
651     fConstraintPZL[i*fgNParCh+0]=0.0;
652     fConstraintPZL[i*fgNParCh+1]=0.0;
653     fConstraintPZL[i*fgNParCh+2]=0.0;
654     fConstraintXYL[i*fgNParCh+0]=0.0;
655     fConstraintXYL[i*fgNParCh+1]=0.0;
656     fConstraintXYL[i*fgNParCh+2]=0.0;
657     fConstraintYYL[i*fgNParCh+0]=0.0;
658     fConstraintYYL[i*fgNParCh+1]=0.0;
659     fConstraintYYL[i*fgNParCh+2]=0.0;
660     fConstraintPYL[i*fgNParCh+0]=0.0;
661     fConstraintPYL[i*fgNParCh+1]=0.0;
662     fConstraintPYL[i*fgNParCh+2]=0.0;
663     fConstraintXB[i*fgNParCh+0]=0.0;
664     fConstraintXB[i*fgNParCh+1]=0.0;
665     fConstraintXB[i*fgNParCh+2]=0.0;
666     fConstraintYB[i*fgNParCh+0]=0.0;
667     fConstraintYB[i*fgNParCh+1]=0.0;
668     fConstraintYB[i*fgNParCh+2]=0.0;
669     fConstraintPB[i*fgNParCh+0]=0.0;
670     fConstraintPB[i*fgNParCh+1]=0.0;
671     fConstraintPB[i*fgNParCh+2]=0.0;
672     fConstraintXZB[i*fgNParCh+0]=0.0;
673     fConstraintXZB[i*fgNParCh+1]=0.0;
674     fConstraintXZB[i*fgNParCh+2]=0.0;
675     fConstraintYZB[i*fgNParCh+0]=0.0;
676     fConstraintYZB[i*fgNParCh+1]=0.0;
677     fConstraintYZB[i*fgNParCh+2]=0.0;
678     fConstraintPZB[i*fgNParCh+0]=0.0;
679     fConstraintPZB[i*fgNParCh+1]=0.0;
680     fConstraintPZB[i*fgNParCh+2]=0.0;
681     fConstraintXYB[i*fgNParCh+0]=0.0;
682     fConstraintXYB[i*fgNParCh+1]=0.0;
683     fConstraintXYB[i*fgNParCh+2]=0.0;
684     fConstraintYYB[i*fgNParCh+0]=0.0;
685     fConstraintYYB[i*fgNParCh+1]=0.0;
686     fConstraintYYB[i*fgNParCh+2]=0.0;
687     fConstraintPYB[i*fgNParCh+0]=0.0;
688     fConstraintPYB[i*fgNParCh+1]=0.0;
689     fConstraintPYB[i*fgNParCh+2]=0.0;
690     fConstraintXR[i*fgNParCh+0]=0.0;
691     fConstraintXR[i*fgNParCh+1]=0.0;
692     fConstraintXR[i*fgNParCh+2]=0.0;
693     fConstraintYR[i*fgNParCh+0]=0.0;
694     fConstraintYR[i*fgNParCh+1]=0.0;
695     fConstraintYR[i*fgNParCh+2]=0.0;
696     fConstraintPR[i*fgNParCh+0]=0.0;
697     fConstraintPR[i*fgNParCh+1]=0.0;
698     fConstraintPR[i*fgNParCh+2]=0.0;
699     fConstraintXZR[i*fgNParCh+0]=0.0;
700     fConstraintXZR[i*fgNParCh+1]=0.0;
701     fConstraintXZR[i*fgNParCh+2]=0.0;
702     fConstraintYZR[i*fgNParCh+0]=0.0;
703     fConstraintYZR[i*fgNParCh+1]=0.0;
704     fConstraintYZR[i*fgNParCh+2]=0.0;
705     fConstraintPZR[i*fgNParCh+0]=0.0;
706     fConstraintPZR[i*fgNParCh+1]=0.0;
707     fConstraintPZR[i*fgNParCh+2]=0.0;
708     fConstraintPZR[i*fgNParCh+0]=0.0;
709     fConstraintPZR[i*fgNParCh+1]=0.0;
710     fConstraintPZR[i*fgNParCh+2]=0.0;
711     fConstraintXYR[i*fgNParCh+0]=0.0;
712     fConstraintXYR[i*fgNParCh+1]=0.0;
713     fConstraintXYR[i*fgNParCh+2]=0.0;
714     fConstraintYYR[i*fgNParCh+0]=0.0;
715     fConstraintYYR[i*fgNParCh+1]=0.0;
716     fConstraintYYR[i*fgNParCh+2]=0.0;
717     fConstraintPYR[i*fgNParCh+0]=0.0;
718     fConstraintPYR[i*fgNParCh+1]=0.0;
719     fConstraintPYR[i*fgNParCh+2]=0.0;
720   }
721 }
722
723 void AliMUONAlignment::AddConstraint(Double_t *par, Double_t value) {
724   /// Constrain equation defined by par to value
725   fMillepede->SetGlobalConstraint(par, value);
726   AliInfo("Adding constraint");
727 }
728
729 void AliMUONAlignment::InitGlobalParameters(Double_t *par) {
730   /// Initialize global parameters with par array
731   fMillepede->SetGlobalParameters(par);
732   AliInfo("Init Global Parameters");
733 }
734  
735 void AliMUONAlignment::FixParameter(Int_t iPar, Double_t value) {
736   /// Parameter iPar is encourage to vary in [-value;value]. 
737   /// If value == 0, parameter is fixed
738   fMillepede->SetParSigma(iPar, value);
739   if (value==0) AliInfo(Form("Parameter %i Fixed", iPar));
740 }
741
742 void AliMUONAlignment::ResetLocalEquation()
743 {
744   /// Reset the derivative vectors
745   for(int i=0; i<fNLocal; i++) {
746     fLocalDerivatives[i] = 0.0;
747   }
748   for(int i=0; i<fNGlobal; i++) {
749     fGlobalDerivatives[i] = 0.0;
750   }
751 }
752
753 void AliMUONAlignment::AllowVariations(Bool_t *bChOnOff) {
754   /// Set allowed variation for selected chambers based on fDoF and fAllowVar
755   for (Int_t iCh=1; iCh<=10; iCh++) {
756     if (bChOnOff[iCh-1]) {
757       Int_t iDetElemFirst = (iCh>1) ? fgSNDetElemCh[iCh-2] : 0; 
758       Int_t iDetElemLast = fgSNDetElemCh[iCh-1]; 
759       for (int i=0; i<fgNParCh; i++) {
760         AliDebug(1,Form("fDoF[%d]= %d",i,fDoF[i]));    
761         if (fDoF[i]) {
762           for (Int_t j=iDetElemFirst; j<iDetElemLast; j++){    
763             FixParameter(j*fgNParCh+i, fAllowVar[i]);
764           }
765         }
766       }
767     }
768   }
769 }
770
771 void AliMUONAlignment::SetNonLinear(Int_t iPar  /* set non linear flag */ ) {
772   /// Set nonlinear flag for parameter iPar
773   fMillepede->SetNonLinear(iPar);
774   AliInfo(Form("Parameter %i set to non linear", iPar));
775 }
776
777 void AliMUONAlignment::LocalEquationX() {
778   /// Define local equation for current track and cluster in x coor. measurement
779   // set local derivatives
780   SetLocalDerivative(0, fCosPhi);
781   SetLocalDerivative(1, fCosPhi * (fTrackPos[2] - fTrackPos0[2]));
782   SetLocalDerivative(2, fSinPhi);
783   SetLocalDerivative(3, fSinPhi * (fTrackPos[2] - fTrackPos0[2]));
784
785   // set global derivatives
786   SetGlobalDerivative(fDetElemNumber*fgNParCh+0, -fCosPhi);
787   SetGlobalDerivative(fDetElemNumber*fgNParCh+1, -fSinPhi);
788   if (fBFieldOn){
789     SetGlobalDerivative(fDetElemNumber*fgNParCh+2,
790                         -fSinPhi*(fTrackPos[0]-fDetElemPos[0]) 
791                         +fCosPhi*(fTrackPos[1]-fDetElemPos[1]));
792   }
793   else {
794     SetGlobalDerivative(fDetElemNumber*fgNParCh+2,
795                         -fSinPhi*(fTrackPos0[0]+fTrackSlope0[0]*
796                                   (fTrackPos[2]-fTrackPos0[2])-fDetElemPos[0]) 
797                         +fCosPhi*(fTrackPos0[1]+fTrackSlope0[1]*
798                                   (fTrackPos[2]-fTrackPos0[2])-fDetElemPos[1]));
799   }
800   SetGlobalDerivative(fDetElemNumber*fgNParCh+3, 
801                       fCosPhi*fTrackSlope0[0]+fSinPhi*fTrackSlope0[1]);
802
803   fMillepede->SetLocalEquation(fGlobalDerivatives, fLocalDerivatives, fMeas[0], fSigma[0]);
804 }
805
806 void AliMUONAlignment::LocalEquationY() {
807   /// Define local equation for current track and cluster in y coor. measurement
808   // set local derivatives
809   SetLocalDerivative(0,-fSinPhi);
810   SetLocalDerivative(1,-fSinPhi * (fTrackPos[2] - fTrackPos0[2]));
811   SetLocalDerivative(2, fCosPhi);
812   SetLocalDerivative(3, fCosPhi * (fTrackPos[2] - fTrackPos0[2]));
813
814   // set global derivatives
815   SetGlobalDerivative(fDetElemNumber*fgNParCh+0,  fSinPhi);
816   SetGlobalDerivative(fDetElemNumber*fgNParCh+1, -fCosPhi);
817   if (fBFieldOn){
818   SetGlobalDerivative(fDetElemNumber*fgNParCh+2,
819                       -fCosPhi*(fTrackPos[0]-fDetElemPos[0])
820                       -fSinPhi*(fTrackPos[1]-fDetElemPos[1]));
821   }
822   else {
823     SetGlobalDerivative(fDetElemNumber*fgNParCh+2,
824                         -fCosPhi*(fTrackPos0[0]+fTrackSlope0[0]*
825                                   (fTrackPos[2]-fTrackPos0[2])-fDetElemPos[0])
826                         -fSinPhi*(fTrackPos0[1]+fTrackSlope0[1]*
827                                   (fTrackPos[2]-fTrackPos0[2])-fDetElemPos[1]));
828   }
829   SetGlobalDerivative(fDetElemNumber*fgNParCh+3,
830                       -fSinPhi*fTrackSlope0[0]+fCosPhi*fTrackSlope0[1]);
831
832   fMillepede->SetLocalEquation(fGlobalDerivatives, fLocalDerivatives, fMeas[1], fSigma[1]);
833 }
834
835 void AliMUONAlignment::FillRecPointData() {
836   /// Get information of current cluster
837   fClustPos[0] = fCluster->GetX();
838   fClustPos[1] = fCluster->GetY();
839   fClustPos[2] = fCluster->GetZ();
840   fTransform->Global2Local(fDetElemId,fClustPos[0],fClustPos[1],fClustPos[2],
841                             fClustPosLoc[0],fClustPosLoc[1],fClustPosLoc[2]);
842 }
843
844 void AliMUONAlignment::FillTrackParamData() {
845   /// Get information of current track at current cluster
846   fTrackPos[0] = fTrackParam->GetNonBendingCoor();
847   fTrackPos[1] = fTrackParam->GetBendingCoor();
848   fTrackPos[2] = fTrackParam->GetZ();
849   fTrackSlope[0] = fTrackParam->GetNonBendingSlope();
850   fTrackSlope[1] = fTrackParam->GetBendingSlope();
851   fTransform->Global2Local(fDetElemId,fTrackPos[0],fTrackPos[1],fTrackPos[2],
852                             fTrackPosLoc[0],fTrackPosLoc[1],fTrackPosLoc[2]);
853 }
854
855 void AliMUONAlignment::FillDetElemData() {
856   /// Get information of current detection element
857   Double_t lDetElemLocX = 0.;
858   Double_t lDetElemLocY = 0.;
859   Double_t lDetElemLocZ = 0.;
860   fDetElemId = fCluster->GetDetElemId();
861   fDetElemNumber = fDetElemId%100;
862   for (int iCh=0; iCh<fDetElemId/100-1; iCh++){
863     fDetElemNumber += fgNDetElemCh[iCh];
864   }
865   fTransform->Local2Global(fDetElemId,lDetElemLocX,lDetElemLocY,lDetElemLocZ,
866                            fDetElemPos[0],fDetElemPos[1],fDetElemPos[2]);
867 }
868
869 void AliMUONAlignment::ProcessTrack(AliMUONTrack * track) {
870   /// Process track; Loop over clusters and set local equations
871   fTrack = track;
872   // get tclones arrays.
873   fTrackParamAtCluster = fTrack->GetTrackParamAtCluster();
874   
875   // get size of arrays
876   Int_t nTrackParam = fTrackParamAtCluster->GetEntries();
877   AliDebug(1,Form("Number of track param entries : %i ", nTrackParam));
878
879   for(Int_t iCluster=0; iCluster<nTrackParam; iCluster++) {
880     fTrackParam = (AliMUONTrackParam *) fTrack->GetTrackParamAtCluster()->At(iCluster);
881     fCluster = fTrackParam->GetClusterPtr();
882     if (!fCluster || !fTrackParam) continue;
883     // fill local variables for this position --> one measurement
884     FillDetElemData();
885     FillRecPointData();
886     FillTrackParamData();         
887 //     if (fDetElemId<500) continue;
888     fTrackPos0[0]      = fTrackPos[0];    
889     fTrackPos0[1]      = fTrackPos[1];    
890     fTrackPos0[2]      = fTrackPos[2];    
891     fTrackSlope0[0] = fTrackSlope[0]; 
892     fTrackSlope0[1] = fTrackSlope[1];   
893     break;
894   }
895
896   for(Int_t iCluster=0; iCluster<nTrackParam; iCluster++) {
897     // and get new pointers
898     fTrackParam = (AliMUONTrackParam *) fTrack->GetTrackParamAtCluster()->At(iCluster);
899     fCluster = fTrackParam->GetClusterPtr();
900     if (!fCluster || !fTrackParam) continue;
901     // fill local variables for this position --> one measurement
902     FillDetElemData();        
903     FillRecPointData();
904     FillTrackParamData();
905 //     if (fDetElemId<500) continue;
906     AliDebug(1,Form("cluster: %i", iCluster));
907     AliDebug(1,Form("x: %f\t y: %f\t z: %f\t DetElemID: %i\t ", fClustPos[0], fClustPos[1], fClustPos[2], fDetElemId));
908     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));
909
910     AliDebug(1,Form("Track Parameter: %i", iCluster));
911     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]));
912     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]));
913     
914     fCosPhi = TMath::Cos(fPhi);
915     fSinPhi = TMath::Sin(fPhi);
916     if (fBFieldOn){
917       fMeas[0] = fTrackPos[0] - fClustPos[0];
918       fMeas[1] = fTrackPos[1] - fClustPos[1];
919     }
920     else {
921       fMeas[0] = - fClustPos[0];
922       fMeas[1] = - fClustPos[1];
923     }
924     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]));    
925     // Set local equations
926     LocalEquationX();
927     LocalEquationY();
928   }
929 }
930
931 void AliMUONAlignment::LocalFit(Int_t iTrack, Double_t *lTrackParam, Int_t lSingleFit) {
932   /// Call local fit for this tracks
933   Int_t iRes = fMillepede->LocalFit(iTrack,lTrackParam,lSingleFit);
934   if (iRes && !lSingleFit) {
935     fMillepede->SetNLocalEquations(fMillepede->GetNLocalEquations()+1);
936   }
937 }
938
939 void AliMUONAlignment::GlobalFit(Double_t *parameters,Double_t *errors,Double_t *pulls) {
940   /// Call global fit; Global parameters are stored in parameters
941   fMillepede->GlobalFit(parameters,errors,pulls);
942
943   AliInfo("Done fitting global parameters!");
944   for (int iGlob=0; iGlob<fgNDetElem; iGlob++){
945     printf("%d\t %f\t %f\t %f \n",iGlob,parameters[iGlob*fgNParCh+0],parameters[iGlob*fgNParCh+1],parameters[iGlob*fgNParCh+2]);
946   }
947 }
948
949 Double_t AliMUONAlignment::GetParError(Int_t iPar) {
950   /// Get error of parameter iPar
951   Double_t lErr = fMillepede->GetParError(iPar);
952   return lErr;
953 }
954
955 void AliMUONAlignment::PrintGlobalParameters() {
956   /// Print global parameters
957   fMillepede->PrintGlobalParameters();
958 }
959
960 //_________________________________________________________________________
961 TGeoCombiTrans AliMUONAlignment::ReAlign(const TGeoCombiTrans & transform, double *lMisAlignment) const
962 {
963   /// Realign given transformation by given misalignment and return the misaligned transformation
964   
965   Double_t cartMisAlig[3] = {0,0,0};
966   Double_t angMisAlig[3] = {0,0,0};
967 //   const Double_t *trans = transform.GetTranslation();
968 //   TGeoRotation *rot;
969 //   // check if the rotation we obtain is not NULL
970 //   if (transform.GetRotation()) {    
971 //     rot = transform.GetRotation();
972 //   }
973 //   else {    
974 //     rot = new TGeoRotation("rot");
975 //   }                  // default constructor.
976
977   cartMisAlig[0] = -TMath::Sign(1.0,transform.GetRotationMatrix()[0])*lMisAlignment[0];
978   cartMisAlig[1] = -TMath::Sign(1.0,transform.GetRotationMatrix()[4])*lMisAlignment[1];
979   cartMisAlig[2] = -TMath::Sign(1.0,transform.GetRotationMatrix()[8])*lMisAlignment[3];
980   angMisAlig[2] = -TMath::Sign(1.0,transform.GetRotationMatrix()[0]*transform.GetRotationMatrix()[4])*lMisAlignment[2]*180./TMath::Pi();
981
982   TGeoTranslation deltaTrans(cartMisAlig[0], cartMisAlig[1], cartMisAlig[2]);
983   TGeoRotation deltaRot;
984   deltaRot.RotateX(angMisAlig[0]);
985   deltaRot.RotateY(angMisAlig[1]);
986   deltaRot.RotateZ(angMisAlig[2]);
987
988   TGeoCombiTrans deltaTransf(deltaTrans,deltaRot);
989   TGeoHMatrix newTransfMat = transform * deltaTransf;
990
991   return TGeoCombiTrans(newTransfMat);
992 }
993
994 //______________________________________________________________________
995 AliMUONGeometryTransformer *
996 AliMUONAlignment::ReAlign(const AliMUONGeometryTransformer * transformer,
997                             double *misAlignments, Bool_t verbose)
998                             
999 {
1000   /// Returns a new AliMUONGeometryTransformer with the found misalignments
1001   /// applied. 
1002
1003   // Takes the internal geometry module transformers, copies them
1004   // and gets the Detection Elements from them.
1005   // Takes misalignment parameters and applies these
1006   // to the local transform of the Detection Element
1007   // Obtains the global transform by multiplying the module transformer
1008   // transformation with the local transformation 
1009   // Applies the global transform to a new detection element
1010   // Adds the new detection element to a new module transformer
1011   // Adds the new module transformer to a new geometry transformer
1012   // Returns the new geometry transformer
1013
1014   Double_t lModuleMisAlignment[3] = {0.,0.,0.};
1015   Double_t lDetElemMisAlignment[4] = {0.,0.,0.,0.};
1016   Int_t iDetElemId = 0;
1017   Int_t iDetElemNumber = 0;
1018
1019   AliMUONGeometryTransformer *newGeometryTransformer =
1020     new AliMUONGeometryTransformer();
1021   for (Int_t iMt = 0; iMt < transformer->GetNofModuleTransformers(); iMt++) {
1022     // module transformers    
1023     const AliMUONGeometryModuleTransformer *kModuleTransformer =
1024       transformer->GetModuleTransformer(iMt, true);
1025       
1026     AliMUONGeometryModuleTransformer *newModuleTransformer =
1027       new AliMUONGeometryModuleTransformer(iMt);
1028     newGeometryTransformer->AddModuleTransformer(newModuleTransformer);
1029     
1030     TGeoCombiTrans moduleTransform =
1031       TGeoCombiTrans(*kModuleTransformer->GetTransformation());
1032     // New module transformation
1033     TGeoCombiTrans newModuleTransform = ReAlign(moduleTransform,lModuleMisAlignment);
1034     newModuleTransformer->SetTransformation(newModuleTransform);
1035     
1036     // Get delta transformation: 
1037     // Tdelta = Tnew * Told.inverse
1038     TGeoHMatrix deltaModuleTransform = 
1039       AliMUONGeometryBuilder::Multiply(newModuleTransform, 
1040                                        kModuleTransformer->GetTransformation()->Inverse());    
1041     // Create module mis alignment matrix
1042     newGeometryTransformer
1043       ->AddMisAlignModule(kModuleTransformer->GetModuleId(), deltaModuleTransform);
1044     
1045     AliMpExMap *detElements = kModuleTransformer->GetDetElementStore();
1046     
1047     if (verbose)
1048       AliInfo(Form("%i DEs in old GeometryStore  %i",detElements->GetSize(), iMt));
1049
1050     TIter next(detElements->CreateIterator());
1051     AliMUONGeometryDetElement* detElement;
1052     Int_t iDe(-1);
1053     while ( ( detElement = static_cast<AliMUONGeometryDetElement*>(next()) ) )
1054     {
1055       ++iDe;
1056       // make a new detection element
1057       AliMUONGeometryDetElement *newDetElement =
1058         new AliMUONGeometryDetElement(detElement->GetId(),
1059                                       detElement->GetVolumePath());
1060       TString lDetElemName(detElement->GetDEName());
1061       lDetElemName.ReplaceAll("DE","");
1062       iDetElemId = lDetElemName.Atoi();
1063       iDetElemNumber = iDetElemId%100;
1064       for (int iCh=0; iCh<iDetElemId/100-1; iCh++){
1065         iDetElemNumber += fgNDetElemCh[iCh];
1066       }
1067       for (int i=0; i<fgNParCh; i++) {
1068         lDetElemMisAlignment[i] = 0.0;
1069         if (iMt<fgNTrkMod) {
1070           AliInfo(Form("iMt %i, iCh %i, iDe %i, iDeId %i, iDeNb %i, iPar %i",iMt, iDetElemId/100, iDe, iDetElemId, iDetElemNumber, iDetElemNumber*fgNParCh+i));
1071           lDetElemMisAlignment[i] =  misAlignments[iDetElemNumber*fgNParCh+i];
1072         }
1073       }
1074       // local transformation of this detection element.
1075       TGeoCombiTrans localTransform
1076         = TGeoCombiTrans(*detElement->GetLocalTransformation());
1077       TGeoCombiTrans newLocalTransform = ReAlign(localTransform,lDetElemMisAlignment);
1078       newDetElement->SetLocalTransformation(newLocalTransform);   
1079
1080       // global transformation
1081       TGeoHMatrix newGlobalTransform =
1082         AliMUONGeometryBuilder::Multiply(newModuleTransform,
1083                                          newLocalTransform);
1084       newDetElement->SetGlobalTransformation(newGlobalTransform);
1085           
1086       // add this det element to module
1087       newModuleTransformer->GetDetElementStore()->Add(newDetElement->GetId(),
1088                                                       newDetElement);
1089
1090       // In the Alice Alignment Framework misalignment objects store
1091       // global delta transformation
1092       // Get detection "intermediate" global transformation
1093       TGeoHMatrix newOldGlobalTransform = newModuleTransform * localTransform;
1094       // Get detection element global delta transformation: 
1095       // Tdelta = Tnew * Told.inverse
1096       TGeoHMatrix  deltaGlobalTransform
1097         = AliMUONGeometryBuilder::Multiply(newGlobalTransform, 
1098                                            newOldGlobalTransform.Inverse());
1099           
1100       // Create mis alignment matrix
1101       newGeometryTransformer
1102         ->AddMisAlignDetElement(detElement->GetId(), deltaGlobalTransform);
1103     }
1104       
1105     if (verbose)
1106       AliInfo(Form("Added module transformer %i to the transformer", iMt));
1107     newGeometryTransformer->AddModuleTransformer(newModuleTransformer);
1108   }
1109   return newGeometryTransformer;
1110 }
1111
1112 //______________________________________________________________________
1113 void AliMUONAlignment::SetAlignmentResolution(const TClonesArray* misAlignArray, Int_t rChId, Double_t rChResX, Double_t rChResY, Double_t rDeResX, Double_t rDeResY){
1114   /// Set alignment resolution to misalign objects to be stored in CDB
1115   Int_t chIdMin = (rChId<0)? 0 : rChId;
1116   Int_t chIdMax = (rChId<0)? 9 : rChId;
1117   Double_t chResX = rChResX;
1118   Double_t chResY = rChResY;
1119   Double_t deResX = rDeResX;
1120   Double_t deResY = rDeResY;
1121
1122   TMatrixDSym mChCorrMatrix(6);
1123   mChCorrMatrix[0][0]=chResX*chResX;
1124   mChCorrMatrix[1][1]=chResY*chResY;
1125   //  mChCorrMatrix.Print();
1126
1127   TMatrixDSym mDECorrMatrix(6);
1128   mDECorrMatrix[0][0]=deResX*deResX;
1129   mDECorrMatrix[1][1]=deResY*deResY;
1130   //  mDECorrMatrix.Print();
1131
1132   AliAlignObjMatrix *alignMat = 0x0;
1133
1134   for(Int_t chId=chIdMin; chId<=chIdMax; chId++) {
1135     TString chName1;
1136     TString chName2;
1137     if (chId<4){
1138       chName1 = Form("GM%d",chId);
1139       chName2 = Form("GM%d",chId);
1140     } else {
1141       chName1 = Form("GM%d",4+(chId-4)*2);
1142       chName2 = Form("GM%d",4+(chId-4)*2+1);
1143     }
1144     
1145     for (int i=0; i<misAlignArray->GetEntries(); i++) {
1146       alignMat = (AliAlignObjMatrix*)misAlignArray->At(i);
1147       TString volName(alignMat->GetSymName());
1148       if((volName.Contains(chName1)&&
1149           ((volName.Last('/')==volName.Index(chName1)+chName1.Length())||
1150            (volName.Length()==volName.Index(chName1)+chName1.Length())))||
1151          (volName.Contains(chName2)&&
1152           ((volName.Last('/')==volName.Index(chName2)+chName2.Length())||
1153            (volName.Length()==volName.Index(chName2)+chName2.Length())))){
1154         volName.Remove(0,volName.Last('/')+1);
1155         if (volName.Contains("GM")) {
1156           //    alignMat->Print("NULL");
1157           alignMat->SetCorrMatrix(mChCorrMatrix);
1158         } else if (volName.Contains("DE")) {
1159           //    alignMat->Print("NULL");
1160           alignMat->SetCorrMatrix(mDECorrMatrix);
1161         }
1162       }
1163     }
1164   }
1165 }