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