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