1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 //-----------------------------------------------------------------------------
19 /// \class AliMUONAlignment
20 /// Alignment class for the ALICE DiMuon spectrometer
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.
27 /// \author Bruce Becker, Javier Castillo
28 //-----------------------------------------------------------------------------
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"
40 #include "AliMpExMap.h"
41 #include "AliMpExMapIterator.h"
43 #include "AliAlignObjMatrix.h"
47 #include "TMatrixDSym.h"
50 ClassImp(AliMUONAlignment)
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;
61 AliMUONAlignment::AliMUONAlignment()
68 fTrackParamAtCluster(0),
72 fNGlobal(fgNDetElem*fgNParCh),
82 /// Default constructor, setting define alignment parameters
86 AliInfo(Form("fSigma[0]: %f\t fSigma[1]: %f",fSigma[0],fSigma[1]));
88 fDoF[0] = kTRUE; fDoF[1] = kTRUE; fDoF[2] = kTRUE; fDoF[3] = kTRUE;
89 fAllowVar[0] = 0.05; fAllowVar[1] = 0.05; fAllowVar[2] = 0.001; fAllowVar[3] = 0.5;
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));
93 fMillepede = new AliMillepede();
95 Init(fNGlobal, fNLocal, fNStdDev);
98 AliInfo("Parameters initialized to zero");
102 AliMUONAlignment::~AliMUONAlignment() {
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 */ )
110 /// Initialization of AliMillepede. Fix parameters, define constraints ...
111 fMillepede->InitMille(nGlobal,nLocal,nStdDev,fResCut,fResCutInitial);
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};
117 // AllowVariations(bChOnOff);
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);
125 // FixHalfSpectrometer(bChOnOff,bSpecLROnOff);
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);
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);
140 bVarXYT[0] = kTRUE; bVarXYT[1] = kTRUE; bVarXYT[2] = kFALSE;
141 // AddConstraints(bStOnOff,bVarXYT);
144 if (fStartFac>1) fMillepede->SetIterations(fStartFac);
147 void AliMUONAlignment::FixStation(Int_t iSt){
148 /// Fix all detection elements of station iSt
149 Int_t iDetElemFirst = (iSt>1) ? fgSNDetElemCh[2*(iSt-1)-1] : 0;
150 Int_t iDetElemLast = fgSNDetElemCh[2*(iSt)-1];
151 for (Int_t i = iDetElemFirst; i < iDetElemLast; i++){
152 FixParameter(i*fgNParCh+0, 0.0);
153 FixParameter(i*fgNParCh+1, 0.0);
154 FixParameter(i*fgNParCh+2, 0.0);
155 FixParameter(i*fgNParCh+3, 0.0);
159 void AliMUONAlignment::FixChamber(Int_t iCh){
160 /// Fix all detection elements of chamber iCh
161 Int_t iDetElemFirst = (iCh>1) ? fgSNDetElemCh[iCh-2] : 0;
162 Int_t iDetElemLast = fgSNDetElemCh[iCh-1];
163 for (Int_t i = iDetElemFirst; i < iDetElemLast; i++){
164 FixParameter(i*fgNParCh+0, 0.0);
165 FixParameter(i*fgNParCh+1, 0.0);
166 FixParameter(i*fgNParCh+2, 0.0);
167 FixParameter(i*fgNParCh+3, 0.0);
171 void AliMUONAlignment::FixDetElem(Int_t iDetElemId, TString sVarXYT){
172 /// Fix a given detection element
173 Int_t iDetElemNumber = iDetElemId%100;
174 for (int iCh=0; iCh<iDetElemId/100-1; iCh++){
175 iDetElemNumber += fgNDetElemCh[iCh];
177 if (sVarXYT.Contains("X")) { // X constraint
178 FixParameter(iDetElemNumber*fgNParCh+0, 0.0);
180 if (sVarXYT.Contains("Y")) { // Y constraint
181 FixParameter(iDetElemNumber*fgNParCh+1, 0.0);
183 if (sVarXYT.Contains("T")) { // T constraint
184 FixParameter(iDetElemNumber*fgNParCh+2, 0.0);
186 if (sVarXYT.Contains("Z")) { // T constraint
187 FixParameter(iDetElemNumber*fgNParCh+3, 0.0);
191 void AliMUONAlignment::FixHalfSpectrometer(const Bool_t *lChOnOff, const Bool_t *lSpecLROnOff){
192 /// Fix left or right detector
193 for (Int_t i = 0; i < fgNDetElem; i++){
195 for (iCh=1; iCh<=fgNCh; iCh++){
196 if (i<fgSNDetElemCh[iCh-1]) break;
198 if (lChOnOff[iCh-1]){
199 Int_t lDetElemNumber = (iCh==1) ? i : i-fgSNDetElemCh[iCh-2];
200 if (iCh>=1 && iCh<=4){
201 if ((lDetElemNumber==1 || lDetElemNumber==2) && !lSpecLROnOff[0]){ // From track crossings
202 FixParameter(i*fgNParCh+0, 0.0);
203 FixParameter(i*fgNParCh+1, 0.0);
204 FixParameter(i*fgNParCh+2, 0.0);
205 FixParameter(i*fgNParCh+3, 0.0);
207 if ((lDetElemNumber==0 || lDetElemNumber==3) && !lSpecLROnOff[1]){ // From track crossings
208 FixParameter(i*fgNParCh+0, 0.0);
209 FixParameter(i*fgNParCh+1, 0.0);
210 FixParameter(i*fgNParCh+2, 0.0);
211 FixParameter(i*fgNParCh+3, 0.0);
214 if (iCh>=5 && iCh<=6){
215 if ((lDetElemNumber>=5&&lDetElemNumber<=13) && !lSpecLROnOff[0]){
216 FixParameter(i*fgNParCh+0, 0.0);
217 FixParameter(i*fgNParCh+1, 0.0);
218 FixParameter(i*fgNParCh+2, 0.0);
219 FixParameter(i*fgNParCh+3, 0.0);
221 if (((lDetElemNumber>=0&&lDetElemNumber<=4) ||
222 (lDetElemNumber>=14&&lDetElemNumber<=17)) && !lSpecLROnOff[1]){
223 FixParameter(i*fgNParCh+0, 0.0);
224 FixParameter(i*fgNParCh+1, 0.0);
225 FixParameter(i*fgNParCh+2, 0.0);
226 FixParameter(i*fgNParCh+3, 0.0);
229 if (iCh>=7 && iCh<=10){
230 if ((lDetElemNumber>=7&&lDetElemNumber<=19) && !lSpecLROnOff[0]){
231 FixParameter(i*fgNParCh+0, 0.0);
232 FixParameter(i*fgNParCh+1, 0.0);
233 FixParameter(i*fgNParCh+2, 0.0);
234 FixParameter(i*fgNParCh+3, 0.0);
236 if (((lDetElemNumber>=0&&lDetElemNumber<=6) ||
237 (lDetElemNumber>=20&&lDetElemNumber<=25)) && !lSpecLROnOff[1]){
238 FixParameter(i*fgNParCh+0, 0.0);
239 FixParameter(i*fgNParCh+1, 0.0);
240 FixParameter(i*fgNParCh+2, 0.0);
241 FixParameter(i*fgNParCh+3, 0.0);
248 void AliMUONAlignment::SetNonLinear(const Bool_t *lChOnOff, const Bool_t *lVarXYT){
249 /// Set non linear parameter flag selected chambers and degrees of freedom
250 for (Int_t i = 0; i < fgNDetElem; i++){
252 for (iCh=1; iCh<=fgNCh; iCh++){
253 if (i<fgSNDetElemCh[iCh-1]) break;
255 if (lChOnOff[iCh-1]){
256 if (lVarXYT[0]) { // X constraint
257 SetNonLinear(i*fgNParCh+0);
259 if (lVarXYT[1]) { // Y constraint
260 SetNonLinear(i*fgNParCh+1);
262 if (lVarXYT[2]) { // T constraint
263 SetNonLinear(i*fgNParCh+2);
265 if (lVarXYT[3]) { // Z constraint
266 SetNonLinear(i*fgNParCh+3);
272 void AliMUONAlignment::AddConstraints(const Bool_t *lChOnOff, const Bool_t *lVarXYT){
273 /// Add constraint equations for selected chambers and degrees of freedom
274 for (Int_t i = 0; i < fgNDetElem; i++){
276 for (iCh=1; iCh<=fgNCh; iCh++){
277 if (i<fgSNDetElemCh[iCh-1]) break;
279 if (lChOnOff[iCh-1]){
280 if (lVarXYT[0]) { // X constraint
281 fConstraintX[i*fgNParCh+0]=1.0;
283 if (lVarXYT[1]) { // Y constraint
284 fConstraintY[i*fgNParCh+1]=1.0;
286 if (lVarXYT[2]) { // T constraint
287 fConstraintP[i*fgNParCh+2]=1.0;
289 // if (lVarXYT[3]) { // Z constraint
290 // fConstraintP[i*fgNParCh+3]=1.0;
294 if (lVarXYT[0]) { // X constraint
295 AddConstraint(fConstraintX,0.0);
297 if (lVarXYT[1]) { // Y constraint
298 AddConstraint(fConstraintY,0.0);
300 if (lVarXYT[2]) { // T constraint
301 AddConstraint(fConstraintP,0.0);
303 // if (lVarXYT[3]) { // Z constraint
304 // AddConstraint(fConstraintP,0.0);
308 void AliMUONAlignment::AddConstraints(const Bool_t *lChOnOff, const Bool_t *lVarXYT, const Bool_t *lDetTLBR, const Bool_t *lSpecLROnOff){
309 /// Add constraint equations for selected chambers, degrees of freedom and detector half
310 Double_t lDetElemLocX = 0.;
311 Double_t lDetElemLocY = 0.;
312 Double_t lDetElemLocZ = 0.;
313 Double_t lDetElemGloX = 0.;
314 Double_t lDetElemGloY = 0.;
315 Double_t lDetElemGloZ = 0.;
316 Double_t lMeanY = 0.;
317 Double_t lSigmaY = 0.;
318 Double_t lMeanZ = 0.;
319 Double_t lSigmaZ = 0.;
321 for (Int_t i = 0; i < fgNDetElem; i++){
323 for (iCh=1; iCh<=fgNCh; iCh++){
324 if (i<fgSNDetElemCh[iCh-1]) break;
326 if (lChOnOff[iCh-1]){
327 Int_t lDetElemNumber = (iCh==1) ? i : i-fgSNDetElemCh[iCh-2];
328 Int_t lDetElemId = iCh*100+lDetElemNumber;
329 fTransform->Local2Global(lDetElemId,lDetElemLocX,lDetElemLocY,lDetElemLocZ,
330 lDetElemGloX,lDetElemGloY,lDetElemGloZ);
331 if (iCh>=1 && iCh<=4){
332 if ((lDetElemNumber==1 || lDetElemNumber==2) && lSpecLROnOff[0]){ // From track crossings
333 lMeanY += lDetElemGloY;
334 lSigmaY += lDetElemGloY*lDetElemGloY;
335 lMeanZ += lDetElemGloZ;
336 lSigmaZ += lDetElemGloZ*lDetElemGloZ;
339 if ((lDetElemNumber==0 || lDetElemNumber==3) && lSpecLROnOff[1]){ // From track crossings
340 lMeanY += lDetElemGloY;
341 lSigmaY += lDetElemGloY*lDetElemGloY;
342 lMeanZ += lDetElemGloZ;
343 lSigmaZ += lDetElemGloZ*lDetElemGloZ;
347 if (iCh>=5 && iCh<=6){
348 if ((lDetElemNumber>=5&&lDetElemNumber<=13) && lSpecLROnOff[0]){
349 lMeanY += lDetElemGloY;
350 lSigmaY += lDetElemGloY*lDetElemGloY;
351 lMeanZ += lDetElemGloZ;
352 lSigmaZ += lDetElemGloZ*lDetElemGloZ;
355 if (((lDetElemNumber>=0&&lDetElemNumber<=4) ||
356 (lDetElemNumber>=14&&lDetElemNumber<=17)) && lSpecLROnOff[1]){
357 lMeanY += lDetElemGloY;
358 lSigmaY += lDetElemGloY*lDetElemGloY;
359 lMeanZ += lDetElemGloZ;
360 lSigmaZ += lDetElemGloZ*lDetElemGloZ;
364 if (iCh>=7 && iCh<=10){
365 if ((lDetElemNumber>=7&&lDetElemNumber<=19) && lSpecLROnOff[0]){
366 lMeanY += lDetElemGloY;
367 lSigmaY += lDetElemGloY*lDetElemGloY;
368 lMeanZ += lDetElemGloZ;
369 lSigmaZ += lDetElemGloZ*lDetElemGloZ;
372 if (((lDetElemNumber>=0&&lDetElemNumber<=6) ||
373 (lDetElemNumber>=20&&lDetElemNumber<=25)) && lSpecLROnOff[1]){
374 lMeanY += lDetElemGloY;
375 lSigmaY += lDetElemGloY*lDetElemGloY;
376 lMeanZ += lDetElemGloZ;
377 lSigmaZ += lDetElemGloZ*lDetElemGloZ;
384 lSigmaY /= lNDetElem;
385 lSigmaY = TMath::Sqrt(lSigmaY-lMeanY*lMeanY);
387 lSigmaZ /= lNDetElem;
388 lSigmaZ = TMath::Sqrt(lSigmaZ-lMeanZ*lMeanZ);
389 AliInfo(Form("Used %i DetElem, MeanZ= %f , SigmaZ= %f", lNDetElem,lMeanZ,lSigmaZ));
391 for (Int_t i = 0; i < fgNDetElem; i++){
393 for (iCh=1; iCh<=fgNCh; iCh++){
394 if (i<fgSNDetElemCh[iCh-1]) break;
396 if (lChOnOff[iCh-1]){
397 Int_t lDetElemNumber = (iCh==1) ? i : i-fgSNDetElemCh[iCh-2];
398 Int_t lDetElemId = iCh*100+lDetElemNumber;
399 fTransform->Local2Global(lDetElemId,lDetElemLocX,lDetElemLocY,lDetElemLocZ,
400 lDetElemGloX,lDetElemGloY,lDetElemGloZ);
401 if (lVarXYT[0]) { // X constraint
402 if (lDetTLBR[0]) ConstrainT(i,iCh,fConstraintXT,0); // Top half
403 if (lDetTLBR[1]) ConstrainL(i,iCh,fConstraintXL,0); // Left half
404 if (lDetTLBR[2]) ConstrainB(i,iCh,fConstraintXB,0); // Bottom half
405 if (lDetTLBR[3]) ConstrainR(i,iCh,fConstraintXR,0); // Right half
407 if (lVarXYT[1]) { // Y constraint
408 if (lDetTLBR[0]) ConstrainT(i,iCh,fConstraintYT,1); // Top half
409 if (lDetTLBR[1]) ConstrainL(i,iCh,fConstraintYL,1); // Left half
410 if (lDetTLBR[2]) ConstrainB(i,iCh,fConstraintYB,1); // Bottom half
411 if (lDetTLBR[3]) ConstrainR(i,iCh,fConstraintYR,1); // Right half
413 if (lVarXYT[2]) { // P constraint
414 if (lDetTLBR[0]) ConstrainT(i,iCh,fConstraintPT,2); // Top half
415 if (lDetTLBR[1]) ConstrainL(i,iCh,fConstraintPL,2); // Left half
416 if (lDetTLBR[2]) ConstrainB(i,iCh,fConstraintPB,2); // Bottom half
417 if (lDetTLBR[3]) ConstrainR(i,iCh,fConstraintPR,2); // Right half
419 if (lVarXYT[3]) { // X-Z shearing
420 if (lDetTLBR[0]) ConstrainT(i,iCh,fConstraintXZT,0,(lDetElemGloZ-lMeanZ)/lSigmaZ); // Top half
421 if (lDetTLBR[1]) ConstrainL(i,iCh,fConstraintXZL,0,(lDetElemGloZ-lMeanZ)/lSigmaZ); // Left half
422 if (lDetTLBR[2]) ConstrainB(i,iCh,fConstraintXZB,0,(lDetElemGloZ-lMeanZ)/lSigmaZ); // Bottom half
423 if (lDetTLBR[3]) ConstrainR(i,iCh,fConstraintXZR,0,(lDetElemGloZ-lMeanZ)/lSigmaZ); // Right half
425 if (lVarXYT[4]) { // Y-Z shearing
426 if (lDetTLBR[0]) ConstrainT(i,iCh,fConstraintYZT,1,(lDetElemGloZ-lMeanZ)/lSigmaZ); // Top half
427 if (lDetTLBR[1]) ConstrainL(i,iCh,fConstraintYZL,1,(lDetElemGloZ-lMeanZ)/lSigmaZ); // Left half
428 if (lDetTLBR[2]) ConstrainB(i,iCh,fConstraintYZB,1,(lDetElemGloZ-lMeanZ)/lSigmaZ); // Bottom half
429 if (lDetTLBR[3]) ConstrainR(i,iCh,fConstraintYZR,1,(lDetElemGloZ-lMeanZ)/lSigmaZ); // Right half
431 if (lVarXYT[5]) { // P-Z rotation
432 if (lDetTLBR[0]) ConstrainT(i,iCh,fConstraintPZT,2,(lDetElemGloZ-lMeanZ)/lSigmaZ); // Top half
433 if (lDetTLBR[1]) ConstrainL(i,iCh,fConstraintPZL,2,(lDetElemGloZ-lMeanZ)/lSigmaZ); // Left half
434 if (lDetTLBR[2]) ConstrainB(i,iCh,fConstraintPZB,2,(lDetElemGloZ-lMeanZ)/lSigmaZ); // Bottom half
435 if (lDetTLBR[3]) ConstrainR(i,iCh,fConstraintPZR,2,(lDetElemGloZ-lMeanZ)/lSigmaZ); // Right half
437 if (lVarXYT[6]) { // X-Y shearing
438 if (lDetTLBR[0]) ConstrainT(i,iCh,fConstraintXYT,0,(lDetElemGloY-lMeanY)/lSigmaY); // Top half
439 if (lDetTLBR[1]) ConstrainL(i,iCh,fConstraintXYL,0,(lDetElemGloY-lMeanY)/lSigmaY); // Left half
440 if (lDetTLBR[2]) ConstrainB(i,iCh,fConstraintXYB,0,(lDetElemGloY-lMeanY)/lSigmaY); // Bottom half
441 if (lDetTLBR[3]) ConstrainR(i,iCh,fConstraintXYR,0,(lDetElemGloY-lMeanY)/lSigmaY); // Right half
443 if (lVarXYT[7]) { // Y-Y scaling
444 if (lDetTLBR[0]) ConstrainT(i,iCh,fConstraintYYT,1,(lDetElemGloY-lMeanY)/lSigmaY); // Top half
445 if (lDetTLBR[1]) ConstrainL(i,iCh,fConstraintYYL,1,(lDetElemGloY-lMeanY)/lSigmaY); // Left half
446 if (lDetTLBR[2]) ConstrainB(i,iCh,fConstraintYYB,1,(lDetElemGloY-lMeanY)/lSigmaY); // Bottom half
447 if (lDetTLBR[3]) ConstrainR(i,iCh,fConstraintYYR,1,(lDetElemGloY-lMeanY)/lSigmaY); // Right half
449 if (lVarXYT[8]) { // P-Y rotation
450 if (lDetTLBR[0]) ConstrainT(i,iCh,fConstraintPYT,2,(lDetElemGloY-lMeanY)/lSigmaY); // Top half
451 if (lDetTLBR[1]) ConstrainL(i,iCh,fConstraintPYL,2,(lDetElemGloY-lMeanY)/lSigmaY); // Left half
452 if (lDetTLBR[2]) ConstrainB(i,iCh,fConstraintPYB,2,(lDetElemGloY-lMeanY)/lSigmaY); // Bottom half
453 if (lDetTLBR[3]) ConstrainR(i,iCh,fConstraintPYR,2,(lDetElemGloY-lMeanY)/lSigmaY); // Right half
457 if (lVarXYT[0]) { // X constraint
458 if (lDetTLBR[0]) AddConstraint(fConstraintXT,0.0); // Top half
459 if (lDetTLBR[1]) AddConstraint(fConstraintXL,0.0); // Left half
460 if (lDetTLBR[2]) AddConstraint(fConstraintXB,0.0); // Bottom half
461 if (lDetTLBR[3]) AddConstraint(fConstraintXR,0.0); // Right half
463 if (lVarXYT[1]) { // Y constraint
464 if (lDetTLBR[0]) AddConstraint(fConstraintYT,0.0); // Top half
465 if (lDetTLBR[1]) AddConstraint(fConstraintYL,0.0); // Left half
466 if (lDetTLBR[2]) AddConstraint(fConstraintYB,0.0); // Bottom half
467 if (lDetTLBR[3]) AddConstraint(fConstraintYR,0.0); // Right half
469 if (lVarXYT[2]) { // T constraint
470 if (lDetTLBR[0]) AddConstraint(fConstraintPT,0.0); // Top half
471 if (lDetTLBR[1]) AddConstraint(fConstraintPL,0.0); // Left half
472 if (lDetTLBR[2]) AddConstraint(fConstraintPB,0.0); // Bottom half
473 if (lDetTLBR[3]) AddConstraint(fConstraintPR,0.0); // Right half
475 if (lVarXYT[3]) { // X-Z constraint
476 if (lDetTLBR[0]) AddConstraint(fConstraintXZT,0.0); // Top half
477 if (lDetTLBR[1]) AddConstraint(fConstraintXZL,0.0); // Left half
478 if (lDetTLBR[2]) AddConstraint(fConstraintXZB,0.0); // Bottom half
479 if (lDetTLBR[3]) AddConstraint(fConstraintXZR,0.0); // Right half
481 if (lVarXYT[4]) { // Y-Z constraint
482 if (lDetTLBR[0]) AddConstraint(fConstraintYZT,0.0); // Top half
483 if (lDetTLBR[1]) AddConstraint(fConstraintYZL,0.0); // Left half
484 if (lDetTLBR[2]) AddConstraint(fConstraintYZB,0.0); // Bottom half
485 if (lDetTLBR[3]) AddConstraint(fConstraintYZR,0.0); // Right half
487 if (lVarXYT[5]) { // P-Z constraint
488 if (lDetTLBR[0]) AddConstraint(fConstraintPZT,0.0); // Top half
489 if (lDetTLBR[1]) AddConstraint(fConstraintPZL,0.0); // Left half
490 if (lDetTLBR[2]) AddConstraint(fConstraintPZB,0.0); // Bottom half
491 if (lDetTLBR[3]) AddConstraint(fConstraintPZR,0.0); // Right half
493 if (lVarXYT[6]) { // X-Y constraint
494 if (lDetTLBR[0]) AddConstraint(fConstraintXYT,0.0); // Top half
495 if (lDetTLBR[1]) AddConstraint(fConstraintXYL,0.0); // Left half
496 if (lDetTLBR[2]) AddConstraint(fConstraintXYB,0.0); // Bottom half
497 if (lDetTLBR[3]) AddConstraint(fConstraintXYR,0.0); // Right half
499 if (lVarXYT[7]) { // Y-Y constraint
500 if (lDetTLBR[0]) AddConstraint(fConstraintYYT,0.0); // Top half
501 if (lDetTLBR[1]) AddConstraint(fConstraintYYL,0.0); // Left half
502 if (lDetTLBR[2]) AddConstraint(fConstraintYYB,0.0); // Bottom half
503 if (lDetTLBR[3]) AddConstraint(fConstraintYYR,0.0); // Right half
505 if (lVarXYT[8]) { // P-Y constraint
506 if (lDetTLBR[0]) AddConstraint(fConstraintPYT,0.0); // Top half
507 if (lDetTLBR[1]) AddConstraint(fConstraintPYL,0.0); // Left half
508 if (lDetTLBR[2]) AddConstraint(fConstraintPYB,0.0); // Bottom half
509 if (lDetTLBR[3]) AddConstraint(fConstraintPYR,0.0); // Right half
513 void AliMUONAlignment::ConstrainT(Int_t lDetElem, Int_t lCh, Double_t *lConstraintT, Int_t iVar, Double_t /*lWeight*/) const{
514 /// Set constrain equation for top half of spectrometer
515 Int_t lDetElemNumber = (lCh==1) ? lDetElem : lDetElem-fgSNDetElemCh[lCh-2];
516 if (lCh>=1 && lCh<=4){
517 if (lDetElemNumber==0 || lDetElemNumber==1){ // From track crossings
518 lConstraintT[lDetElem*fgNParCh+iVar]=1.0;
521 if (lCh>=5 && lCh<=6){
522 if (lDetElemNumber>=0&&lDetElemNumber<=9){
523 lConstraintT[lDetElem*fgNParCh+iVar]=1.0;
526 if (lCh>=7 && lCh<=10){
527 if (lDetElemNumber>=0&&lDetElemNumber<=13){
528 lConstraintT[lDetElem*fgNParCh+iVar]=1.0;
533 void AliMUONAlignment::ConstrainL(Int_t lDetElem, Int_t lCh, Double_t *lConstraintL, Int_t iVar, Double_t lWeight) const{
534 /// Set constrain equation for left half of spectrometer
535 Int_t lDetElemNumber = (lCh==1) ? lDetElem : lDetElem-fgSNDetElemCh[lCh-2];
536 if (lCh>=1 && lCh<=4){
537 if (lDetElemNumber==1 || lDetElemNumber==2){ // From track crossings
538 lConstraintL[lDetElem*fgNParCh+iVar]=lWeight;
541 if (lCh>=5 && lCh<=6){
542 if (lDetElemNumber>=5&&lDetElemNumber<=13){
543 lConstraintL[lDetElem*fgNParCh+iVar]=lWeight;
546 if (lCh>=7 && lCh<=10){
547 if (lDetElemNumber>=7&&lDetElemNumber<=19){
548 lConstraintL[lDetElem*fgNParCh+iVar]=lWeight;
553 void AliMUONAlignment::ConstrainB(Int_t lDetElem, Int_t lCh, Double_t *lConstraintB, Int_t iVar, Double_t /*lWeight*/) const{
554 /// Set constrain equation for bottom half of spectrometer
555 Int_t lDetElemNumber = (lCh==1) ? lDetElem : lDetElem-fgSNDetElemCh[lCh-2];
556 if (lCh>=1 && lCh<=4){
557 if (lDetElemNumber==2 || lDetElemNumber==3){ // From track crossings
558 lConstraintB[lDetElem*fgNParCh+iVar]=1.0;
561 if (lCh>=5 && lCh<=6){
562 if ((lDetElemNumber>=9&&lDetElemNumber<=17) ||
563 (lDetElemNumber==0)){
564 lConstraintB[lDetElem*fgNParCh+iVar]=1.0;
567 if (lCh>=7 && lCh<=10){
568 if ((lDetElemNumber>=13&&lDetElemNumber<=25) ||
569 (lDetElemNumber==0)){
570 lConstraintB[lDetElem*fgNParCh+iVar]=1.0;
575 void AliMUONAlignment::ConstrainR(Int_t lDetElem, Int_t lCh, Double_t *lConstraintR, Int_t iVar, Double_t lWeight) const{
576 /// Set constrain equation for right half of spectrometer
577 Int_t lDetElemNumber = (lCh==1) ? lDetElem : lDetElem-fgSNDetElemCh[lCh-2];
578 if (lCh>=1 && lCh<=4){
579 if (lDetElemNumber==0 || lDetElemNumber==3){ // From track crossings
580 lConstraintR[lDetElem*fgNParCh+iVar]=lWeight;
583 if (lCh>=5 && lCh<=6){
584 if ((lDetElemNumber>=0&&lDetElemNumber<=4) ||
585 (lDetElemNumber>=14&&lDetElemNumber<=17)){
586 lConstraintR[lDetElem*fgNParCh+iVar]=lWeight;
589 if (lCh>=7 && lCh<=10){
590 if ((lDetElemNumber>=0&&lDetElemNumber<=6) ||
591 (lDetElemNumber>=20&&lDetElemNumber<=25)){
592 lConstraintR[lDetElem*fgNParCh+iVar]=lWeight;
597 void AliMUONAlignment::ResetConstraints(){
598 /// Reset all constraint equations
599 for (Int_t i = 0; i < fgNDetElem; i++){
600 fConstraintX[i*fgNParCh+0]=0.0;
601 fConstraintX[i*fgNParCh+1]=0.0;
602 fConstraintX[i*fgNParCh+2]=0.0;
603 fConstraintY[i*fgNParCh+0]=0.0;
604 fConstraintY[i*fgNParCh+1]=0.0;
605 fConstraintY[i*fgNParCh+2]=0.0;
606 fConstraintP[i*fgNParCh+0]=0.0;
607 fConstraintP[i*fgNParCh+1]=0.0;
608 fConstraintP[i*fgNParCh+2]=0.0;
609 fConstraintXT[i*fgNParCh+0]=0.0;
610 fConstraintXT[i*fgNParCh+1]=0.0;
611 fConstraintXT[i*fgNParCh+2]=0.0;
612 fConstraintYT[i*fgNParCh+0]=0.0;
613 fConstraintYT[i*fgNParCh+1]=0.0;
614 fConstraintYT[i*fgNParCh+2]=0.0;
615 fConstraintPT[i*fgNParCh+0]=0.0;
616 fConstraintPT[i*fgNParCh+1]=0.0;
617 fConstraintPT[i*fgNParCh+2]=0.0;
618 fConstraintXZT[i*fgNParCh+0]=0.0;
619 fConstraintXZT[i*fgNParCh+1]=0.0;
620 fConstraintXZT[i*fgNParCh+2]=0.0;
621 fConstraintYZT[i*fgNParCh+0]=0.0;
622 fConstraintYZT[i*fgNParCh+1]=0.0;
623 fConstraintYZT[i*fgNParCh+2]=0.0;
624 fConstraintPZT[i*fgNParCh+0]=0.0;
625 fConstraintPZT[i*fgNParCh+1]=0.0;
626 fConstraintPZT[i*fgNParCh+2]=0.0;
627 fConstraintXYT[i*fgNParCh+0]=0.0;
628 fConstraintXYT[i*fgNParCh+1]=0.0;
629 fConstraintXYT[i*fgNParCh+2]=0.0;
630 fConstraintYYT[i*fgNParCh+0]=0.0;
631 fConstraintYYT[i*fgNParCh+1]=0.0;
632 fConstraintYYT[i*fgNParCh+2]=0.0;
633 fConstraintPYT[i*fgNParCh+0]=0.0;
634 fConstraintPYT[i*fgNParCh+1]=0.0;
635 fConstraintPYT[i*fgNParCh+2]=0.0;
636 fConstraintXL[i*fgNParCh+0]=0.0;
637 fConstraintXL[i*fgNParCh+1]=0.0;
638 fConstraintXL[i*fgNParCh+2]=0.0;
639 fConstraintYL[i*fgNParCh+0]=0.0;
640 fConstraintYL[i*fgNParCh+1]=0.0;
641 fConstraintYL[i*fgNParCh+2]=0.0;
642 fConstraintPL[i*fgNParCh+0]=0.0;
643 fConstraintPL[i*fgNParCh+1]=0.0;
644 fConstraintPL[i*fgNParCh+2]=0.0;
645 fConstraintXZL[i*fgNParCh+0]=0.0;
646 fConstraintXZL[i*fgNParCh+1]=0.0;
647 fConstraintXZL[i*fgNParCh+2]=0.0;
648 fConstraintYZL[i*fgNParCh+0]=0.0;
649 fConstraintYZL[i*fgNParCh+1]=0.0;
650 fConstraintYZL[i*fgNParCh+2]=0.0;
651 fConstraintPZL[i*fgNParCh+0]=0.0;
652 fConstraintPZL[i*fgNParCh+1]=0.0;
653 fConstraintPZL[i*fgNParCh+2]=0.0;
654 fConstraintXYL[i*fgNParCh+0]=0.0;
655 fConstraintXYL[i*fgNParCh+1]=0.0;
656 fConstraintXYL[i*fgNParCh+2]=0.0;
657 fConstraintYYL[i*fgNParCh+0]=0.0;
658 fConstraintYYL[i*fgNParCh+1]=0.0;
659 fConstraintYYL[i*fgNParCh+2]=0.0;
660 fConstraintPYL[i*fgNParCh+0]=0.0;
661 fConstraintPYL[i*fgNParCh+1]=0.0;
662 fConstraintPYL[i*fgNParCh+2]=0.0;
663 fConstraintXB[i*fgNParCh+0]=0.0;
664 fConstraintXB[i*fgNParCh+1]=0.0;
665 fConstraintXB[i*fgNParCh+2]=0.0;
666 fConstraintYB[i*fgNParCh+0]=0.0;
667 fConstraintYB[i*fgNParCh+1]=0.0;
668 fConstraintYB[i*fgNParCh+2]=0.0;
669 fConstraintPB[i*fgNParCh+0]=0.0;
670 fConstraintPB[i*fgNParCh+1]=0.0;
671 fConstraintPB[i*fgNParCh+2]=0.0;
672 fConstraintXZB[i*fgNParCh+0]=0.0;
673 fConstraintXZB[i*fgNParCh+1]=0.0;
674 fConstraintXZB[i*fgNParCh+2]=0.0;
675 fConstraintYZB[i*fgNParCh+0]=0.0;
676 fConstraintYZB[i*fgNParCh+1]=0.0;
677 fConstraintYZB[i*fgNParCh+2]=0.0;
678 fConstraintPZB[i*fgNParCh+0]=0.0;
679 fConstraintPZB[i*fgNParCh+1]=0.0;
680 fConstraintPZB[i*fgNParCh+2]=0.0;
681 fConstraintXYB[i*fgNParCh+0]=0.0;
682 fConstraintXYB[i*fgNParCh+1]=0.0;
683 fConstraintXYB[i*fgNParCh+2]=0.0;
684 fConstraintYYB[i*fgNParCh+0]=0.0;
685 fConstraintYYB[i*fgNParCh+1]=0.0;
686 fConstraintYYB[i*fgNParCh+2]=0.0;
687 fConstraintPYB[i*fgNParCh+0]=0.0;
688 fConstraintPYB[i*fgNParCh+1]=0.0;
689 fConstraintPYB[i*fgNParCh+2]=0.0;
690 fConstraintXR[i*fgNParCh+0]=0.0;
691 fConstraintXR[i*fgNParCh+1]=0.0;
692 fConstraintXR[i*fgNParCh+2]=0.0;
693 fConstraintYR[i*fgNParCh+0]=0.0;
694 fConstraintYR[i*fgNParCh+1]=0.0;
695 fConstraintYR[i*fgNParCh+2]=0.0;
696 fConstraintPR[i*fgNParCh+0]=0.0;
697 fConstraintPR[i*fgNParCh+1]=0.0;
698 fConstraintPR[i*fgNParCh+2]=0.0;
699 fConstraintXZR[i*fgNParCh+0]=0.0;
700 fConstraintXZR[i*fgNParCh+1]=0.0;
701 fConstraintXZR[i*fgNParCh+2]=0.0;
702 fConstraintYZR[i*fgNParCh+0]=0.0;
703 fConstraintYZR[i*fgNParCh+1]=0.0;
704 fConstraintYZR[i*fgNParCh+2]=0.0;
705 fConstraintPZR[i*fgNParCh+0]=0.0;
706 fConstraintPZR[i*fgNParCh+1]=0.0;
707 fConstraintPZR[i*fgNParCh+2]=0.0;
708 fConstraintPZR[i*fgNParCh+0]=0.0;
709 fConstraintPZR[i*fgNParCh+1]=0.0;
710 fConstraintPZR[i*fgNParCh+2]=0.0;
711 fConstraintXYR[i*fgNParCh+0]=0.0;
712 fConstraintXYR[i*fgNParCh+1]=0.0;
713 fConstraintXYR[i*fgNParCh+2]=0.0;
714 fConstraintYYR[i*fgNParCh+0]=0.0;
715 fConstraintYYR[i*fgNParCh+1]=0.0;
716 fConstraintYYR[i*fgNParCh+2]=0.0;
717 fConstraintPYR[i*fgNParCh+0]=0.0;
718 fConstraintPYR[i*fgNParCh+1]=0.0;
719 fConstraintPYR[i*fgNParCh+2]=0.0;
723 void AliMUONAlignment::AddConstraint(Double_t *par, Double_t value) {
724 /// Constrain equation defined by par to value
725 fMillepede->SetGlobalConstraint(par, value);
726 AliInfo("Adding constraint");
729 void AliMUONAlignment::InitGlobalParameters(Double_t *par) {
730 /// Initialize global parameters with par array
731 fMillepede->SetGlobalParameters(par);
732 AliInfo("Init Global Parameters");
735 void AliMUONAlignment::FixParameter(Int_t iPar, Double_t value) {
736 /// Parameter iPar is encourage to vary in [-value;value].
737 /// If value == 0, parameter is fixed
738 fMillepede->SetParSigma(iPar, value);
739 if (TMath::Abs(value)<1e-4) AliInfo(Form("Parameter %i Fixed", iPar));
742 void AliMUONAlignment::ResetLocalEquation()
744 /// Reset the derivative vectors
745 for(int i=0; i<fNLocal; i++) {
746 fLocalDerivatives[i] = 0.0;
748 for(int i=0; i<fNGlobal; i++) {
749 fGlobalDerivatives[i] = 0.0;
753 void AliMUONAlignment::AllowVariations(const Bool_t *bChOnOff) {
754 /// Set allowed variation for selected chambers based on fDoF and fAllowVar
755 for (Int_t iCh=1; iCh<=10; iCh++) {
756 if (bChOnOff[iCh-1]) {
757 Int_t iDetElemFirst = (iCh>1) ? fgSNDetElemCh[iCh-2] : 0;
758 Int_t iDetElemLast = fgSNDetElemCh[iCh-1];
759 for (int i=0; i<fgNParCh; i++) {
760 AliDebug(1,Form("fDoF[%d]= %d",i,fDoF[i]));
762 for (Int_t j=iDetElemFirst; j<iDetElemLast; j++){
763 FixParameter(j*fgNParCh+i, fAllowVar[i]);
771 void AliMUONAlignment::SetNonLinear(Int_t iPar /* set non linear flag */ ) {
772 /// Set nonlinear flag for parameter iPar
773 fMillepede->SetNonLinear(iPar);
774 AliInfo(Form("Parameter %i set to non linear", iPar));
778 void AliMUONAlignment::SetSigmaXY(Double_t sigmaX, Double_t sigmaY) {
779 /// Set expected measurement resolution
780 fSigma[0] = sigmaX; fSigma[1] = sigmaY;
781 AliInfo(Form("Using fSigma[0]=%f and fSigma[1]=%f",fSigma[0],fSigma[1]));
785 void AliMUONAlignment::LocalEquationX() {
786 /// Define local equation for current track and cluster in x coor. measurement
787 // set local derivatives
788 SetLocalDerivative(0, fCosPhi);
789 SetLocalDerivative(1, fCosPhi * (fTrackPos[2] - fTrackPos0[2]));
790 SetLocalDerivative(2, fSinPhi);
791 SetLocalDerivative(3, fSinPhi * (fTrackPos[2] - fTrackPos0[2]));
793 // set global derivatives
794 SetGlobalDerivative(fDetElemNumber*fgNParCh+0, -fCosPhi);
795 SetGlobalDerivative(fDetElemNumber*fgNParCh+1, -fSinPhi);
797 SetGlobalDerivative(fDetElemNumber*fgNParCh+2,
798 -fSinPhi*(fTrackPos[0]-fDetElemPos[0])
799 +fCosPhi*(fTrackPos[1]-fDetElemPos[1]));
802 SetGlobalDerivative(fDetElemNumber*fgNParCh+2,
803 -fSinPhi*(fTrackPos0[0]+fTrackSlope0[0]*
804 (fTrackPos[2]-fTrackPos0[2])-fDetElemPos[0])
805 +fCosPhi*(fTrackPos0[1]+fTrackSlope0[1]*
806 (fTrackPos[2]-fTrackPos0[2])-fDetElemPos[1]));
808 SetGlobalDerivative(fDetElemNumber*fgNParCh+3,
809 fCosPhi*fTrackSlope0[0]+fSinPhi*fTrackSlope0[1]);
811 fMillepede->SetLocalEquation(fGlobalDerivatives, fLocalDerivatives, fMeas[0], fSigma[0]);
814 void AliMUONAlignment::LocalEquationY() {
815 /// Define local equation for current track and cluster in y coor. measurement
816 // set local derivatives
817 SetLocalDerivative(0,-fSinPhi);
818 SetLocalDerivative(1,-fSinPhi * (fTrackPos[2] - fTrackPos0[2]));
819 SetLocalDerivative(2, fCosPhi);
820 SetLocalDerivative(3, fCosPhi * (fTrackPos[2] - fTrackPos0[2]));
822 // set global derivatives
823 SetGlobalDerivative(fDetElemNumber*fgNParCh+0, fSinPhi);
824 SetGlobalDerivative(fDetElemNumber*fgNParCh+1, -fCosPhi);
826 SetGlobalDerivative(fDetElemNumber*fgNParCh+2,
827 -fCosPhi*(fTrackPos[0]-fDetElemPos[0])
828 -fSinPhi*(fTrackPos[1]-fDetElemPos[1]));
831 SetGlobalDerivative(fDetElemNumber*fgNParCh+2,
832 -fCosPhi*(fTrackPos0[0]+fTrackSlope0[0]*
833 (fTrackPos[2]-fTrackPos0[2])-fDetElemPos[0])
834 -fSinPhi*(fTrackPos0[1]+fTrackSlope0[1]*
835 (fTrackPos[2]-fTrackPos0[2])-fDetElemPos[1]));
837 SetGlobalDerivative(fDetElemNumber*fgNParCh+3,
838 -fSinPhi*fTrackSlope0[0]+fCosPhi*fTrackSlope0[1]);
840 fMillepede->SetLocalEquation(fGlobalDerivatives, fLocalDerivatives, fMeas[1], fSigma[1]);
843 void AliMUONAlignment::FillRecPointData() {
844 /// Get information of current cluster
845 fClustPos[0] = fCluster->GetX();
846 fClustPos[1] = fCluster->GetY();
847 fClustPos[2] = fCluster->GetZ();
848 fTransform->Global2Local(fDetElemId,fClustPos[0],fClustPos[1],fClustPos[2],
849 fClustPosLoc[0],fClustPosLoc[1],fClustPosLoc[2]);
852 void AliMUONAlignment::FillTrackParamData() {
853 /// Get information of current track at current cluster
854 fTrackPos[0] = fTrackParam->GetNonBendingCoor();
855 fTrackPos[1] = fTrackParam->GetBendingCoor();
856 fTrackPos[2] = fTrackParam->GetZ();
857 fTrackSlope[0] = fTrackParam->GetNonBendingSlope();
858 fTrackSlope[1] = fTrackParam->GetBendingSlope();
859 fTransform->Global2Local(fDetElemId,fTrackPos[0],fTrackPos[1],fTrackPos[2],
860 fTrackPosLoc[0],fTrackPosLoc[1],fTrackPosLoc[2]);
863 void AliMUONAlignment::FillDetElemData() {
864 /// Get information of current detection element
865 Double_t lDetElemLocX = 0.;
866 Double_t lDetElemLocY = 0.;
867 Double_t lDetElemLocZ = 0.;
868 fDetElemId = fCluster->GetDetElemId();
869 fDetElemNumber = fDetElemId%100;
870 for (int iCh=0; iCh<fDetElemId/100-1; iCh++){
871 fDetElemNumber += fgNDetElemCh[iCh];
873 fTransform->Local2Global(fDetElemId,lDetElemLocX,lDetElemLocY,lDetElemLocZ,
874 fDetElemPos[0],fDetElemPos[1],fDetElemPos[2]);
877 void AliMUONAlignment::ProcessTrack(AliMUONTrack * track) {
878 /// Process track; Loop over clusters and set local equations
880 // get tclones arrays.
881 fTrackParamAtCluster = fTrack->GetTrackParamAtCluster();
883 // get size of arrays
884 Int_t nTrackParam = fTrackParamAtCluster->GetEntries();
885 AliDebug(1,Form("Number of track param entries : %i ", nTrackParam));
887 for(Int_t iCluster=0; iCluster<nTrackParam; iCluster++) {
888 fTrackParam = (AliMUONTrackParam *) fTrack->GetTrackParamAtCluster()->At(iCluster);
889 fCluster = fTrackParam->GetClusterPtr();
890 if (!fCluster || !fTrackParam) continue;
891 // fill local variables for this position --> one measurement
894 FillTrackParamData();
895 // if (fDetElemId<500) continue;
896 fTrackPos0[0] = fTrackPos[0];
897 fTrackPos0[1] = fTrackPos[1];
898 fTrackPos0[2] = fTrackPos[2];
899 fTrackSlope0[0] = fTrackSlope[0];
900 fTrackSlope0[1] = fTrackSlope[1];
904 for(Int_t iCluster=0; iCluster<nTrackParam; iCluster++) {
905 // and get new pointers
906 fTrackParam = (AliMUONTrackParam *) fTrack->GetTrackParamAtCluster()->At(iCluster);
907 fCluster = fTrackParam->GetClusterPtr();
908 if (!fCluster || !fTrackParam) continue;
909 // fill local variables for this position --> one measurement
912 FillTrackParamData();
913 // if (fDetElemId<500) continue;
914 AliDebug(1,Form("cluster: %i", iCluster));
915 AliDebug(1,Form("x: %f\t y: %f\t z: %f\t DetElemID: %i\t ", fClustPos[0], fClustPos[1], fClustPos[2], fDetElemId));
916 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));
918 AliDebug(1,Form("Track Parameter: %i", iCluster));
919 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]));
920 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]));
922 fCosPhi = TMath::Cos(fPhi);
923 fSinPhi = TMath::Sin(fPhi);
925 fMeas[0] = fTrackPos[0] - fClustPos[0];
926 fMeas[1] = fTrackPos[1] - fClustPos[1];
929 fMeas[0] = - fClustPos[0];
930 fMeas[1] = - fClustPos[1];
932 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]));
933 // Set local equations
939 void AliMUONAlignment::LocalFit(Int_t iTrack, Double_t *lTrackParam, Int_t lSingleFit) {
940 /// Call local fit for this tracks
941 Int_t iRes = fMillepede->LocalFit(iTrack,lTrackParam,lSingleFit);
942 if (iRes && !lSingleFit) {
943 fMillepede->SetNLocalEquations(fMillepede->GetNLocalEquations()+1);
947 void AliMUONAlignment::GlobalFit(Double_t *parameters,Double_t *errors,Double_t *pulls) {
948 /// Call global fit; Global parameters are stored in parameters
949 fMillepede->GlobalFit(parameters,errors,pulls);
951 AliInfo("Done fitting global parameters!");
952 for (int iGlob=0; iGlob<fgNDetElem; iGlob++){
953 printf("%d\t %f\t %f\t %f \n",iGlob,parameters[iGlob*fgNParCh+0],parameters[iGlob*fgNParCh+1],parameters[iGlob*fgNParCh+2]);
957 Double_t AliMUONAlignment::GetParError(Int_t iPar) {
958 /// Get error of parameter iPar
959 Double_t lErr = fMillepede->GetParError(iPar);
963 void AliMUONAlignment::PrintGlobalParameters() {
964 /// Print global parameters
965 fMillepede->PrintGlobalParameters();
968 //_________________________________________________________________________
969 TGeoCombiTrans AliMUONAlignment::ReAlign(const TGeoCombiTrans & transform, const double *lMisAlignment) const
971 /// Realign given transformation by given misalignment and return the misaligned transformation
973 Double_t cartMisAlig[3] = {0,0,0};
974 Double_t angMisAlig[3] = {0,0,0};
975 // const Double_t *trans = transform.GetTranslation();
976 // TGeoRotation *rot;
977 // // check if the rotation we obtain is not NULL
978 // if (transform.GetRotation()) {
979 // rot = transform.GetRotation();
982 // rot = new TGeoRotation("rot");
983 // } // default constructor.
985 cartMisAlig[0] = -TMath::Sign(1.0,transform.GetRotationMatrix()[0])*lMisAlignment[0];
986 cartMisAlig[1] = -TMath::Sign(1.0,transform.GetRotationMatrix()[4])*lMisAlignment[1];
987 cartMisAlig[2] = -TMath::Sign(1.0,transform.GetRotationMatrix()[8])*lMisAlignment[3];
988 angMisAlig[2] = -TMath::Sign(1.0,transform.GetRotationMatrix()[0]*transform.GetRotationMatrix()[4])*lMisAlignment[2]*180./TMath::Pi();
990 TGeoTranslation deltaTrans(cartMisAlig[0], cartMisAlig[1], cartMisAlig[2]);
991 TGeoRotation deltaRot;
992 deltaRot.RotateX(angMisAlig[0]);
993 deltaRot.RotateY(angMisAlig[1]);
994 deltaRot.RotateZ(angMisAlig[2]);
996 TGeoCombiTrans deltaTransf(deltaTrans,deltaRot);
997 TGeoHMatrix newTransfMat = transform * deltaTransf;
999 return TGeoCombiTrans(newTransfMat);
1002 //______________________________________________________________________
1003 AliMUONGeometryTransformer *
1004 AliMUONAlignment::ReAlign(const AliMUONGeometryTransformer * transformer,
1005 const double *misAlignments, Bool_t verbose)
1008 /// Returns a new AliMUONGeometryTransformer with the found misalignments
1011 // Takes the internal geometry module transformers, copies them
1012 // and gets the Detection Elements from them.
1013 // Takes misalignment parameters and applies these
1014 // to the local transform of the Detection Element
1015 // Obtains the global transform by multiplying the module transformer
1016 // transformation with the local transformation
1017 // Applies the global transform to a new detection element
1018 // Adds the new detection element to a new module transformer
1019 // Adds the new module transformer to a new geometry transformer
1020 // Returns the new geometry transformer
1022 Double_t lModuleMisAlignment[4] = {0.,0.,0.,0.};
1023 Double_t lDetElemMisAlignment[4] = {0.,0.,0.,0.};
1024 Int_t iDetElemId = 0;
1025 Int_t iDetElemNumber = 0;
1027 AliMUONGeometryTransformer *newGeometryTransformer =
1028 new AliMUONGeometryTransformer();
1029 for (Int_t iMt = 0; iMt < transformer->GetNofModuleTransformers(); iMt++) {
1030 // module transformers
1031 const AliMUONGeometryModuleTransformer *kModuleTransformer =
1032 transformer->GetModuleTransformer(iMt, true);
1034 AliMUONGeometryModuleTransformer *newModuleTransformer =
1035 new AliMUONGeometryModuleTransformer(iMt);
1036 newGeometryTransformer->AddModuleTransformer(newModuleTransformer);
1038 TGeoCombiTrans moduleTransform =
1039 TGeoCombiTrans(*kModuleTransformer->GetTransformation());
1040 // New module transformation
1041 TGeoCombiTrans newModuleTransform = ReAlign(moduleTransform,lModuleMisAlignment);
1042 newModuleTransformer->SetTransformation(newModuleTransform);
1044 // Get delta transformation:
1045 // Tdelta = Tnew * Told.inverse
1046 TGeoHMatrix deltaModuleTransform =
1047 AliMUONGeometryBuilder::Multiply(newModuleTransform,
1048 kModuleTransformer->GetTransformation()->Inverse());
1049 // Create module mis alignment matrix
1050 newGeometryTransformer
1051 ->AddMisAlignModule(kModuleTransformer->GetModuleId(), deltaModuleTransform);
1053 AliMpExMap *detElements = kModuleTransformer->GetDetElementStore();
1056 AliInfo(Form("%i DEs in old GeometryStore %i",detElements->GetSize(), iMt));
1058 TIter next(detElements->CreateIterator());
1059 AliMUONGeometryDetElement* detElement;
1061 while ( ( detElement = static_cast<AliMUONGeometryDetElement*>(next()) ) )
1064 // make a new detection element
1065 AliMUONGeometryDetElement *newDetElement =
1066 new AliMUONGeometryDetElement(detElement->GetId(),
1067 detElement->GetVolumePath());
1068 TString lDetElemName(detElement->GetDEName());
1069 lDetElemName.ReplaceAll("DE","");
1070 iDetElemId = lDetElemName.Atoi();
1071 iDetElemNumber = iDetElemId%100;
1072 for (int iCh=0; iCh<iDetElemId/100-1; iCh++){
1073 iDetElemNumber += fgNDetElemCh[iCh];
1075 for (int i=0; i<fgNParCh; i++) {
1076 lDetElemMisAlignment[i] = 0.0;
1077 if (iMt<fgNTrkMod) {
1078 AliInfo(Form("iMt %i, iCh %i, iDe %i, iDeId %i, iDeNb %i, iPar %i",iMt, iDetElemId/100, iDe, iDetElemId, iDetElemNumber, iDetElemNumber*fgNParCh+i));
1079 lDetElemMisAlignment[i] = misAlignments[iDetElemNumber*fgNParCh+i];
1082 // local transformation of this detection element.
1083 TGeoCombiTrans localTransform
1084 = TGeoCombiTrans(*detElement->GetLocalTransformation());
1085 TGeoCombiTrans newLocalTransform = ReAlign(localTransform,lDetElemMisAlignment);
1086 newDetElement->SetLocalTransformation(newLocalTransform);
1088 // global transformation
1089 TGeoHMatrix newGlobalTransform =
1090 AliMUONGeometryBuilder::Multiply(newModuleTransform,
1092 newDetElement->SetGlobalTransformation(newGlobalTransform);
1094 // add this det element to module
1095 newModuleTransformer->GetDetElementStore()->Add(newDetElement->GetId(),
1098 // In the Alice Alignment Framework misalignment objects store
1099 // global delta transformation
1100 // Get detection "intermediate" global transformation
1101 TGeoHMatrix newOldGlobalTransform = newModuleTransform * localTransform;
1102 // Get detection element global delta transformation:
1103 // Tdelta = Tnew * Told.inverse
1104 TGeoHMatrix deltaGlobalTransform
1105 = AliMUONGeometryBuilder::Multiply(newGlobalTransform,
1106 newOldGlobalTransform.Inverse());
1108 // Create mis alignment matrix
1109 newGeometryTransformer
1110 ->AddMisAlignDetElement(detElement->GetId(), deltaGlobalTransform);
1114 AliInfo(Form("Added module transformer %i to the transformer", iMt));
1115 newGeometryTransformer->AddModuleTransformer(newModuleTransformer);
1117 return newGeometryTransformer;
1120 //______________________________________________________________________
1121 void AliMUONAlignment::SetAlignmentResolution(const TClonesArray* misAlignArray, Int_t rChId, Double_t rChResX, Double_t rChResY, Double_t rDeResX, Double_t rDeResY){
1122 /// Set alignment resolution to misalign objects to be stored in CDB
1123 Int_t chIdMin = (rChId<0)? 0 : rChId;
1124 Int_t chIdMax = (rChId<0)? 9 : rChId;
1125 Double_t chResX = rChResX;
1126 Double_t chResY = rChResY;
1127 Double_t deResX = rDeResX;
1128 Double_t deResY = rDeResY;
1130 TMatrixDSym mChCorrMatrix(6);
1131 mChCorrMatrix[0][0]=chResX*chResX;
1132 mChCorrMatrix[1][1]=chResY*chResY;
1133 // mChCorrMatrix.Print();
1135 TMatrixDSym mDECorrMatrix(6);
1136 mDECorrMatrix[0][0]=deResX*deResX;
1137 mDECorrMatrix[1][1]=deResY*deResY;
1138 // mDECorrMatrix.Print();
1140 AliAlignObjMatrix *alignMat = 0x0;
1142 for(Int_t chId=chIdMin; chId<=chIdMax; chId++) {
1146 chName1 = Form("GM%d",chId);
1147 chName2 = Form("GM%d",chId);
1149 chName1 = Form("GM%d",4+(chId-4)*2);
1150 chName2 = Form("GM%d",4+(chId-4)*2+1);
1153 for (int i=0; i<misAlignArray->GetEntries(); i++) {
1154 alignMat = (AliAlignObjMatrix*)misAlignArray->At(i);
1155 TString volName(alignMat->GetSymName());
1156 if((volName.Contains(chName1)&&
1157 ((volName.Last('/')==volName.Index(chName1)+chName1.Length())||
1158 (volName.Length()==volName.Index(chName1)+chName1.Length())))||
1159 (volName.Contains(chName2)&&
1160 ((volName.Last('/')==volName.Index(chName2)+chName2.Length())||
1161 (volName.Length()==volName.Index(chName2)+chName2.Length())))){
1162 volName.Remove(0,volName.Last('/')+1);
1163 if (volName.Contains("GM")) {
1164 // alignMat->Print("NULL");
1165 alignMat->SetCorrMatrix(mChCorrMatrix);
1166 } else if (volName.Contains("DE")) {
1167 // alignMat->Print("NULL");
1168 alignMat->SetCorrMatrix(mDECorrMatrix);