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"
48 #include "TClonesArray.h"
51 ClassImp(AliMUONAlignment)
54 //_____________________________________________________________________
56 Int_t AliMUONAlignment::fgNDetElem = 4*2+4*2+18*2+26*2+26*2;
57 Int_t AliMUONAlignment::fgNDetElemCh[10] = {4,4,4,4,18,18,26,26,26,26};
58 Int_t AliMUONAlignment::fgSNDetElemCh[10] = {4,8,12,16,34,52,78,104,130,156};
59 Int_t AliMUONAlignment::fgNParCh = 4;
60 Int_t AliMUONAlignment::fgNTrkMod = 16;
61 Int_t AliMUONAlignment::fgNCh = 10;
62 Int_t AliMUONAlignment::fgNSt = 5;
64 //_____________________________________________________________________
65 AliMUONAlignment::AliMUONAlignment()
75 fNGlobal(fgNDetElem*fgNParCh),
90 AliInfo(Form("fSigma[0]: %f\t fSigma[1]: %f",fSigma[0],fSigma[1]));
92 fDoF[0] = kTRUE; fDoF[1] = kTRUE; fDoF[2] = kTRUE; fDoF[3] = kTRUE;
93 fAllowVar[0] = 0.05; fAllowVar[1] = 0.05; fAllowVar[2] = 0.001; fAllowVar[3] = 0.5;
95 AliInfo(Form("fAllowVar[0]: %f\t fAllowVar[1]: %f\t fPhi: %f\t fgNDetElem: %i\t fNGlobal: %i\t fNLocal: %i",fAllowVar[0],fAllowVar[1],fPhi,fgNDetElem,fNGlobal,fNLocal));
97 fMillepede = new AliMillepede();
99 Init(fNGlobal, fNLocal, fNStdDev);
101 ResetLocalEquation();
102 AliInfo("Parameters initialized to zero");
106 //_____________________________________________________________________
107 AliMUONAlignment::~AliMUONAlignment()
110 //_____________________________________________________________________
111 void AliMUONAlignment::Init(
112 Int_t nGlobal, /* number of global paramers */
113 Int_t nLocal, /* number of local parameters */
114 Int_t nStdDev /* std dev cut */ )
116 /// Initialization of AliMillepede. Fix parameters, define constraints ...
117 fMillepede->InitMille(nGlobal,nLocal,nStdDev,fResCut,fResCutInitial);
119 // Bool_t bStOnOff[5] = {kTRUE,kTRUE,kTRUE,kTRUE,kTRUE};
120 // Bool_t bChOnOff[10] = {kTRUE,kTRUE,kTRUE,kTRUE,kTRUE,kTRUE,kTRUE,kTRUE,kTRUE,kTRUE};
121 // Bool_t bSpecLROnOff[2] = {kTRUE,kTRUE};
123 // AllowVariations(bChOnOff);
125 // Fix parameters or add constraints here
126 // for (Int_t iSt=0; iSt<5; iSt++)
127 // { if (!bStOnOff[iSt]) FixStation(iSt+1); }
129 // for (Int_t iCh=0; iCh<10; iCh++)
130 // { if (!bChOnOff[iCh]) FixChamber(iCh+1); }
132 // FixHalfSpectrometer(bChOnOff,bSpecLROnOff);
136 // Define global constrains to be applied
137 // X, Y, P, XvsZ, YvsZ, PvsZ, XvsY, YvsY, PvsY
138 Bool_t bVarXYT[9] = {kTRUE,kTRUE,kTRUE,kTRUE,kTRUE,kTRUE,kTRUE,kTRUE,kTRUE};
139 Bool_t bDetTLBR[4] = {kFALSE,kTRUE,kFALSE,kTRUE};
140 // AddConstraints(bStOnOff,bVarXYT,bDetTLBR,bSpecLROnOff);
142 // Other possible way to add constrains
143 bVarXYT[0] = kFALSE; bVarXYT[1] = kFALSE; bVarXYT[2] = kTRUE;
144 bDetTLBR[0] = kFALSE; bDetTLBR[1] = kTRUE; bDetTLBR[2] = kFALSE; bDetTLBR[3] = kFALSE;
145 // AddConstraints(bStOnOff,bVarXYT,bDetTLBR);
147 bVarXYT[0] = kTRUE; bVarXYT[1] = kTRUE; bVarXYT[2] = kFALSE;
148 // AddConstraints(bStOnOff,bVarXYT);
151 if (fStartFac>1) fMillepede->SetIterations(fStartFac);
154 //_____________________________________________________________________
155 void AliMUONAlignment::FixStation(Int_t iSt)
157 /// Fix all detection elements of station iSt
158 Int_t iDetElemFirst = (iSt>1) ? fgSNDetElemCh[2*(iSt-1)-1] : 0;
159 Int_t iDetElemLast = fgSNDetElemCh[2*(iSt)-1];
160 for (Int_t i = iDetElemFirst; i < iDetElemLast; i++)
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);
170 //_____________________________________________________________________
171 void AliMUONAlignment::FixChamber(Int_t iCh)
173 /// Fix all detection elements of chamber iCh
174 Int_t iDetElemFirst = (iCh>1) ? fgSNDetElemCh[iCh-2] : 0;
175 Int_t iDetElemLast = fgSNDetElemCh[iCh-1];
176 for (Int_t i = iDetElemFirst; i < iDetElemLast; i++)
178 FixParameter(i*fgNParCh+0, 0.0);
179 FixParameter(i*fgNParCh+1, 0.0);
180 FixParameter(i*fgNParCh+2, 0.0);
181 FixParameter(i*fgNParCh+3, 0.0);
185 //_____________________________________________________________________
186 void AliMUONAlignment::FixDetElem(Int_t iDetElemId, TString sVarXYT)
189 /// Fix a given detection element
190 Int_t iDetElemNumber = iDetElemId%100;
191 for (int iCh=0; iCh<iDetElemId/100-1; iCh++)
193 iDetElemNumber += fgNDetElemCh[iCh];
196 if (sVarXYT.Contains("X"))
199 FixParameter(iDetElemNumber*fgNParCh+0, 0.0);
202 if (sVarXYT.Contains("Y"))
205 FixParameter(iDetElemNumber*fgNParCh+1, 0.0);
208 if (sVarXYT.Contains("T"))
211 FixParameter(iDetElemNumber*fgNParCh+2, 0.0);
214 if (sVarXYT.Contains("Z"))
217 FixParameter(iDetElemNumber*fgNParCh+3, 0.0);
222 //_____________________________________________________________________
223 void AliMUONAlignment::FixHalfSpectrometer(const Bool_t *lChOnOff, const Bool_t *lSpecLROnOff)
226 /// Fix left or right detector
227 for (Int_t i = 0; i < fgNDetElem; i++)
231 for (iCh=1; iCh<=fgNCh; iCh++)
232 { if (i<fgSNDetElemCh[iCh-1]) break; }
237 Int_t lDetElemNumber = (iCh==1) ? i : i-fgSNDetElemCh[iCh-2];
238 if (iCh>=1 && iCh<=4)
241 if ((lDetElemNumber==1 || lDetElemNumber==2) && !lSpecLROnOff[0])
243 // From track crossings
244 FixParameter(i*fgNParCh+0, 0.0);
245 FixParameter(i*fgNParCh+1, 0.0);
246 FixParameter(i*fgNParCh+2, 0.0);
247 FixParameter(i*fgNParCh+3, 0.0);
250 if ((lDetElemNumber==0 || lDetElemNumber==3) && !lSpecLROnOff[1])
252 // From track crossings
253 FixParameter(i*fgNParCh+0, 0.0);
254 FixParameter(i*fgNParCh+1, 0.0);
255 FixParameter(i*fgNParCh+2, 0.0);
256 FixParameter(i*fgNParCh+3, 0.0);
261 if (iCh>=5 && iCh<=6)
264 if ((lDetElemNumber>=5&&lDetElemNumber<=13) && !lSpecLROnOff[0])
266 FixParameter(i*fgNParCh+0, 0.0);
267 FixParameter(i*fgNParCh+1, 0.0);
268 FixParameter(i*fgNParCh+2, 0.0);
269 FixParameter(i*fgNParCh+3, 0.0);
272 if (((lDetElemNumber>=0&&lDetElemNumber<=4) ||
273 (lDetElemNumber>=14&&lDetElemNumber<=17)) && !lSpecLROnOff[1])
276 FixParameter(i*fgNParCh+0, 0.0);
277 FixParameter(i*fgNParCh+1, 0.0);
278 FixParameter(i*fgNParCh+2, 0.0);
279 FixParameter(i*fgNParCh+3, 0.0);
284 if (iCh>=7 && iCh<=10)
287 if ((lDetElemNumber>=7&&lDetElemNumber<=19) && !lSpecLROnOff[0])
289 FixParameter(i*fgNParCh+0, 0.0);
290 FixParameter(i*fgNParCh+1, 0.0);
291 FixParameter(i*fgNParCh+2, 0.0);
292 FixParameter(i*fgNParCh+3, 0.0);
295 if (((lDetElemNumber>=0&&lDetElemNumber<=6) ||
296 (lDetElemNumber>=20&&lDetElemNumber<=25)) && !lSpecLROnOff[1])
298 FixParameter(i*fgNParCh+0, 0.0);
299 FixParameter(i*fgNParCh+1, 0.0);
300 FixParameter(i*fgNParCh+2, 0.0);
301 FixParameter(i*fgNParCh+3, 0.0);
312 //______________________________________________________________________
313 void AliMUONAlignment::SetNonLinear(const Bool_t *lChOnOff, const Bool_t *lVarXYT)
316 /// Set non linear parameter flag selected chambers and degrees of freedom
317 for (Int_t i = 0; i < fgNDetElem; i++)
321 for (iCh=1; iCh<=fgNCh; iCh++)
322 { if (i<fgSNDetElemCh[iCh-1]) break; }
330 SetNonLinear(i*fgNParCh+0);
336 SetNonLinear(i*fgNParCh+1);
342 SetNonLinear(i*fgNParCh+2);
348 SetNonLinear(i*fgNParCh+3);
357 //______________________________________________________________________
358 void AliMUONAlignment::AddConstraints(const Bool_t *lChOnOff, const Bool_t *lVarXYT)
361 /// Add constraint equations for selected chambers and degrees of freedom
362 for (Int_t i = 0; i < fgNDetElem; i++)
366 for (iCh=1; iCh<=fgNCh; iCh++)
368 if (i<fgSNDetElemCh[iCh-1]) break;
376 fConstraintX[i*fgNParCh+0]=1.0;
382 fConstraintY[i*fgNParCh+1]=1.0;
388 fConstraintP[i*fgNParCh+2]=1.0;
390 // if (lVarXYT[3]) { // Z constraint
391 // fConstraintP[i*fgNParCh+3]=1.0;
400 AddConstraint(fConstraintX,0.0);
406 AddConstraint(fConstraintY,0.0);
412 AddConstraint(fConstraintP,0.0);
415 // if (lVarXYT[3]) { // Z constraint
416 // AddConstraint(fConstraintP,0.0);
420 //______________________________________________________________________
421 void AliMUONAlignment::AddConstraints(const Bool_t *lChOnOff, const Bool_t *lVarXYT, const Bool_t *lDetTLBR, const Bool_t *lSpecLROnOff)
423 /// Add constraint equations for selected chambers, degrees of freedom and detector half
424 Double_t lDetElemLocX = 0.;
425 Double_t lDetElemLocY = 0.;
426 Double_t lDetElemLocZ = 0.;
427 Double_t lDetElemGloX = 0.;
428 Double_t lDetElemGloY = 0.;
429 Double_t lDetElemGloZ = 0.;
430 Double_t lMeanY = 0.;
431 Double_t lSigmaY = 0.;
432 Double_t lMeanZ = 0.;
433 Double_t lSigmaZ = 0.;
435 for (Int_t i = 0; i < fgNDetElem; i++)
439 for (iCh=1; iCh<=fgNCh; iCh++){
440 if (i<fgSNDetElemCh[iCh-1]) break;
442 if (lChOnOff[iCh-1]){
443 Int_t lDetElemNumber = (iCh==1) ? i : i-fgSNDetElemCh[iCh-2];
444 Int_t lDetElemId = iCh*100+lDetElemNumber;
445 fTransform->Local2Global(lDetElemId,lDetElemLocX,lDetElemLocY,lDetElemLocZ,
446 lDetElemGloX,lDetElemGloY,lDetElemGloZ);
447 if (iCh>=1 && iCh<=4){
448 if ((lDetElemNumber==1 || lDetElemNumber==2) && lSpecLROnOff[0]){ // From track crossings
449 lMeanY += lDetElemGloY;
450 lSigmaY += lDetElemGloY*lDetElemGloY;
451 lMeanZ += lDetElemGloZ;
452 lSigmaZ += lDetElemGloZ*lDetElemGloZ;
455 if ((lDetElemNumber==0 || lDetElemNumber==3) && lSpecLROnOff[1]){ // From track crossings
456 lMeanY += lDetElemGloY;
457 lSigmaY += lDetElemGloY*lDetElemGloY;
458 lMeanZ += lDetElemGloZ;
459 lSigmaZ += lDetElemGloZ*lDetElemGloZ;
463 if (iCh>=5 && iCh<=6){
464 if ((lDetElemNumber>=5&&lDetElemNumber<=13) && lSpecLROnOff[0]){
465 lMeanY += lDetElemGloY;
466 lSigmaY += lDetElemGloY*lDetElemGloY;
467 lMeanZ += lDetElemGloZ;
468 lSigmaZ += lDetElemGloZ*lDetElemGloZ;
471 if (((lDetElemNumber>=0&&lDetElemNumber<=4) ||
472 (lDetElemNumber>=14&&lDetElemNumber<=17)) && lSpecLROnOff[1]){
473 lMeanY += lDetElemGloY;
474 lSigmaY += lDetElemGloY*lDetElemGloY;
475 lMeanZ += lDetElemGloZ;
476 lSigmaZ += lDetElemGloZ*lDetElemGloZ;
480 if (iCh>=7 && iCh<=10){
481 if ((lDetElemNumber>=7&&lDetElemNumber<=19) && lSpecLROnOff[0]){
482 lMeanY += lDetElemGloY;
483 lSigmaY += lDetElemGloY*lDetElemGloY;
484 lMeanZ += lDetElemGloZ;
485 lSigmaZ += lDetElemGloZ*lDetElemGloZ;
488 if (((lDetElemNumber>=0&&lDetElemNumber<=6) ||
489 (lDetElemNumber>=20&&lDetElemNumber<=25)) && lSpecLROnOff[1]){
490 lMeanY += lDetElemGloY;
491 lSigmaY += lDetElemGloY*lDetElemGloY;
492 lMeanZ += lDetElemGloZ;
493 lSigmaZ += lDetElemGloZ*lDetElemGloZ;
501 lSigmaY /= lNDetElem;
502 lSigmaY = TMath::Sqrt(lSigmaY-lMeanY*lMeanY);
504 lSigmaZ /= lNDetElem;
505 lSigmaZ = TMath::Sqrt(lSigmaZ-lMeanZ*lMeanZ);
506 AliInfo(Form("Used %i DetElem, MeanZ= %f , SigmaZ= %f", lNDetElem,lMeanZ,lSigmaZ));
508 AliError("No detection elements to constrain!!!");
512 for (Int_t i = 0; i < fgNDetElem; i++){
514 for (iCh=1; iCh<=fgNCh; iCh++){
515 if (i<fgSNDetElemCh[iCh-1]) break;
517 if (lChOnOff[iCh-1]){
518 Int_t lDetElemNumber = (iCh==1) ? i : i-fgSNDetElemCh[iCh-2];
519 Int_t lDetElemId = iCh*100+lDetElemNumber;
520 fTransform->Local2Global(lDetElemId,lDetElemLocX,lDetElemLocY,lDetElemLocZ,
521 lDetElemGloX,lDetElemGloY,lDetElemGloZ);
522 if (lVarXYT[0]) { // X constraint
523 if (lDetTLBR[0]) ConstrainT(i,iCh,fConstraintXT,0); // Top half
524 if (lDetTLBR[1]) ConstrainL(i,iCh,fConstraintXL,0); // Left half
525 if (lDetTLBR[2]) ConstrainB(i,iCh,fConstraintXB,0); // Bottom half
526 if (lDetTLBR[3]) ConstrainR(i,iCh,fConstraintXR,0); // Right half
528 if (lVarXYT[1]) { // Y constraint
529 if (lDetTLBR[0]) ConstrainT(i,iCh,fConstraintYT,1); // Top half
530 if (lDetTLBR[1]) ConstrainL(i,iCh,fConstraintYL,1); // Left half
531 if (lDetTLBR[2]) ConstrainB(i,iCh,fConstraintYB,1); // Bottom half
532 if (lDetTLBR[3]) ConstrainR(i,iCh,fConstraintYR,1); // Right half
534 if (lVarXYT[2]) { // P constraint
535 if (lDetTLBR[0]) ConstrainT(i,iCh,fConstraintPT,2); // Top half
536 if (lDetTLBR[1]) ConstrainL(i,iCh,fConstraintPL,2); // Left half
537 if (lDetTLBR[2]) ConstrainB(i,iCh,fConstraintPB,2); // Bottom half
538 if (lDetTLBR[3]) ConstrainR(i,iCh,fConstraintPR,2); // Right half
540 if (lVarXYT[3]) { // X-Z shearing
541 if (lDetTLBR[0]) ConstrainT(i,iCh,fConstraintXZT,0,(lDetElemGloZ-lMeanZ)/lSigmaZ); // Top half
542 if (lDetTLBR[1]) ConstrainL(i,iCh,fConstraintXZL,0,(lDetElemGloZ-lMeanZ)/lSigmaZ); // Left half
543 if (lDetTLBR[2]) ConstrainB(i,iCh,fConstraintXZB,0,(lDetElemGloZ-lMeanZ)/lSigmaZ); // Bottom half
544 if (lDetTLBR[3]) ConstrainR(i,iCh,fConstraintXZR,0,(lDetElemGloZ-lMeanZ)/lSigmaZ); // Right half
546 if (lVarXYT[4]) { // Y-Z shearing
547 if (lDetTLBR[0]) ConstrainT(i,iCh,fConstraintYZT,1,(lDetElemGloZ-lMeanZ)/lSigmaZ); // Top half
548 if (lDetTLBR[1]) ConstrainL(i,iCh,fConstraintYZL,1,(lDetElemGloZ-lMeanZ)/lSigmaZ); // Left half
549 if (lDetTLBR[2]) ConstrainB(i,iCh,fConstraintYZB,1,(lDetElemGloZ-lMeanZ)/lSigmaZ); // Bottom half
550 if (lDetTLBR[3]) ConstrainR(i,iCh,fConstraintYZR,1,(lDetElemGloZ-lMeanZ)/lSigmaZ); // Right half
552 if (lVarXYT[5]) { // P-Z rotation
553 if (lDetTLBR[0]) ConstrainT(i,iCh,fConstraintPZT,2,(lDetElemGloZ-lMeanZ)/lSigmaZ); // Top half
554 if (lDetTLBR[1]) ConstrainL(i,iCh,fConstraintPZL,2,(lDetElemGloZ-lMeanZ)/lSigmaZ); // Left half
555 if (lDetTLBR[2]) ConstrainB(i,iCh,fConstraintPZB,2,(lDetElemGloZ-lMeanZ)/lSigmaZ); // Bottom half
556 if (lDetTLBR[3]) ConstrainR(i,iCh,fConstraintPZR,2,(lDetElemGloZ-lMeanZ)/lSigmaZ); // Right half
558 if (lVarXYT[6]) { // X-Y shearing
559 if (lDetTLBR[0]) ConstrainT(i,iCh,fConstraintXYT,0,(lDetElemGloY-lMeanY)/lSigmaY); // Top half
560 if (lDetTLBR[1]) ConstrainL(i,iCh,fConstraintXYL,0,(lDetElemGloY-lMeanY)/lSigmaY); // Left half
561 if (lDetTLBR[2]) ConstrainB(i,iCh,fConstraintXYB,0,(lDetElemGloY-lMeanY)/lSigmaY); // Bottom half
562 if (lDetTLBR[3]) ConstrainR(i,iCh,fConstraintXYR,0,(lDetElemGloY-lMeanY)/lSigmaY); // Right half
564 if (lVarXYT[7]) { // Y-Y scaling
565 if (lDetTLBR[0]) ConstrainT(i,iCh,fConstraintYYT,1,(lDetElemGloY-lMeanY)/lSigmaY); // Top half
566 if (lDetTLBR[1]) ConstrainL(i,iCh,fConstraintYYL,1,(lDetElemGloY-lMeanY)/lSigmaY); // Left half
567 if (lDetTLBR[2]) ConstrainB(i,iCh,fConstraintYYB,1,(lDetElemGloY-lMeanY)/lSigmaY); // Bottom half
568 if (lDetTLBR[3]) ConstrainR(i,iCh,fConstraintYYR,1,(lDetElemGloY-lMeanY)/lSigmaY); // Right half
570 if (lVarXYT[8]) { // P-Y rotation
571 if (lDetTLBR[0]) ConstrainT(i,iCh,fConstraintPYT,2,(lDetElemGloY-lMeanY)/lSigmaY); // Top half
572 if (lDetTLBR[1]) ConstrainL(i,iCh,fConstraintPYL,2,(lDetElemGloY-lMeanY)/lSigmaY); // Left half
573 if (lDetTLBR[2]) ConstrainB(i,iCh,fConstraintPYB,2,(lDetElemGloY-lMeanY)/lSigmaY); // Bottom half
574 if (lDetTLBR[3]) ConstrainR(i,iCh,fConstraintPYR,2,(lDetElemGloY-lMeanY)/lSigmaY); // Right half
578 if (lVarXYT[0]) { // X constraint
579 if (lDetTLBR[0]) AddConstraint(fConstraintXT,0.0); // Top half
580 if (lDetTLBR[1]) AddConstraint(fConstraintXL,0.0); // Left half
581 if (lDetTLBR[2]) AddConstraint(fConstraintXB,0.0); // Bottom half
582 if (lDetTLBR[3]) AddConstraint(fConstraintXR,0.0); // Right half
584 if (lVarXYT[1]) { // Y constraint
585 if (lDetTLBR[0]) AddConstraint(fConstraintYT,0.0); // Top half
586 if (lDetTLBR[1]) AddConstraint(fConstraintYL,0.0); // Left half
587 if (lDetTLBR[2]) AddConstraint(fConstraintYB,0.0); // Bottom half
588 if (lDetTLBR[3]) AddConstraint(fConstraintYR,0.0); // Right half
590 if (lVarXYT[2]) { // T constraint
591 if (lDetTLBR[0]) AddConstraint(fConstraintPT,0.0); // Top half
592 if (lDetTLBR[1]) AddConstraint(fConstraintPL,0.0); // Left half
593 if (lDetTLBR[2]) AddConstraint(fConstraintPB,0.0); // Bottom half
594 if (lDetTLBR[3]) AddConstraint(fConstraintPR,0.0); // Right half
596 if (lVarXYT[3]) { // X-Z constraint
597 if (lDetTLBR[0]) AddConstraint(fConstraintXZT,0.0); // Top half
598 if (lDetTLBR[1]) AddConstraint(fConstraintXZL,0.0); // Left half
599 if (lDetTLBR[2]) AddConstraint(fConstraintXZB,0.0); // Bottom half
600 if (lDetTLBR[3]) AddConstraint(fConstraintXZR,0.0); // Right half
602 if (lVarXYT[4]) { // Y-Z constraint
603 if (lDetTLBR[0]) AddConstraint(fConstraintYZT,0.0); // Top half
604 if (lDetTLBR[1]) AddConstraint(fConstraintYZL,0.0); // Left half
605 if (lDetTLBR[2]) AddConstraint(fConstraintYZB,0.0); // Bottom half
606 if (lDetTLBR[3]) AddConstraint(fConstraintYZR,0.0); // Right half
608 if (lVarXYT[5]) { // P-Z constraint
609 if (lDetTLBR[0]) AddConstraint(fConstraintPZT,0.0); // Top half
610 if (lDetTLBR[1]) AddConstraint(fConstraintPZL,0.0); // Left half
611 if (lDetTLBR[2]) AddConstraint(fConstraintPZB,0.0); // Bottom half
612 if (lDetTLBR[3]) AddConstraint(fConstraintPZR,0.0); // Right half
614 if (lVarXYT[6]) { // X-Y constraint
615 if (lDetTLBR[0]) AddConstraint(fConstraintXYT,0.0); // Top half
616 if (lDetTLBR[1]) AddConstraint(fConstraintXYL,0.0); // Left half
617 if (lDetTLBR[2]) AddConstraint(fConstraintXYB,0.0); // Bottom half
618 if (lDetTLBR[3]) AddConstraint(fConstraintXYR,0.0); // Right half
620 if (lVarXYT[7]) { // Y-Y constraint
621 if (lDetTLBR[0]) AddConstraint(fConstraintYYT,0.0); // Top half
622 if (lDetTLBR[1]) AddConstraint(fConstraintYYL,0.0); // Left half
623 if (lDetTLBR[2]) AddConstraint(fConstraintYYB,0.0); // Bottom half
624 if (lDetTLBR[3]) AddConstraint(fConstraintYYR,0.0); // Right half
626 if (lVarXYT[8]) { // P-Y constraint
627 if (lDetTLBR[0]) AddConstraint(fConstraintPYT,0.0); // Top half
628 if (lDetTLBR[1]) AddConstraint(fConstraintPYL,0.0); // Left half
629 if (lDetTLBR[2]) AddConstraint(fConstraintPYB,0.0); // Bottom half
630 if (lDetTLBR[3]) AddConstraint(fConstraintPYR,0.0); // Right half
634 void AliMUONAlignment::ConstrainT(Int_t lDetElem, Int_t lCh, Double_t *lConstraintT, Int_t iVar, Double_t /*lWeight*/) const{
635 /// Set constrain equation for top half of spectrometer
636 Int_t lDetElemNumber = (lCh==1) ? lDetElem : lDetElem-fgSNDetElemCh[lCh-2];
637 if (lCh>=1 && lCh<=4){
638 if (lDetElemNumber==0 || lDetElemNumber==1){ // From track crossings
639 lConstraintT[lDetElem*fgNParCh+iVar]=1.0;
642 if (lCh>=5 && lCh<=6){
643 if (lDetElemNumber>=0&&lDetElemNumber<=9){
644 lConstraintT[lDetElem*fgNParCh+iVar]=1.0;
647 if (lCh>=7 && lCh<=10){
648 if (lDetElemNumber>=0&&lDetElemNumber<=13){
649 lConstraintT[lDetElem*fgNParCh+iVar]=1.0;
654 //______________________________________________________________________
655 void AliMUONAlignment::ConstrainL(Int_t lDetElem, Int_t lCh, Double_t *lConstraintL, Int_t iVar, Double_t lWeight) const{
656 /// Set constrain equation for left half of spectrometer
657 Int_t lDetElemNumber = (lCh==1) ? lDetElem : lDetElem-fgSNDetElemCh[lCh-2];
658 if (lCh>=1 && lCh<=4){
659 if (lDetElemNumber==1 || lDetElemNumber==2){ // From track crossings
660 lConstraintL[lDetElem*fgNParCh+iVar]=lWeight;
663 if (lCh>=5 && lCh<=6){
664 if (lDetElemNumber>=5&&lDetElemNumber<=13){
665 lConstraintL[lDetElem*fgNParCh+iVar]=lWeight;
668 if (lCh>=7 && lCh<=10){
669 if (lDetElemNumber>=7&&lDetElemNumber<=19){
670 lConstraintL[lDetElem*fgNParCh+iVar]=lWeight;
675 //______________________________________________________________________
676 void AliMUONAlignment::ConstrainB(Int_t lDetElem, Int_t lCh, Double_t *lConstraintB, Int_t iVar, Double_t /*lWeight*/) const{
677 /// Set constrain equation for bottom half of spectrometer
678 Int_t lDetElemNumber = (lCh==1) ? lDetElem : lDetElem-fgSNDetElemCh[lCh-2];
679 if (lCh>=1 && lCh<=4){
680 if (lDetElemNumber==2 || lDetElemNumber==3){ // From track crossings
681 lConstraintB[lDetElem*fgNParCh+iVar]=1.0;
684 if (lCh>=5 && lCh<=6){
685 if ((lDetElemNumber>=9&&lDetElemNumber<=17) ||
686 (lDetElemNumber==0)){
687 lConstraintB[lDetElem*fgNParCh+iVar]=1.0;
690 if (lCh>=7 && lCh<=10){
691 if ((lDetElemNumber>=13&&lDetElemNumber<=25) ||
692 (lDetElemNumber==0)){
693 lConstraintB[lDetElem*fgNParCh+iVar]=1.0;
698 //______________________________________________________________________
699 void AliMUONAlignment::ConstrainR(Int_t lDetElem, Int_t lCh, Double_t *lConstraintR, Int_t iVar, Double_t lWeight) const{
700 /// Set constrain equation for right half of spectrometer
701 Int_t lDetElemNumber = (lCh==1) ? lDetElem : lDetElem-fgSNDetElemCh[lCh-2];
702 if (lCh>=1 && lCh<=4){
703 if (lDetElemNumber==0 || lDetElemNumber==3){ // From track crossings
704 lConstraintR[lDetElem*fgNParCh+iVar]=lWeight;
707 if (lCh>=5 && lCh<=6){
708 if ((lDetElemNumber>=0&&lDetElemNumber<=4) ||
709 (lDetElemNumber>=14&&lDetElemNumber<=17)){
710 lConstraintR[lDetElem*fgNParCh+iVar]=lWeight;
713 if (lCh>=7 && lCh<=10){
714 if ((lDetElemNumber>=0&&lDetElemNumber<=6) ||
715 (lDetElemNumber>=20&&lDetElemNumber<=25)){
716 lConstraintR[lDetElem*fgNParCh+iVar]=lWeight;
721 //______________________________________________________________________
722 void AliMUONAlignment::ResetConstraints(){
723 /// Reset all constraint equations
724 for (Int_t i = 0; i < fgNDetElem; i++){
725 fConstraintX[i*fgNParCh+0]=0.0;
726 fConstraintX[i*fgNParCh+1]=0.0;
727 fConstraintX[i*fgNParCh+2]=0.0;
728 fConstraintY[i*fgNParCh+0]=0.0;
729 fConstraintY[i*fgNParCh+1]=0.0;
730 fConstraintY[i*fgNParCh+2]=0.0;
731 fConstraintP[i*fgNParCh+0]=0.0;
732 fConstraintP[i*fgNParCh+1]=0.0;
733 fConstraintP[i*fgNParCh+2]=0.0;
734 fConstraintXT[i*fgNParCh+0]=0.0;
735 fConstraintXT[i*fgNParCh+1]=0.0;
736 fConstraintXT[i*fgNParCh+2]=0.0;
737 fConstraintYT[i*fgNParCh+0]=0.0;
738 fConstraintYT[i*fgNParCh+1]=0.0;
739 fConstraintYT[i*fgNParCh+2]=0.0;
740 fConstraintPT[i*fgNParCh+0]=0.0;
741 fConstraintPT[i*fgNParCh+1]=0.0;
742 fConstraintPT[i*fgNParCh+2]=0.0;
743 fConstraintXZT[i*fgNParCh+0]=0.0;
744 fConstraintXZT[i*fgNParCh+1]=0.0;
745 fConstraintXZT[i*fgNParCh+2]=0.0;
746 fConstraintYZT[i*fgNParCh+0]=0.0;
747 fConstraintYZT[i*fgNParCh+1]=0.0;
748 fConstraintYZT[i*fgNParCh+2]=0.0;
749 fConstraintPZT[i*fgNParCh+0]=0.0;
750 fConstraintPZT[i*fgNParCh+1]=0.0;
751 fConstraintPZT[i*fgNParCh+2]=0.0;
752 fConstraintXYT[i*fgNParCh+0]=0.0;
753 fConstraintXYT[i*fgNParCh+1]=0.0;
754 fConstraintXYT[i*fgNParCh+2]=0.0;
755 fConstraintYYT[i*fgNParCh+0]=0.0;
756 fConstraintYYT[i*fgNParCh+1]=0.0;
757 fConstraintYYT[i*fgNParCh+2]=0.0;
758 fConstraintPYT[i*fgNParCh+0]=0.0;
759 fConstraintPYT[i*fgNParCh+1]=0.0;
760 fConstraintPYT[i*fgNParCh+2]=0.0;
761 fConstraintXL[i*fgNParCh+0]=0.0;
762 fConstraintXL[i*fgNParCh+1]=0.0;
763 fConstraintXL[i*fgNParCh+2]=0.0;
764 fConstraintYL[i*fgNParCh+0]=0.0;
765 fConstraintYL[i*fgNParCh+1]=0.0;
766 fConstraintYL[i*fgNParCh+2]=0.0;
767 fConstraintPL[i*fgNParCh+0]=0.0;
768 fConstraintPL[i*fgNParCh+1]=0.0;
769 fConstraintPL[i*fgNParCh+2]=0.0;
770 fConstraintXZL[i*fgNParCh+0]=0.0;
771 fConstraintXZL[i*fgNParCh+1]=0.0;
772 fConstraintXZL[i*fgNParCh+2]=0.0;
773 fConstraintYZL[i*fgNParCh+0]=0.0;
774 fConstraintYZL[i*fgNParCh+1]=0.0;
775 fConstraintYZL[i*fgNParCh+2]=0.0;
776 fConstraintPZL[i*fgNParCh+0]=0.0;
777 fConstraintPZL[i*fgNParCh+1]=0.0;
778 fConstraintPZL[i*fgNParCh+2]=0.0;
779 fConstraintXYL[i*fgNParCh+0]=0.0;
780 fConstraintXYL[i*fgNParCh+1]=0.0;
781 fConstraintXYL[i*fgNParCh+2]=0.0;
782 fConstraintYYL[i*fgNParCh+0]=0.0;
783 fConstraintYYL[i*fgNParCh+1]=0.0;
784 fConstraintYYL[i*fgNParCh+2]=0.0;
785 fConstraintPYL[i*fgNParCh+0]=0.0;
786 fConstraintPYL[i*fgNParCh+1]=0.0;
787 fConstraintPYL[i*fgNParCh+2]=0.0;
788 fConstraintXB[i*fgNParCh+0]=0.0;
789 fConstraintXB[i*fgNParCh+1]=0.0;
790 fConstraintXB[i*fgNParCh+2]=0.0;
791 fConstraintYB[i*fgNParCh+0]=0.0;
792 fConstraintYB[i*fgNParCh+1]=0.0;
793 fConstraintYB[i*fgNParCh+2]=0.0;
794 fConstraintPB[i*fgNParCh+0]=0.0;
795 fConstraintPB[i*fgNParCh+1]=0.0;
796 fConstraintPB[i*fgNParCh+2]=0.0;
797 fConstraintXZB[i*fgNParCh+0]=0.0;
798 fConstraintXZB[i*fgNParCh+1]=0.0;
799 fConstraintXZB[i*fgNParCh+2]=0.0;
800 fConstraintYZB[i*fgNParCh+0]=0.0;
801 fConstraintYZB[i*fgNParCh+1]=0.0;
802 fConstraintYZB[i*fgNParCh+2]=0.0;
803 fConstraintPZB[i*fgNParCh+0]=0.0;
804 fConstraintPZB[i*fgNParCh+1]=0.0;
805 fConstraintPZB[i*fgNParCh+2]=0.0;
806 fConstraintXYB[i*fgNParCh+0]=0.0;
807 fConstraintXYB[i*fgNParCh+1]=0.0;
808 fConstraintXYB[i*fgNParCh+2]=0.0;
809 fConstraintYYB[i*fgNParCh+0]=0.0;
810 fConstraintYYB[i*fgNParCh+1]=0.0;
811 fConstraintYYB[i*fgNParCh+2]=0.0;
812 fConstraintPYB[i*fgNParCh+0]=0.0;
813 fConstraintPYB[i*fgNParCh+1]=0.0;
814 fConstraintPYB[i*fgNParCh+2]=0.0;
815 fConstraintXR[i*fgNParCh+0]=0.0;
816 fConstraintXR[i*fgNParCh+1]=0.0;
817 fConstraintXR[i*fgNParCh+2]=0.0;
818 fConstraintYR[i*fgNParCh+0]=0.0;
819 fConstraintYR[i*fgNParCh+1]=0.0;
820 fConstraintYR[i*fgNParCh+2]=0.0;
821 fConstraintPR[i*fgNParCh+0]=0.0;
822 fConstraintPR[i*fgNParCh+1]=0.0;
823 fConstraintPR[i*fgNParCh+2]=0.0;
824 fConstraintXZR[i*fgNParCh+0]=0.0;
825 fConstraintXZR[i*fgNParCh+1]=0.0;
826 fConstraintXZR[i*fgNParCh+2]=0.0;
827 fConstraintYZR[i*fgNParCh+0]=0.0;
828 fConstraintYZR[i*fgNParCh+1]=0.0;
829 fConstraintYZR[i*fgNParCh+2]=0.0;
830 fConstraintPZR[i*fgNParCh+0]=0.0;
831 fConstraintPZR[i*fgNParCh+1]=0.0;
832 fConstraintPZR[i*fgNParCh+2]=0.0;
833 fConstraintPZR[i*fgNParCh+0]=0.0;
834 fConstraintPZR[i*fgNParCh+1]=0.0;
835 fConstraintPZR[i*fgNParCh+2]=0.0;
836 fConstraintXYR[i*fgNParCh+0]=0.0;
837 fConstraintXYR[i*fgNParCh+1]=0.0;
838 fConstraintXYR[i*fgNParCh+2]=0.0;
839 fConstraintYYR[i*fgNParCh+0]=0.0;
840 fConstraintYYR[i*fgNParCh+1]=0.0;
841 fConstraintYYR[i*fgNParCh+2]=0.0;
842 fConstraintPYR[i*fgNParCh+0]=0.0;
843 fConstraintPYR[i*fgNParCh+1]=0.0;
844 fConstraintPYR[i*fgNParCh+2]=0.0;
848 //______________________________________________________________________
849 void AliMUONAlignment::AddConstraint(Double_t *par, Double_t value) {
850 /// Constrain equation defined by par to value
851 fMillepede->SetGlobalConstraint(par, value);
852 AliInfo("Adding constraint");
855 //______________________________________________________________________
856 void AliMUONAlignment::InitGlobalParameters(Double_t *par) {
857 /// Initialize global parameters with par array
858 fMillepede->SetGlobalParameters(par);
859 AliInfo("Init Global Parameters");
862 //______________________________________________________________________
863 void AliMUONAlignment::FixParameter(Int_t iPar, Double_t value) {
864 /// Parameter iPar is encourage to vary in [-value;value].
865 /// If value == 0, parameter is fixed
866 fMillepede->SetParSigma(iPar, value);
867 if (TMath::Abs(value)<1e-4) AliInfo(Form("Parameter %i Fixed", iPar));
870 //______________________________________________________________________
871 void AliMUONAlignment::ResetLocalEquation()
873 /// Reset the derivative vectors
874 for(int i=0; i<fNLocal; i++) {
875 fLocalDerivatives[i] = 0.0;
877 for(int i=0; i<fNGlobal; i++) {
878 fGlobalDerivatives[i] = 0.0;
882 //______________________________________________________________________
883 void AliMUONAlignment::AllowVariations(const Bool_t *bChOnOff)
886 /// Set allowed variation for selected chambers based on fDoF and fAllowVar
887 for (Int_t iCh=1; iCh<=10; iCh++)
892 Int_t iDetElemFirst = (iCh>1) ? fgSNDetElemCh[iCh-2] : 0;
893 Int_t iDetElemLast = fgSNDetElemCh[iCh-1];
894 for (int i=0; i<fgNParCh; i++)
896 AliDebug(1,Form("fDoF[%d]= %d",i,fDoF[i]));
899 for (Int_t j=iDetElemFirst; j<iDetElemLast; j++){
900 FixParameter(j*fgNParCh+i, fAllowVar[i]);
912 //______________________________________________________________________
913 void AliMUONAlignment::SetNonLinear(Int_t iPar /* set non linear flag */ )
915 /// Set nonlinear flag for parameter iPar
916 fMillepede->SetNonLinear(iPar);
917 AliInfo(Form("Parameter %i set to non linear", iPar));
920 //______________________________________________________________________
921 void AliMUONAlignment::SetSigmaXY(Double_t sigmaX, Double_t sigmaY)
924 /// Set expected measurement resolution
925 fSigma[0] = sigmaX; fSigma[1] = sigmaY;
926 AliInfo(Form("Using fSigma[0]=%f and fSigma[1]=%f",fSigma[0],fSigma[1]));
930 //______________________________________________________________________
931 void AliMUONAlignment::LocalEquationX( Bool_t doAlignment )
935 // local cluster record
936 AliMUONAlignmentClusterRecord clusterRecord;
938 // store detector and measurement
939 clusterRecord.SetDetElemId( fDetElemId );
940 clusterRecord.SetDetElemNumber( fDetElemNumber );
941 clusterRecord.SetMeas( fMeas[0] );
942 clusterRecord.SetSigma( fSigma[0] );
944 // store local derivatives
945 clusterRecord.SetLocalDerivative( 0, fCosPhi );
946 clusterRecord.SetLocalDerivative( 1, fCosPhi*(fTrackPos[2] - fTrackPos0[2]) );
947 clusterRecord.SetLocalDerivative( 2, fSinPhi );
948 clusterRecord.SetLocalDerivative( 3, fSinPhi*(fTrackPos[2] - fTrackPos0[2]) );
950 // store global derivatives
951 clusterRecord.SetGlobalDerivative( 0, -fCosPhi );
952 clusterRecord.SetGlobalDerivative( 1, -fSinPhi );
957 clusterRecord.SetGlobalDerivative(
959 -fSinPhi*(fTrackPos[0]-fDetElemPos[0])
960 +fCosPhi*(fTrackPos[1]-fDetElemPos[1]) );
964 clusterRecord.SetGlobalDerivative(
966 -fSinPhi*(fTrackPos0[0]+fTrackSlope0[0]*(fTrackPos[2]-fTrackPos0[2])-fDetElemPos[0])
967 +fCosPhi*(fTrackPos0[1]+fTrackSlope0[1]*(fTrackPos[2]-fTrackPos0[2])-fDetElemPos[1]) );
971 clusterRecord.SetGlobalDerivative( 3, fCosPhi*fTrackSlope0[0]+fSinPhi*fTrackSlope0[1] );
973 // append to trackRecord
974 fTrackRecord.AddClusterRecord( clusterRecord );
976 // store local equation
977 if( doAlignment ) LocalEquation( clusterRecord );
981 //______________________________________________________________________
982 void AliMUONAlignment::LocalEquationY(Bool_t doAlignment )
985 // local cluster record
986 AliMUONAlignmentClusterRecord clusterRecord;
988 // store detector and measurement
989 clusterRecord.SetDetElemId( fDetElemId );
990 clusterRecord.SetDetElemNumber( fDetElemNumber );
991 clusterRecord.SetMeas( fMeas[1] );
992 clusterRecord.SetSigma( fSigma[1] );
994 // store local derivatives
995 clusterRecord.SetLocalDerivative( 0, -fSinPhi );
996 clusterRecord.SetLocalDerivative( 1, -fSinPhi*(fTrackPos[2] - fTrackPos0[2] ) );
997 clusterRecord.SetLocalDerivative( 2, fCosPhi );
998 clusterRecord.SetLocalDerivative( 3, fCosPhi*(fTrackPos[2] - fTrackPos0[2] ) );
1000 // set global derivatives
1001 clusterRecord.SetGlobalDerivative( 0, fSinPhi);
1002 clusterRecord.SetGlobalDerivative( 1, -fCosPhi);
1007 clusterRecord.SetGlobalDerivative(
1009 -fCosPhi*(fTrackPos[0]-fDetElemPos[0])
1010 -fSinPhi*(fTrackPos[1]-fDetElemPos[1]));
1014 clusterRecord.SetGlobalDerivative(
1016 -fCosPhi*(fTrackPos0[0]+fTrackSlope0[0]*(fTrackPos[2]-fTrackPos0[2])-fDetElemPos[0])
1017 -fSinPhi*(fTrackPos0[1]+fTrackSlope0[1]*(fTrackPos[2]-fTrackPos0[2])-fDetElemPos[1]));
1020 clusterRecord.SetGlobalDerivative( 3, -fSinPhi*fTrackSlope0[0]+fCosPhi*fTrackSlope0[1]);
1022 // append to trackRecord
1023 fTrackRecord.AddClusterRecord( clusterRecord );
1025 // store local equation
1026 if( doAlignment ) LocalEquation( clusterRecord );
1030 //_____________________________________________________
1031 void AliMUONAlignment::LocalEquation( const AliMUONAlignmentClusterRecord& clusterRecord )
1034 // copy to local variables
1035 for( Int_t index = 0; index < 4; ++index )
1037 SetLocalDerivative( index, clusterRecord.GetLocalDerivative( index ) );
1038 SetGlobalDerivative( clusterRecord.GetDetElemNumber()*fgNParCh + index, clusterRecord.GetGlobalDerivative( index ) );
1041 // pass equation parameters to millepede
1042 fMillepede->SetLocalEquation( fGlobalDerivatives, fLocalDerivatives, clusterRecord.GetMeas(), clusterRecord.GetSigma() );
1046 //______________________________________________________________________
1047 void AliMUONAlignment::FillRecPointData()
1050 /// Get information of current cluster
1051 fClustPos[0] = fCluster->GetX();
1052 fClustPos[1] = fCluster->GetY();
1053 fClustPos[2] = fCluster->GetZ();
1054 fTransform->Global2Local(
1055 fDetElemId,fClustPos[0],fClustPos[1],fClustPos[2],
1056 fClustPosLoc[0],fClustPosLoc[1],fClustPosLoc[2]);
1060 //______________________________________________________________________
1061 void AliMUONAlignment::FillTrackParamData()
1064 /// Get information of current track at current cluster
1065 fTrackPos[0] = fTrackParam->GetNonBendingCoor();
1066 fTrackPos[1] = fTrackParam->GetBendingCoor();
1067 fTrackPos[2] = fTrackParam->GetZ();
1068 fTrackSlope[0] = fTrackParam->GetNonBendingSlope();
1069 fTrackSlope[1] = fTrackParam->GetBendingSlope();
1070 fTransform->Global2Local(
1071 fDetElemId,fTrackPos[0],fTrackPos[1],fTrackPos[2],
1072 fTrackPosLoc[0],fTrackPosLoc[1],fTrackPosLoc[2]);
1076 //_____________________________________________________
1077 void AliMUONAlignment::FillDetElemData()
1080 /// Get information of current detection element
1081 Double_t lDetElemLocX = 0.;
1082 Double_t lDetElemLocY = 0.;
1083 Double_t lDetElemLocZ = 0.;
1084 fDetElemId = fCluster->GetDetElemId();
1085 fDetElemNumber = fDetElemId%100;
1086 for (int iCh=0; iCh<fDetElemId/100-1; iCh++)
1087 { fDetElemNumber += fgNDetElemCh[iCh]; }
1089 fTransform->Local2Global(
1090 fDetElemId,lDetElemLocX,lDetElemLocY,lDetElemLocZ,
1091 fDetElemPos[0],fDetElemPos[1],fDetElemPos[2]);
1095 //_____________________________________________________
1096 AliMUONAlignmentTrackRecord* AliMUONAlignment::ProcessTrack( AliMUONTrack* track, Bool_t doAlignment )
1099 // store current track in running member.
1102 // clear track record
1103 fTrackRecord.Clear();
1105 // get number of tracks
1106 Int_t nTrackParam = fTrack->GetTrackParamAtCluster()->GetEntries();
1107 AliDebug(1,Form("Number of track param entries : %i ", nTrackParam));
1109 Bool_t first( kTRUE );
1110 for(Int_t iCluster=0; iCluster<nTrackParam; iCluster++)
1113 // and get new pointers
1114 fTrackParam = (AliMUONTrackParam *) fTrack->GetTrackParamAtCluster()->At(iCluster);
1115 if ( ! fTrackParam ) continue;
1116 fCluster = fTrackParam->GetClusterPtr();
1117 if ( ! fCluster ) continue;
1118 // if (fDetElemId<500) continue;
1120 // fill local variables for this position --> one measurement
1123 FillTrackParamData();
1128 // for first valid cluster, save track position as "starting" values
1131 fTrackPos0[0] = fTrackPos[0];
1132 fTrackPos0[1] = fTrackPos[1];
1133 fTrackPos0[2] = fTrackPos[2];
1134 fTrackSlope0[0] = fTrackSlope[0];
1135 fTrackSlope0[1] = fTrackSlope[1];
1139 // calculate measurements
1140 fCosPhi = TMath::Cos(fPhi);
1141 fSinPhi = TMath::Sin(fPhi);
1145 fMeas[0] = fTrackPos[0] - fClustPos[0];
1146 fMeas[1] = fTrackPos[1] - fClustPos[1];
1150 fMeas[0] = - fClustPos[0];
1151 fMeas[1] = - fClustPos[1];
1156 AliDebug(1,Form("cluster: %i", iCluster));
1157 AliDebug(1,Form("x: %f\t y: %f\t z: %f\t DetElemID: %i\t ", fClustPos[0], fClustPos[1], fClustPos[2], fDetElemId));
1158 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));
1160 AliDebug(1,Form("Track Parameter: %i", iCluster));
1161 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]));
1162 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]));
1164 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]));
1166 // Set local equations
1167 LocalEquationX( doAlignment );
1168 LocalEquationY( doAlignment );
1172 return &fTrackRecord;
1176 //______________________________________________________________________________
1177 void AliMUONAlignment::ProcessTrack( AliMUONAlignmentTrackRecord* track, Bool_t doAlignment )
1180 if( !( track && doAlignment ) ) return;
1182 // loop over clusters
1183 for( Int_t index = 0; index < track->GetNRecords(); ++index )
1184 { if( AliMUONAlignmentClusterRecord* record = track->GetRecord( index ) ) LocalEquation( *record ); }
1190 //_____________________________________________________
1191 void AliMUONAlignment::LocalFit(Int_t iTrack, Double_t *lTrackParam, Int_t lSingleFit)
1194 /// Call local fit for this tracks
1195 Int_t iRes = fMillepede->LocalFit(iTrack,lTrackParam,lSingleFit);
1196 if (iRes && !lSingleFit)
1197 { fMillepede->SetNLocalEquations(fMillepede->GetNLocalEquations()+1); }
1201 //_____________________________________________________
1202 void AliMUONAlignment::GlobalFit(Double_t *parameters,Double_t *errors,Double_t *pulls)
1205 /// Call global fit; Global parameters are stored in parameters
1206 fMillepede->GlobalFit(parameters,errors,pulls);
1208 AliInfo("Done fitting global parameters!");
1209 for (int iGlob=0; iGlob<fgNDetElem; iGlob++)
1210 { printf("%d\t %f\t %f\t %f\t %f \n",iGlob,parameters[iGlob*fgNParCh+0],parameters[iGlob*fgNParCh+1],parameters[iGlob*fgNParCh+3],parameters[iGlob*fgNParCh+2]); }
1215 //_____________________________________________________
1216 Double_t AliMUONAlignment::GetParError(Int_t iPar)
1218 /// Get error of parameter iPar
1219 Double_t lErr = fMillepede->GetParError(iPar);
1223 //_____________________________________________________
1224 void AliMUONAlignment::PrintGlobalParameters()
1226 /// Print global parameters
1227 fMillepede->PrintGlobalParameters();
1230 //_________________________________________________________________________
1231 TGeoCombiTrans AliMUONAlignment::ReAlign(const TGeoCombiTrans & transform, const double *lMisAlignment) const
1233 /// Realign given transformation by given misalignment and return the misaligned transformation
1235 Double_t cartMisAlig[3] = {0,0,0};
1236 Double_t angMisAlig[3] = {0,0,0};
1237 // const Double_t *trans = transform.GetTranslation();
1238 // TGeoRotation *rot;
1239 // // check if the rotation we obtain is not NULL
1240 // if (transform.GetRotation()) {
1241 // rot = transform.GetRotation();
1244 // rot = new TGeoRotation("rot");
1245 // } // default constructor.
1247 cartMisAlig[0] = -TMath::Sign(1.0,transform.GetRotationMatrix()[0])*lMisAlignment[0];
1248 cartMisAlig[1] = -TMath::Sign(1.0,transform.GetRotationMatrix()[4])*lMisAlignment[1];
1249 cartMisAlig[2] = -TMath::Sign(1.0,transform.GetRotationMatrix()[8])*lMisAlignment[3];
1250 angMisAlig[2] = -TMath::Sign(1.0,transform.GetRotationMatrix()[0]*transform.GetRotationMatrix()[4])*lMisAlignment[2]*180./TMath::Pi();
1252 TGeoTranslation deltaTrans(cartMisAlig[0], cartMisAlig[1], cartMisAlig[2]);
1253 TGeoRotation deltaRot;
1254 deltaRot.RotateX(angMisAlig[0]);
1255 deltaRot.RotateY(angMisAlig[1]);
1256 deltaRot.RotateZ(angMisAlig[2]);
1258 TGeoCombiTrans deltaTransf(deltaTrans,deltaRot);
1259 TGeoHMatrix newTransfMat = transform * deltaTransf;
1261 return TGeoCombiTrans(newTransfMat);
1264 //______________________________________________________________________
1265 AliMUONGeometryTransformer *
1266 AliMUONAlignment::ReAlign(const AliMUONGeometryTransformer * transformer,
1267 const double *misAlignments, Bool_t verbose)
1270 /// Returns a new AliMUONGeometryTransformer with the found misalignments
1273 // Takes the internal geometry module transformers, copies them
1274 // and gets the Detection Elements from them.
1275 // Takes misalignment parameters and applies these
1276 // to the local transform of the Detection Element
1277 // Obtains the global transform by multiplying the module transformer
1278 // transformation with the local transformation
1279 // Applies the global transform to a new detection element
1280 // Adds the new detection element to a new module transformer
1281 // Adds the new module transformer to a new geometry transformer
1282 // Returns the new geometry transformer
1284 Double_t lModuleMisAlignment[4] = {0.,0.,0.,0.};
1285 Double_t lDetElemMisAlignment[4] = {0.,0.,0.,0.};
1286 Int_t iDetElemId = 0;
1287 Int_t iDetElemNumber = 0;
1289 AliMUONGeometryTransformer *newGeometryTransformer =
1290 new AliMUONGeometryTransformer();
1291 for (Int_t iMt = 0; iMt < transformer->GetNofModuleTransformers(); iMt++) {
1292 // module transformers
1293 const AliMUONGeometryModuleTransformer *kModuleTransformer =
1294 transformer->GetModuleTransformer(iMt, true);
1296 AliMUONGeometryModuleTransformer *newModuleTransformer =
1297 new AliMUONGeometryModuleTransformer(iMt);
1298 newGeometryTransformer->AddModuleTransformer(newModuleTransformer);
1300 TGeoCombiTrans moduleTransform =
1301 TGeoCombiTrans(*kModuleTransformer->GetTransformation());
1302 // New module transformation
1303 TGeoCombiTrans newModuleTransform = ReAlign(moduleTransform,lModuleMisAlignment);
1304 newModuleTransformer->SetTransformation(newModuleTransform);
1306 // Get delta transformation:
1307 // Tdelta = Tnew * Told.inverse
1308 TGeoHMatrix deltaModuleTransform =
1309 AliMUONGeometryBuilder::Multiply(newModuleTransform,
1310 kModuleTransformer->GetTransformation()->Inverse());
1311 // Create module mis alignment matrix
1312 newGeometryTransformer
1313 ->AddMisAlignModule(kModuleTransformer->GetModuleId(), deltaModuleTransform);
1315 AliMpExMap *detElements = kModuleTransformer->GetDetElementStore();
1318 AliInfo(Form("%i DEs in old GeometryStore %i",detElements->GetSize(), iMt));
1320 TIter next(detElements->CreateIterator());
1321 AliMUONGeometryDetElement* detElement;
1323 while ( ( detElement = static_cast<AliMUONGeometryDetElement*>(next()) ) )
1326 // make a new detection element
1327 AliMUONGeometryDetElement *newDetElement =
1328 new AliMUONGeometryDetElement(detElement->GetId(),
1329 detElement->GetVolumePath());
1330 TString lDetElemName(detElement->GetDEName());
1331 lDetElemName.ReplaceAll("DE","");
1332 iDetElemId = lDetElemName.Atoi();
1333 iDetElemNumber = iDetElemId%100;
1334 for (int iCh=0; iCh<iDetElemId/100-1; iCh++){
1335 iDetElemNumber += fgNDetElemCh[iCh];
1337 for (int i=0; i<fgNParCh; i++) {
1338 lDetElemMisAlignment[i] = 0.0;
1339 if (iMt<fgNTrkMod) {
1340 AliInfo(Form("iMt %i, iCh %i, iDe %i, iDeId %i, iDeNb %i, iPar %i",iMt, iDetElemId/100, iDe, iDetElemId, iDetElemNumber, iDetElemNumber*fgNParCh+i));
1341 lDetElemMisAlignment[i] = misAlignments[iDetElemNumber*fgNParCh+i];
1344 // local transformation of this detection element.
1345 TGeoCombiTrans localTransform
1346 = TGeoCombiTrans(*detElement->GetLocalTransformation());
1347 TGeoCombiTrans newLocalTransform = ReAlign(localTransform,lDetElemMisAlignment);
1348 newDetElement->SetLocalTransformation(newLocalTransform);
1350 // global transformation
1351 TGeoHMatrix newGlobalTransform =
1352 AliMUONGeometryBuilder::Multiply(newModuleTransform,
1354 newDetElement->SetGlobalTransformation(newGlobalTransform);
1356 // add this det element to module
1357 newModuleTransformer->GetDetElementStore()->Add(newDetElement->GetId(),
1360 // In the Alice Alignment Framework misalignment objects store
1361 // global delta transformation
1362 // Get detection "intermediate" global transformation
1363 TGeoHMatrix newOldGlobalTransform = newModuleTransform * localTransform;
1364 // Get detection element global delta transformation:
1365 // Tdelta = Tnew * Told.inverse
1366 TGeoHMatrix deltaGlobalTransform
1367 = AliMUONGeometryBuilder::Multiply(newGlobalTransform,
1368 newOldGlobalTransform.Inverse());
1370 // Create mis alignment matrix
1371 newGeometryTransformer
1372 ->AddMisAlignDetElement(detElement->GetId(), deltaGlobalTransform);
1376 AliInfo(Form("Added module transformer %i to the transformer", iMt));
1377 newGeometryTransformer->AddModuleTransformer(newModuleTransformer);
1379 return newGeometryTransformer;
1382 //______________________________________________________________________
1383 void AliMUONAlignment::SetAlignmentResolution(const TClonesArray* misAlignArray, Int_t rChId, Double_t rChResX, Double_t rChResY, Double_t rDeResX, Double_t rDeResY)
1386 /// Set alignment resolution to misalign objects to be stored in CDB
1387 Int_t chIdMin = (rChId<0)? 0 : rChId;
1388 Int_t chIdMax = (rChId<0)? 9 : rChId;
1389 Double_t chResX = rChResX;
1390 Double_t chResY = rChResY;
1391 Double_t deResX = rDeResX;
1392 Double_t deResY = rDeResY;
1394 TMatrixDSym mChCorrMatrix(6);
1395 mChCorrMatrix[0][0]=chResX*chResX;
1396 mChCorrMatrix[1][1]=chResY*chResY;
1397 // mChCorrMatrix.Print();
1399 TMatrixDSym mDECorrMatrix(6);
1400 mDECorrMatrix[0][0]=deResX*deResX;
1401 mDECorrMatrix[1][1]=deResY*deResY;
1402 // mDECorrMatrix.Print();
1404 AliAlignObjMatrix *alignMat = 0x0;
1406 for(Int_t chId=chIdMin; chId<=chIdMax; chId++) {
1410 chName1 = Form("GM%d",chId);
1411 chName2 = Form("GM%d",chId);
1413 chName1 = Form("GM%d",4+(chId-4)*2);
1414 chName2 = Form("GM%d",4+(chId-4)*2+1);
1417 for (int i=0; i<misAlignArray->GetEntries(); i++) {
1418 alignMat = (AliAlignObjMatrix*)misAlignArray->At(i);
1419 TString volName(alignMat->GetSymName());
1420 if((volName.Contains(chName1)&&
1421 ((volName.Last('/')==volName.Index(chName1)+chName1.Length())||
1422 (volName.Length()==volName.Index(chName1)+chName1.Length())))||
1423 (volName.Contains(chName2)&&
1424 ((volName.Last('/')==volName.Index(chName2)+chName2.Length())||
1425 (volName.Length()==volName.Index(chName2)+chName2.Length())))){
1426 volName.Remove(0,volName.Last('/')+1);
1427 if (volName.Contains("GM")) {
1428 // alignMat->Print("NULL");
1429 alignMat->SetCorrMatrix(mChCorrMatrix);
1430 } else if (volName.Contains("DE")) {
1431 // alignMat->Print("NULL");
1432 alignMat->SetCorrMatrix(mDECorrMatrix);