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(TRootIOCtor* /*dummy*/)
128 /// Root IO constructor
130 for (Int_t iCh=0; iCh<10; iCh++) {
131 fChOnOff[iCh] = kFALSE;
133 fSpecLROnOff[0] = kFALSE; fSpecLROnOff[1] = kFALSE;
134 for (Int_t iDoF=0; iDoF<4; iDoF++) {
136 fAllowVar[iDoF] = kFALSE;
138 for (Int_t i=0; i<3; i++) {
145 fTrackSlope0[0]=0; fTrackSlope0[1]=0;
146 fTrackSlope[0]=0; fTrackSlope[1]=0;
147 fMeas[0]=0; fMeas[1]=0;
148 fSigma[0]=0; fSigma[1]=0;
149 for (Int_t iLPar=0; iLPar<4; iLPar++) {
150 fLocalDerivatives[iLPar]=0;
152 for (Int_t iGPar=0; iGPar<624; iGPar++) {
153 fLocalDerivatives[iGPar]=0;
157 //_____________________________________________________________________
158 AliMUONAlignment::~AliMUONAlignment()
161 //_____________________________________________________________________
162 void AliMUONAlignment::Init(
163 Int_t nGlobal, /* number of global paramers */
164 Int_t nLocal, /* number of local parameters */
165 Int_t nStdDev /* std dev cut */ )
167 /// Initialization of AliMillepede. Fix parameters, define constraints ...
168 fMillepede->InitMille(nGlobal,nLocal,nStdDev,fResCut,fResCutInitial);
170 // Bool_t bStOnOff[5] = {kTRUE,kTRUE,kTRUE,kTRUE,kTRUE};
171 // Bool_t bChOnOff[10] = {kTRUE,kTRUE,kTRUE,kTRUE,kTRUE,kTRUE,kTRUE,kTRUE,kTRUE,kTRUE};
172 // Bool_t bSpecLROnOff[2] = {kTRUE,kTRUE};
174 // AllowVariations(bChOnOff);
176 // Fix parameters or add constraints here
177 // for (Int_t iSt=0; iSt<5; iSt++)
178 // { if (!bStOnOff[iSt]) FixStation(iSt+1); }
180 // for (Int_t iCh=0; iCh<10; iCh++)
181 // { if (!bChOnOff[iCh]) FixChamber(iCh+1); }
183 // FixHalfSpectrometer(bChOnOff,bSpecLROnOff);
187 // Define global constrains to be applied
188 // X, Y, P, XvsZ, YvsZ, PvsZ, XvsY, YvsY, PvsY
189 Bool_t bVarXYT[9] = {kTRUE,kTRUE,kTRUE,kTRUE,kTRUE,kTRUE,kTRUE,kTRUE,kTRUE};
190 Bool_t bDetTLBR[4] = {kFALSE,kTRUE,kFALSE,kTRUE};
191 // AddConstraints(bStOnOff,bVarXYT,bDetTLBR,bSpecLROnOff);
193 // Other possible way to add constrains
194 bVarXYT[0] = kFALSE; bVarXYT[1] = kFALSE; bVarXYT[2] = kTRUE;
195 bDetTLBR[0] = kFALSE; bDetTLBR[1] = kTRUE; bDetTLBR[2] = kFALSE; bDetTLBR[3] = kFALSE;
196 // AddConstraints(bStOnOff,bVarXYT,bDetTLBR);
198 bVarXYT[0] = kTRUE; bVarXYT[1] = kTRUE; bVarXYT[2] = kFALSE;
199 // AddConstraints(bStOnOff,bVarXYT);
202 if (fStartFac>1) fMillepede->SetIterations(fStartFac);
205 //_____________________________________________________________________
206 void AliMUONAlignment::FixStation(Int_t iSt)
208 /// Fix all detection elements of station iSt
209 Int_t iDetElemFirst = (iSt>1) ? fgSNDetElemCh[2*(iSt-1)-1] : 0;
210 Int_t iDetElemLast = fgSNDetElemCh[2*(iSt)-1];
211 for (Int_t i = iDetElemFirst; i < iDetElemLast; i++)
213 FixParameter(i*fgNParCh+0, 0.0);
214 FixParameter(i*fgNParCh+1, 0.0);
215 FixParameter(i*fgNParCh+2, 0.0);
216 FixParameter(i*fgNParCh+3, 0.0);
221 //_____________________________________________________________________
222 void AliMUONAlignment::FixChamber(Int_t iCh)
224 /// Fix all detection elements of chamber iCh
225 Int_t iDetElemFirst = (iCh>1) ? fgSNDetElemCh[iCh-2] : 0;
226 Int_t iDetElemLast = fgSNDetElemCh[iCh-1];
227 for (Int_t i = iDetElemFirst; i < iDetElemLast; i++)
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);
236 //_____________________________________________________________________
237 void AliMUONAlignment::FixDetElem(Int_t iDetElemId, TString sVarXYT)
240 /// Fix a given detection element
241 Int_t iDetElemNumber = iDetElemId%100;
242 for (int iCh=0; iCh<iDetElemId/100-1; iCh++)
244 iDetElemNumber += fgNDetElemCh[iCh];
247 if (sVarXYT.Contains("X"))
250 FixParameter(iDetElemNumber*fgNParCh+0, 0.0);
253 if (sVarXYT.Contains("Y"))
256 FixParameter(iDetElemNumber*fgNParCh+1, 0.0);
259 if (sVarXYT.Contains("T"))
262 FixParameter(iDetElemNumber*fgNParCh+2, 0.0);
265 if (sVarXYT.Contains("Z"))
268 FixParameter(iDetElemNumber*fgNParCh+3, 0.0);
273 //_____________________________________________________________________
274 void AliMUONAlignment::FixHalfSpectrometer(const Bool_t *lChOnOff, const Bool_t *lSpecLROnOff)
277 /// Fix left or right detector
278 for (Int_t i = 0; i < fgNDetElem; i++)
282 for (iCh=1; iCh<=fgNCh; iCh++)
283 { if (i<fgSNDetElemCh[iCh-1]) break; }
288 Int_t lDetElemNumber = (iCh==1) ? i : i-fgSNDetElemCh[iCh-2];
289 if (iCh>=1 && iCh<=4)
292 if ((lDetElemNumber==1 || lDetElemNumber==2) && !lSpecLROnOff[0])
294 // From track crossings
295 FixParameter(i*fgNParCh+0, 0.0);
296 FixParameter(i*fgNParCh+1, 0.0);
297 FixParameter(i*fgNParCh+2, 0.0);
298 FixParameter(i*fgNParCh+3, 0.0);
301 if ((lDetElemNumber==0 || lDetElemNumber==3) && !lSpecLROnOff[1])
303 // From track crossings
304 FixParameter(i*fgNParCh+0, 0.0);
305 FixParameter(i*fgNParCh+1, 0.0);
306 FixParameter(i*fgNParCh+2, 0.0);
307 FixParameter(i*fgNParCh+3, 0.0);
312 if (iCh>=5 && iCh<=6)
315 if ((lDetElemNumber>=5&&lDetElemNumber<=13) && !lSpecLROnOff[0])
317 FixParameter(i*fgNParCh+0, 0.0);
318 FixParameter(i*fgNParCh+1, 0.0);
319 FixParameter(i*fgNParCh+2, 0.0);
320 FixParameter(i*fgNParCh+3, 0.0);
323 if (((lDetElemNumber>=0&&lDetElemNumber<=4) ||
324 (lDetElemNumber>=14&&lDetElemNumber<=17)) && !lSpecLROnOff[1])
327 FixParameter(i*fgNParCh+0, 0.0);
328 FixParameter(i*fgNParCh+1, 0.0);
329 FixParameter(i*fgNParCh+2, 0.0);
330 FixParameter(i*fgNParCh+3, 0.0);
335 if (iCh>=7 && iCh<=10)
338 if ((lDetElemNumber>=7&&lDetElemNumber<=19) && !lSpecLROnOff[0])
340 FixParameter(i*fgNParCh+0, 0.0);
341 FixParameter(i*fgNParCh+1, 0.0);
342 FixParameter(i*fgNParCh+2, 0.0);
343 FixParameter(i*fgNParCh+3, 0.0);
346 if (((lDetElemNumber>=0&&lDetElemNumber<=6) ||
347 (lDetElemNumber>=20&&lDetElemNumber<=25)) && !lSpecLROnOff[1])
349 FixParameter(i*fgNParCh+0, 0.0);
350 FixParameter(i*fgNParCh+1, 0.0);
351 FixParameter(i*fgNParCh+2, 0.0);
352 FixParameter(i*fgNParCh+3, 0.0);
363 //______________________________________________________________________
364 void AliMUONAlignment::SetNonLinear(const Bool_t *lChOnOff, const Bool_t *lVarXYT)
367 /// Set non linear parameter flag selected chambers and degrees of freedom
368 for (Int_t i = 0; i < fgNDetElem; i++)
372 for (iCh=1; iCh<=fgNCh; iCh++)
373 { if (i<fgSNDetElemCh[iCh-1]) break; }
381 SetNonLinear(i*fgNParCh+0);
387 SetNonLinear(i*fgNParCh+1);
393 SetNonLinear(i*fgNParCh+2);
399 SetNonLinear(i*fgNParCh+3);
408 //______________________________________________________________________
409 void AliMUONAlignment::AddConstraints(const Bool_t *lChOnOff, const Bool_t *lVarXYT)
412 /// Add constraint equations for selected chambers and degrees of freedom
413 for (Int_t i = 0; i < fgNDetElem; i++)
417 for (iCh=1; iCh<=fgNCh; iCh++)
419 if (i<fgSNDetElemCh[iCh-1]) break;
427 fConstraintX[i*fgNParCh+0]=1.0;
433 fConstraintY[i*fgNParCh+1]=1.0;
439 fConstraintP[i*fgNParCh+2]=1.0;
441 // if (lVarXYT[3]) { // Z constraint
442 // fConstraintP[i*fgNParCh+3]=1.0;
451 AddConstraint(fConstraintX,0.0);
457 AddConstraint(fConstraintY,0.0);
463 AddConstraint(fConstraintP,0.0);
466 // if (lVarXYT[3]) { // Z constraint
467 // AddConstraint(fConstraintP,0.0);
471 //______________________________________________________________________
472 void AliMUONAlignment::AddConstraints(const Bool_t *lChOnOff, const Bool_t *lVarXYT, const Bool_t *lDetTLBR, const Bool_t *lSpecLROnOff)
474 /// Add constraint equations for selected chambers, degrees of freedom and detector half
475 Double_t lDetElemLocX = 0.;
476 Double_t lDetElemLocY = 0.;
477 Double_t lDetElemLocZ = 0.;
478 Double_t lDetElemGloX = 0.;
479 Double_t lDetElemGloY = 0.;
480 Double_t lDetElemGloZ = 0.;
481 Double_t lMeanY = 0.;
482 Double_t lSigmaY = 0.;
483 Double_t lMeanZ = 0.;
484 Double_t lSigmaZ = 0.;
486 for (Int_t i = 0; i < fgNDetElem; i++)
490 for (iCh=1; iCh<=fgNCh; iCh++){
491 if (i<fgSNDetElemCh[iCh-1]) break;
493 if (lChOnOff[iCh-1]){
494 Int_t lDetElemNumber = (iCh==1) ? i : i-fgSNDetElemCh[iCh-2];
495 Int_t lDetElemId = iCh*100+lDetElemNumber;
496 fTransform->Local2Global(lDetElemId,lDetElemLocX,lDetElemLocY,lDetElemLocZ,
497 lDetElemGloX,lDetElemGloY,lDetElemGloZ);
498 if (iCh>=1 && iCh<=4){
499 if ((lDetElemNumber==1 || lDetElemNumber==2) && lSpecLROnOff[0]){ // From track crossings
500 lMeanY += lDetElemGloY;
501 lSigmaY += lDetElemGloY*lDetElemGloY;
502 lMeanZ += lDetElemGloZ;
503 lSigmaZ += lDetElemGloZ*lDetElemGloZ;
506 if ((lDetElemNumber==0 || lDetElemNumber==3) && lSpecLROnOff[1]){ // From track crossings
507 lMeanY += lDetElemGloY;
508 lSigmaY += lDetElemGloY*lDetElemGloY;
509 lMeanZ += lDetElemGloZ;
510 lSigmaZ += lDetElemGloZ*lDetElemGloZ;
514 if (iCh>=5 && iCh<=6){
515 if ((lDetElemNumber>=5&&lDetElemNumber<=13) && lSpecLROnOff[0]){
516 lMeanY += lDetElemGloY;
517 lSigmaY += lDetElemGloY*lDetElemGloY;
518 lMeanZ += lDetElemGloZ;
519 lSigmaZ += lDetElemGloZ*lDetElemGloZ;
522 if (((lDetElemNumber>=0&&lDetElemNumber<=4) ||
523 (lDetElemNumber>=14&&lDetElemNumber<=17)) && lSpecLROnOff[1]){
524 lMeanY += lDetElemGloY;
525 lSigmaY += lDetElemGloY*lDetElemGloY;
526 lMeanZ += lDetElemGloZ;
527 lSigmaZ += lDetElemGloZ*lDetElemGloZ;
531 if (iCh>=7 && iCh<=10){
532 if ((lDetElemNumber>=7&&lDetElemNumber<=19) && lSpecLROnOff[0]){
533 lMeanY += lDetElemGloY;
534 lSigmaY += lDetElemGloY*lDetElemGloY;
535 lMeanZ += lDetElemGloZ;
536 lSigmaZ += lDetElemGloZ*lDetElemGloZ;
539 if (((lDetElemNumber>=0&&lDetElemNumber<=6) ||
540 (lDetElemNumber>=20&&lDetElemNumber<=25)) && lSpecLROnOff[1]){
541 lMeanY += lDetElemGloY;
542 lSigmaY += lDetElemGloY*lDetElemGloY;
543 lMeanZ += lDetElemGloZ;
544 lSigmaZ += lDetElemGloZ*lDetElemGloZ;
552 lSigmaY /= lNDetElem;
553 lSigmaY = TMath::Sqrt(lSigmaY-lMeanY*lMeanY);
555 lSigmaZ /= lNDetElem;
556 lSigmaZ = TMath::Sqrt(lSigmaZ-lMeanZ*lMeanZ);
557 AliInfo(Form("Used %i DetElem, MeanZ= %f , SigmaZ= %f", lNDetElem,lMeanZ,lSigmaZ));
559 AliError("No detection elements to constrain!!!");
563 for (Int_t i = 0; i < fgNDetElem; i++){
565 for (iCh=1; iCh<=fgNCh; iCh++){
566 if (i<fgSNDetElemCh[iCh-1]) break;
568 if (lChOnOff[iCh-1]){
569 Int_t lDetElemNumber = (iCh==1) ? i : i-fgSNDetElemCh[iCh-2];
570 Int_t lDetElemId = iCh*100+lDetElemNumber;
571 fTransform->Local2Global(lDetElemId,lDetElemLocX,lDetElemLocY,lDetElemLocZ,
572 lDetElemGloX,lDetElemGloY,lDetElemGloZ);
573 if (lVarXYT[0]) { // X constraint
574 if (lDetTLBR[0]) ConstrainT(i,iCh,fConstraintXT,0); // Top half
575 if (lDetTLBR[1]) ConstrainL(i,iCh,fConstraintXL,0); // Left half
576 if (lDetTLBR[2]) ConstrainB(i,iCh,fConstraintXB,0); // Bottom half
577 if (lDetTLBR[3]) ConstrainR(i,iCh,fConstraintXR,0); // Right half
579 if (lVarXYT[1]) { // Y constraint
580 if (lDetTLBR[0]) ConstrainT(i,iCh,fConstraintYT,1); // Top half
581 if (lDetTLBR[1]) ConstrainL(i,iCh,fConstraintYL,1); // Left half
582 if (lDetTLBR[2]) ConstrainB(i,iCh,fConstraintYB,1); // Bottom half
583 if (lDetTLBR[3]) ConstrainR(i,iCh,fConstraintYR,1); // Right half
585 if (lVarXYT[2]) { // P constraint
586 if (lDetTLBR[0]) ConstrainT(i,iCh,fConstraintPT,2); // Top half
587 if (lDetTLBR[1]) ConstrainL(i,iCh,fConstraintPL,2); // Left half
588 if (lDetTLBR[2]) ConstrainB(i,iCh,fConstraintPB,2); // Bottom half
589 if (lDetTLBR[3]) ConstrainR(i,iCh,fConstraintPR,2); // Right half
591 if (lVarXYT[3]) { // X-Z shearing
592 if (lDetTLBR[0]) ConstrainT(i,iCh,fConstraintXZT,0,(lDetElemGloZ-lMeanZ)/lSigmaZ); // Top half
593 if (lDetTLBR[1]) ConstrainL(i,iCh,fConstraintXZL,0,(lDetElemGloZ-lMeanZ)/lSigmaZ); // Left half
594 if (lDetTLBR[2]) ConstrainB(i,iCh,fConstraintXZB,0,(lDetElemGloZ-lMeanZ)/lSigmaZ); // Bottom half
595 if (lDetTLBR[3]) ConstrainR(i,iCh,fConstraintXZR,0,(lDetElemGloZ-lMeanZ)/lSigmaZ); // Right half
597 if (lVarXYT[4]) { // Y-Z shearing
598 if (lDetTLBR[0]) ConstrainT(i,iCh,fConstraintYZT,1,(lDetElemGloZ-lMeanZ)/lSigmaZ); // Top half
599 if (lDetTLBR[1]) ConstrainL(i,iCh,fConstraintYZL,1,(lDetElemGloZ-lMeanZ)/lSigmaZ); // Left half
600 if (lDetTLBR[2]) ConstrainB(i,iCh,fConstraintYZB,1,(lDetElemGloZ-lMeanZ)/lSigmaZ); // Bottom half
601 if (lDetTLBR[3]) ConstrainR(i,iCh,fConstraintYZR,1,(lDetElemGloZ-lMeanZ)/lSigmaZ); // Right half
603 if (lVarXYT[5]) { // P-Z rotation
604 if (lDetTLBR[0]) ConstrainT(i,iCh,fConstraintPZT,2,(lDetElemGloZ-lMeanZ)/lSigmaZ); // Top half
605 if (lDetTLBR[1]) ConstrainL(i,iCh,fConstraintPZL,2,(lDetElemGloZ-lMeanZ)/lSigmaZ); // Left half
606 if (lDetTLBR[2]) ConstrainB(i,iCh,fConstraintPZB,2,(lDetElemGloZ-lMeanZ)/lSigmaZ); // Bottom half
607 if (lDetTLBR[3]) ConstrainR(i,iCh,fConstraintPZR,2,(lDetElemGloZ-lMeanZ)/lSigmaZ); // Right half
609 if (lVarXYT[6]) { // X-Y shearing
610 if (lDetTLBR[0]) ConstrainT(i,iCh,fConstraintXYT,0,(lDetElemGloY-lMeanY)/lSigmaY); // Top half
611 if (lDetTLBR[1]) ConstrainL(i,iCh,fConstraintXYL,0,(lDetElemGloY-lMeanY)/lSigmaY); // Left half
612 if (lDetTLBR[2]) ConstrainB(i,iCh,fConstraintXYB,0,(lDetElemGloY-lMeanY)/lSigmaY); // Bottom half
613 if (lDetTLBR[3]) ConstrainR(i,iCh,fConstraintXYR,0,(lDetElemGloY-lMeanY)/lSigmaY); // Right half
615 if (lVarXYT[7]) { // Y-Y scaling
616 if (lDetTLBR[0]) ConstrainT(i,iCh,fConstraintYYT,1,(lDetElemGloY-lMeanY)/lSigmaY); // Top half
617 if (lDetTLBR[1]) ConstrainL(i,iCh,fConstraintYYL,1,(lDetElemGloY-lMeanY)/lSigmaY); // Left half
618 if (lDetTLBR[2]) ConstrainB(i,iCh,fConstraintYYB,1,(lDetElemGloY-lMeanY)/lSigmaY); // Bottom half
619 if (lDetTLBR[3]) ConstrainR(i,iCh,fConstraintYYR,1,(lDetElemGloY-lMeanY)/lSigmaY); // Right half
621 if (lVarXYT[8]) { // P-Y rotation
622 if (lDetTLBR[0]) ConstrainT(i,iCh,fConstraintPYT,2,(lDetElemGloY-lMeanY)/lSigmaY); // Top half
623 if (lDetTLBR[1]) ConstrainL(i,iCh,fConstraintPYL,2,(lDetElemGloY-lMeanY)/lSigmaY); // Left half
624 if (lDetTLBR[2]) ConstrainB(i,iCh,fConstraintPYB,2,(lDetElemGloY-lMeanY)/lSigmaY); // Bottom half
625 if (lDetTLBR[3]) ConstrainR(i,iCh,fConstraintPYR,2,(lDetElemGloY-lMeanY)/lSigmaY); // Right half
629 if (lVarXYT[0]) { // X constraint
630 if (lDetTLBR[0]) AddConstraint(fConstraintXT,0.0); // Top half
631 if (lDetTLBR[1]) AddConstraint(fConstraintXL,0.0); // Left half
632 if (lDetTLBR[2]) AddConstraint(fConstraintXB,0.0); // Bottom half
633 if (lDetTLBR[3]) AddConstraint(fConstraintXR,0.0); // Right half
635 if (lVarXYT[1]) { // Y constraint
636 if (lDetTLBR[0]) AddConstraint(fConstraintYT,0.0); // Top half
637 if (lDetTLBR[1]) AddConstraint(fConstraintYL,0.0); // Left half
638 if (lDetTLBR[2]) AddConstraint(fConstraintYB,0.0); // Bottom half
639 if (lDetTLBR[3]) AddConstraint(fConstraintYR,0.0); // Right half
641 if (lVarXYT[2]) { // T constraint
642 if (lDetTLBR[0]) AddConstraint(fConstraintPT,0.0); // Top half
643 if (lDetTLBR[1]) AddConstraint(fConstraintPL,0.0); // Left half
644 if (lDetTLBR[2]) AddConstraint(fConstraintPB,0.0); // Bottom half
645 if (lDetTLBR[3]) AddConstraint(fConstraintPR,0.0); // Right half
647 if (lVarXYT[3]) { // X-Z constraint
648 if (lDetTLBR[0]) AddConstraint(fConstraintXZT,0.0); // Top half
649 if (lDetTLBR[1]) AddConstraint(fConstraintXZL,0.0); // Left half
650 if (lDetTLBR[2]) AddConstraint(fConstraintXZB,0.0); // Bottom half
651 if (lDetTLBR[3]) AddConstraint(fConstraintXZR,0.0); // Right half
653 if (lVarXYT[4]) { // Y-Z constraint
654 if (lDetTLBR[0]) AddConstraint(fConstraintYZT,0.0); // Top half
655 if (lDetTLBR[1]) AddConstraint(fConstraintYZL,0.0); // Left half
656 if (lDetTLBR[2]) AddConstraint(fConstraintYZB,0.0); // Bottom half
657 if (lDetTLBR[3]) AddConstraint(fConstraintYZR,0.0); // Right half
659 if (lVarXYT[5]) { // P-Z constraint
660 if (lDetTLBR[0]) AddConstraint(fConstraintPZT,0.0); // Top half
661 if (lDetTLBR[1]) AddConstraint(fConstraintPZL,0.0); // Left half
662 if (lDetTLBR[2]) AddConstraint(fConstraintPZB,0.0); // Bottom half
663 if (lDetTLBR[3]) AddConstraint(fConstraintPZR,0.0); // Right half
665 if (lVarXYT[6]) { // X-Y constraint
666 if (lDetTLBR[0]) AddConstraint(fConstraintXYT,0.0); // Top half
667 if (lDetTLBR[1]) AddConstraint(fConstraintXYL,0.0); // Left half
668 if (lDetTLBR[2]) AddConstraint(fConstraintXYB,0.0); // Bottom half
669 if (lDetTLBR[3]) AddConstraint(fConstraintXYR,0.0); // Right half
671 if (lVarXYT[7]) { // Y-Y constraint
672 if (lDetTLBR[0]) AddConstraint(fConstraintYYT,0.0); // Top half
673 if (lDetTLBR[1]) AddConstraint(fConstraintYYL,0.0); // Left half
674 if (lDetTLBR[2]) AddConstraint(fConstraintYYB,0.0); // Bottom half
675 if (lDetTLBR[3]) AddConstraint(fConstraintYYR,0.0); // Right half
677 if (lVarXYT[8]) { // P-Y constraint
678 if (lDetTLBR[0]) AddConstraint(fConstraintPYT,0.0); // Top half
679 if (lDetTLBR[1]) AddConstraint(fConstraintPYL,0.0); // Left half
680 if (lDetTLBR[2]) AddConstraint(fConstraintPYB,0.0); // Bottom half
681 if (lDetTLBR[3]) AddConstraint(fConstraintPYR,0.0); // Right half
685 void AliMUONAlignment::ConstrainT(Int_t lDetElem, Int_t lCh, Double_t *lConstraintT, Int_t iVar, Double_t /*lWeight*/) const{
686 /// Set constrain equation for top half of spectrometer
687 Int_t lDetElemNumber = (lCh==1) ? lDetElem : lDetElem-fgSNDetElemCh[lCh-2];
688 if (lCh>=1 && lCh<=4){
689 if (lDetElemNumber==0 || lDetElemNumber==1){ // From track crossings
690 lConstraintT[lDetElem*fgNParCh+iVar]=1.0;
693 if (lCh>=5 && lCh<=6){
694 if (lDetElemNumber>=0&&lDetElemNumber<=9){
695 lConstraintT[lDetElem*fgNParCh+iVar]=1.0;
698 if (lCh>=7 && lCh<=10){
699 if (lDetElemNumber>=0&&lDetElemNumber<=13){
700 lConstraintT[lDetElem*fgNParCh+iVar]=1.0;
705 //______________________________________________________________________
706 void AliMUONAlignment::ConstrainL(Int_t lDetElem, Int_t lCh, Double_t *lConstraintL, Int_t iVar, Double_t lWeight) const{
707 /// Set constrain equation for left half of spectrometer
708 Int_t lDetElemNumber = (lCh==1) ? lDetElem : lDetElem-fgSNDetElemCh[lCh-2];
709 if (lCh>=1 && lCh<=4){
710 if (lDetElemNumber==1 || lDetElemNumber==2){ // From track crossings
711 lConstraintL[lDetElem*fgNParCh+iVar]=lWeight;
714 if (lCh>=5 && lCh<=6){
715 if (lDetElemNumber>=5&&lDetElemNumber<=13){
716 lConstraintL[lDetElem*fgNParCh+iVar]=lWeight;
719 if (lCh>=7 && lCh<=10){
720 if (lDetElemNumber>=7&&lDetElemNumber<=19){
721 lConstraintL[lDetElem*fgNParCh+iVar]=lWeight;
726 //______________________________________________________________________
727 void AliMUONAlignment::ConstrainB(Int_t lDetElem, Int_t lCh, Double_t *lConstraintB, Int_t iVar, Double_t /*lWeight*/) const{
728 /// Set constrain equation for bottom half of spectrometer
729 Int_t lDetElemNumber = (lCh==1) ? lDetElem : lDetElem-fgSNDetElemCh[lCh-2];
730 if (lCh>=1 && lCh<=4){
731 if (lDetElemNumber==2 || lDetElemNumber==3){ // From track crossings
732 lConstraintB[lDetElem*fgNParCh+iVar]=1.0;
735 if (lCh>=5 && lCh<=6){
736 if ((lDetElemNumber>=9&&lDetElemNumber<=17) ||
737 (lDetElemNumber==0)){
738 lConstraintB[lDetElem*fgNParCh+iVar]=1.0;
741 if (lCh>=7 && lCh<=10){
742 if ((lDetElemNumber>=13&&lDetElemNumber<=25) ||
743 (lDetElemNumber==0)){
744 lConstraintB[lDetElem*fgNParCh+iVar]=1.0;
749 //______________________________________________________________________
750 void AliMUONAlignment::ConstrainR(Int_t lDetElem, Int_t lCh, Double_t *lConstraintR, Int_t iVar, Double_t lWeight) const{
751 /// Set constrain equation for right half of spectrometer
752 Int_t lDetElemNumber = (lCh==1) ? lDetElem : lDetElem-fgSNDetElemCh[lCh-2];
753 if (lCh>=1 && lCh<=4){
754 if (lDetElemNumber==0 || lDetElemNumber==3){ // From track crossings
755 lConstraintR[lDetElem*fgNParCh+iVar]=lWeight;
758 if (lCh>=5 && lCh<=6){
759 if ((lDetElemNumber>=0&&lDetElemNumber<=4) ||
760 (lDetElemNumber>=14&&lDetElemNumber<=17)){
761 lConstraintR[lDetElem*fgNParCh+iVar]=lWeight;
764 if (lCh>=7 && lCh<=10){
765 if ((lDetElemNumber>=0&&lDetElemNumber<=6) ||
766 (lDetElemNumber>=20&&lDetElemNumber<=25)){
767 lConstraintR[lDetElem*fgNParCh+iVar]=lWeight;
772 //______________________________________________________________________
773 void AliMUONAlignment::ResetConstraints(){
774 /// Reset all constraint equations
775 for (Int_t i = 0; i < fgNDetElem; i++){
776 fConstraintX[i*fgNParCh+0]=0.0;
777 fConstraintX[i*fgNParCh+1]=0.0;
778 fConstraintX[i*fgNParCh+2]=0.0;
779 fConstraintY[i*fgNParCh+0]=0.0;
780 fConstraintY[i*fgNParCh+1]=0.0;
781 fConstraintY[i*fgNParCh+2]=0.0;
782 fConstraintP[i*fgNParCh+0]=0.0;
783 fConstraintP[i*fgNParCh+1]=0.0;
784 fConstraintP[i*fgNParCh+2]=0.0;
785 fConstraintXT[i*fgNParCh+0]=0.0;
786 fConstraintXT[i*fgNParCh+1]=0.0;
787 fConstraintXT[i*fgNParCh+2]=0.0;
788 fConstraintYT[i*fgNParCh+0]=0.0;
789 fConstraintYT[i*fgNParCh+1]=0.0;
790 fConstraintYT[i*fgNParCh+2]=0.0;
791 fConstraintPT[i*fgNParCh+0]=0.0;
792 fConstraintPT[i*fgNParCh+1]=0.0;
793 fConstraintPT[i*fgNParCh+2]=0.0;
794 fConstraintXZT[i*fgNParCh+0]=0.0;
795 fConstraintXZT[i*fgNParCh+1]=0.0;
796 fConstraintXZT[i*fgNParCh+2]=0.0;
797 fConstraintYZT[i*fgNParCh+0]=0.0;
798 fConstraintYZT[i*fgNParCh+1]=0.0;
799 fConstraintYZT[i*fgNParCh+2]=0.0;
800 fConstraintPZT[i*fgNParCh+0]=0.0;
801 fConstraintPZT[i*fgNParCh+1]=0.0;
802 fConstraintPZT[i*fgNParCh+2]=0.0;
803 fConstraintXYT[i*fgNParCh+0]=0.0;
804 fConstraintXYT[i*fgNParCh+1]=0.0;
805 fConstraintXYT[i*fgNParCh+2]=0.0;
806 fConstraintYYT[i*fgNParCh+0]=0.0;
807 fConstraintYYT[i*fgNParCh+1]=0.0;
808 fConstraintYYT[i*fgNParCh+2]=0.0;
809 fConstraintPYT[i*fgNParCh+0]=0.0;
810 fConstraintPYT[i*fgNParCh+1]=0.0;
811 fConstraintPYT[i*fgNParCh+2]=0.0;
812 fConstraintXL[i*fgNParCh+0]=0.0;
813 fConstraintXL[i*fgNParCh+1]=0.0;
814 fConstraintXL[i*fgNParCh+2]=0.0;
815 fConstraintYL[i*fgNParCh+0]=0.0;
816 fConstraintYL[i*fgNParCh+1]=0.0;
817 fConstraintYL[i*fgNParCh+2]=0.0;
818 fConstraintPL[i*fgNParCh+0]=0.0;
819 fConstraintPL[i*fgNParCh+1]=0.0;
820 fConstraintPL[i*fgNParCh+2]=0.0;
821 fConstraintXZL[i*fgNParCh+0]=0.0;
822 fConstraintXZL[i*fgNParCh+1]=0.0;
823 fConstraintXZL[i*fgNParCh+2]=0.0;
824 fConstraintYZL[i*fgNParCh+0]=0.0;
825 fConstraintYZL[i*fgNParCh+1]=0.0;
826 fConstraintYZL[i*fgNParCh+2]=0.0;
827 fConstraintPZL[i*fgNParCh+0]=0.0;
828 fConstraintPZL[i*fgNParCh+1]=0.0;
829 fConstraintPZL[i*fgNParCh+2]=0.0;
830 fConstraintXYL[i*fgNParCh+0]=0.0;
831 fConstraintXYL[i*fgNParCh+1]=0.0;
832 fConstraintXYL[i*fgNParCh+2]=0.0;
833 fConstraintYYL[i*fgNParCh+0]=0.0;
834 fConstraintYYL[i*fgNParCh+1]=0.0;
835 fConstraintYYL[i*fgNParCh+2]=0.0;
836 fConstraintPYL[i*fgNParCh+0]=0.0;
837 fConstraintPYL[i*fgNParCh+1]=0.0;
838 fConstraintPYL[i*fgNParCh+2]=0.0;
839 fConstraintXB[i*fgNParCh+0]=0.0;
840 fConstraintXB[i*fgNParCh+1]=0.0;
841 fConstraintXB[i*fgNParCh+2]=0.0;
842 fConstraintYB[i*fgNParCh+0]=0.0;
843 fConstraintYB[i*fgNParCh+1]=0.0;
844 fConstraintYB[i*fgNParCh+2]=0.0;
845 fConstraintPB[i*fgNParCh+0]=0.0;
846 fConstraintPB[i*fgNParCh+1]=0.0;
847 fConstraintPB[i*fgNParCh+2]=0.0;
848 fConstraintXZB[i*fgNParCh+0]=0.0;
849 fConstraintXZB[i*fgNParCh+1]=0.0;
850 fConstraintXZB[i*fgNParCh+2]=0.0;
851 fConstraintYZB[i*fgNParCh+0]=0.0;
852 fConstraintYZB[i*fgNParCh+1]=0.0;
853 fConstraintYZB[i*fgNParCh+2]=0.0;
854 fConstraintPZB[i*fgNParCh+0]=0.0;
855 fConstraintPZB[i*fgNParCh+1]=0.0;
856 fConstraintPZB[i*fgNParCh+2]=0.0;
857 fConstraintXYB[i*fgNParCh+0]=0.0;
858 fConstraintXYB[i*fgNParCh+1]=0.0;
859 fConstraintXYB[i*fgNParCh+2]=0.0;
860 fConstraintYYB[i*fgNParCh+0]=0.0;
861 fConstraintYYB[i*fgNParCh+1]=0.0;
862 fConstraintYYB[i*fgNParCh+2]=0.0;
863 fConstraintPYB[i*fgNParCh+0]=0.0;
864 fConstraintPYB[i*fgNParCh+1]=0.0;
865 fConstraintPYB[i*fgNParCh+2]=0.0;
866 fConstraintXR[i*fgNParCh+0]=0.0;
867 fConstraintXR[i*fgNParCh+1]=0.0;
868 fConstraintXR[i*fgNParCh+2]=0.0;
869 fConstraintYR[i*fgNParCh+0]=0.0;
870 fConstraintYR[i*fgNParCh+1]=0.0;
871 fConstraintYR[i*fgNParCh+2]=0.0;
872 fConstraintPR[i*fgNParCh+0]=0.0;
873 fConstraintPR[i*fgNParCh+1]=0.0;
874 fConstraintPR[i*fgNParCh+2]=0.0;
875 fConstraintXZR[i*fgNParCh+0]=0.0;
876 fConstraintXZR[i*fgNParCh+1]=0.0;
877 fConstraintXZR[i*fgNParCh+2]=0.0;
878 fConstraintYZR[i*fgNParCh+0]=0.0;
879 fConstraintYZR[i*fgNParCh+1]=0.0;
880 fConstraintYZR[i*fgNParCh+2]=0.0;
881 fConstraintPZR[i*fgNParCh+0]=0.0;
882 fConstraintPZR[i*fgNParCh+1]=0.0;
883 fConstraintPZR[i*fgNParCh+2]=0.0;
884 fConstraintPZR[i*fgNParCh+0]=0.0;
885 fConstraintPZR[i*fgNParCh+1]=0.0;
886 fConstraintPZR[i*fgNParCh+2]=0.0;
887 fConstraintXYR[i*fgNParCh+0]=0.0;
888 fConstraintXYR[i*fgNParCh+1]=0.0;
889 fConstraintXYR[i*fgNParCh+2]=0.0;
890 fConstraintYYR[i*fgNParCh+0]=0.0;
891 fConstraintYYR[i*fgNParCh+1]=0.0;
892 fConstraintYYR[i*fgNParCh+2]=0.0;
893 fConstraintPYR[i*fgNParCh+0]=0.0;
894 fConstraintPYR[i*fgNParCh+1]=0.0;
895 fConstraintPYR[i*fgNParCh+2]=0.0;
899 //______________________________________________________________________
900 void AliMUONAlignment::AddConstraint(Double_t *par, Double_t value) {
901 /// Constrain equation defined by par to value
902 fMillepede->SetGlobalConstraint(par, value);
903 AliInfo("Adding constraint");
906 //______________________________________________________________________
907 void AliMUONAlignment::InitGlobalParameters(Double_t *par) {
908 /// Initialize global parameters with par array
909 fMillepede->SetGlobalParameters(par);
910 AliInfo("Init Global Parameters");
913 //______________________________________________________________________
914 void AliMUONAlignment::FixParameter(Int_t iPar, Double_t value) {
915 /// Parameter iPar is encourage to vary in [-value;value].
916 /// If value == 0, parameter is fixed
917 fMillepede->SetParSigma(iPar, value);
918 if (TMath::Abs(value)<1e-4) AliInfo(Form("Parameter %i Fixed", iPar));
921 //______________________________________________________________________
922 void AliMUONAlignment::ResetLocalEquation()
924 /// Reset the derivative vectors
925 for(int i=0; i<fNLocal; i++) {
926 fLocalDerivatives[i] = 0.0;
928 for(int i=0; i<fNGlobal; i++) {
929 fGlobalDerivatives[i] = 0.0;
933 //______________________________________________________________________
934 void AliMUONAlignment::AllowVariations(const Bool_t *bChOnOff)
937 /// Set allowed variation for selected chambers based on fDoF and fAllowVar
938 for (Int_t iCh=1; iCh<=10; iCh++)
943 Int_t iDetElemFirst = (iCh>1) ? fgSNDetElemCh[iCh-2] : 0;
944 Int_t iDetElemLast = fgSNDetElemCh[iCh-1];
945 for (int i=0; i<fgNParCh; i++)
947 AliDebug(1,Form("fDoF[%d]= %d",i,fDoF[i]));
950 for (Int_t j=iDetElemFirst; j<iDetElemLast; j++){
951 FixParameter(j*fgNParCh+i, fAllowVar[i]);
963 //______________________________________________________________________
964 void AliMUONAlignment::SetNonLinear(Int_t iPar /* set non linear flag */ )
966 /// Set nonlinear flag for parameter iPar
967 fMillepede->SetNonLinear(iPar);
968 AliInfo(Form("Parameter %i set to non linear", iPar));
971 //______________________________________________________________________
972 void AliMUONAlignment::SetSigmaXY(Double_t sigmaX, Double_t sigmaY)
975 /// Set expected measurement resolution
976 fSigma[0] = sigmaX; fSigma[1] = sigmaY;
977 AliInfo(Form("Using fSigma[0]=%f and fSigma[1]=%f",fSigma[0],fSigma[1]));
981 //______________________________________________________________________
982 void AliMUONAlignment::LocalEquationX( Bool_t doAlignment )
986 // local cluster record
987 AliMUONAlignmentClusterRecord clusterRecord;
989 // store detector and measurement
990 clusterRecord.SetDetElemId( fDetElemId );
991 clusterRecord.SetDetElemNumber( fDetElemNumber );
992 clusterRecord.SetMeas( fMeas[0] );
993 clusterRecord.SetSigma( fSigma[0] );
995 // store local derivatives
996 clusterRecord.SetLocalDerivative( 0, fCosPhi );
997 clusterRecord.SetLocalDerivative( 1, fCosPhi*(fTrackPos[2] - fTrackPos0[2]) );
998 clusterRecord.SetLocalDerivative( 2, fSinPhi );
999 clusterRecord.SetLocalDerivative( 3, fSinPhi*(fTrackPos[2] - fTrackPos0[2]) );
1001 // store global derivatives
1002 clusterRecord.SetGlobalDerivative( 0, -fCosPhi );
1003 clusterRecord.SetGlobalDerivative( 1, -fSinPhi );
1008 clusterRecord.SetGlobalDerivative(
1010 -fSinPhi*(fTrackPos[0]-fDetElemPos[0])
1011 +fCosPhi*(fTrackPos[1]-fDetElemPos[1]) );
1015 clusterRecord.SetGlobalDerivative(
1017 -fSinPhi*(fTrackPos0[0]+fTrackSlope0[0]*(fTrackPos[2]-fTrackPos0[2])-fDetElemPos[0])
1018 +fCosPhi*(fTrackPos0[1]+fTrackSlope0[1]*(fTrackPos[2]-fTrackPos0[2])-fDetElemPos[1]) );
1022 clusterRecord.SetGlobalDerivative( 3, fCosPhi*fTrackSlope0[0]+fSinPhi*fTrackSlope0[1] );
1024 // append to trackRecord
1025 fTrackRecord.AddClusterRecord( clusterRecord );
1027 // store local equation
1028 if( doAlignment ) LocalEquation( clusterRecord );
1032 //______________________________________________________________________
1033 void AliMUONAlignment::LocalEquationY(Bool_t doAlignment )
1036 // local cluster record
1037 AliMUONAlignmentClusterRecord clusterRecord;
1039 // store detector and measurement
1040 clusterRecord.SetDetElemId( fDetElemId );
1041 clusterRecord.SetDetElemNumber( fDetElemNumber );
1042 clusterRecord.SetMeas( fMeas[1] );
1043 clusterRecord.SetSigma( fSigma[1] );
1045 // store local derivatives
1046 clusterRecord.SetLocalDerivative( 0, -fSinPhi );
1047 clusterRecord.SetLocalDerivative( 1, -fSinPhi*(fTrackPos[2] - fTrackPos0[2] ) );
1048 clusterRecord.SetLocalDerivative( 2, fCosPhi );
1049 clusterRecord.SetLocalDerivative( 3, fCosPhi*(fTrackPos[2] - fTrackPos0[2] ) );
1051 // set global derivatives
1052 clusterRecord.SetGlobalDerivative( 0, fSinPhi);
1053 clusterRecord.SetGlobalDerivative( 1, -fCosPhi);
1058 clusterRecord.SetGlobalDerivative(
1060 -fCosPhi*(fTrackPos[0]-fDetElemPos[0])
1061 -fSinPhi*(fTrackPos[1]-fDetElemPos[1]));
1065 clusterRecord.SetGlobalDerivative(
1067 -fCosPhi*(fTrackPos0[0]+fTrackSlope0[0]*(fTrackPos[2]-fTrackPos0[2])-fDetElemPos[0])
1068 -fSinPhi*(fTrackPos0[1]+fTrackSlope0[1]*(fTrackPos[2]-fTrackPos0[2])-fDetElemPos[1]));
1071 clusterRecord.SetGlobalDerivative( 3, -fSinPhi*fTrackSlope0[0]+fCosPhi*fTrackSlope0[1]);
1073 // append to trackRecord
1074 fTrackRecord.AddClusterRecord( clusterRecord );
1076 // store local equation
1077 if( doAlignment ) LocalEquation( clusterRecord );
1081 //_____________________________________________________
1082 void AliMUONAlignment::LocalEquation( const AliMUONAlignmentClusterRecord& clusterRecord )
1085 // copy to local variables
1086 for( Int_t index = 0; index < 4; ++index )
1088 SetLocalDerivative( index, clusterRecord.GetLocalDerivative( index ) );
1089 SetGlobalDerivative( clusterRecord.GetDetElemNumber()*fgNParCh + index, clusterRecord.GetGlobalDerivative( index ) );
1092 // pass equation parameters to millepede
1093 fMillepede->SetLocalEquation( fGlobalDerivatives, fLocalDerivatives, clusterRecord.GetMeas(), clusterRecord.GetSigma() );
1097 //______________________________________________________________________
1098 void AliMUONAlignment::FillRecPointData()
1101 /// Get information of current cluster
1102 fClustPos[0] = fCluster->GetX();
1103 fClustPos[1] = fCluster->GetY();
1104 fClustPos[2] = fCluster->GetZ();
1105 fTransform->Global2Local(
1106 fDetElemId,fClustPos[0],fClustPos[1],fClustPos[2],
1107 fClustPosLoc[0],fClustPosLoc[1],fClustPosLoc[2]);
1111 //______________________________________________________________________
1112 void AliMUONAlignment::FillTrackParamData()
1115 /// Get information of current track at current cluster
1116 fTrackPos[0] = fTrackParam->GetNonBendingCoor();
1117 fTrackPos[1] = fTrackParam->GetBendingCoor();
1118 fTrackPos[2] = fTrackParam->GetZ();
1119 fTrackSlope[0] = fTrackParam->GetNonBendingSlope();
1120 fTrackSlope[1] = fTrackParam->GetBendingSlope();
1121 fTransform->Global2Local(
1122 fDetElemId,fTrackPos[0],fTrackPos[1],fTrackPos[2],
1123 fTrackPosLoc[0],fTrackPosLoc[1],fTrackPosLoc[2]);
1127 //_____________________________________________________
1128 void AliMUONAlignment::FillDetElemData()
1131 /// Get information of current detection element
1132 Double_t lDetElemLocX = 0.;
1133 Double_t lDetElemLocY = 0.;
1134 Double_t lDetElemLocZ = 0.;
1135 fDetElemId = fCluster->GetDetElemId();
1136 fDetElemNumber = fDetElemId%100;
1137 for (int iCh=0; iCh<fDetElemId/100-1; iCh++)
1138 { fDetElemNumber += fgNDetElemCh[iCh]; }
1140 fTransform->Local2Global(
1141 fDetElemId,lDetElemLocX,lDetElemLocY,lDetElemLocZ,
1142 fDetElemPos[0],fDetElemPos[1],fDetElemPos[2]);
1146 //_____________________________________________________
1147 AliMUONAlignmentTrackRecord* AliMUONAlignment::ProcessTrack( AliMUONTrack* track, Bool_t doAlignment )
1150 // store current track in running member.
1153 // clear track record
1154 fTrackRecord.Clear();
1156 // get number of tracks
1157 Int_t nTrackParam = fTrack->GetTrackParamAtCluster()->GetEntries();
1158 AliDebug(1,Form("Number of track param entries : %i ", nTrackParam));
1160 Bool_t first( kTRUE );
1161 for(Int_t iCluster=0; iCluster<nTrackParam; iCluster++)
1164 // and get new pointers
1165 fTrackParam = (AliMUONTrackParam *) fTrack->GetTrackParamAtCluster()->At(iCluster);
1166 if ( ! fTrackParam ) continue;
1167 fCluster = fTrackParam->GetClusterPtr();
1168 if ( ! fCluster ) continue;
1169 // if (fDetElemId<500) continue;
1171 // fill local variables for this position --> one measurement
1174 FillTrackParamData();
1179 // for first valid cluster, save track position as "starting" values
1182 fTrackPos0[0] = fTrackPos[0];
1183 fTrackPos0[1] = fTrackPos[1];
1184 fTrackPos0[2] = fTrackPos[2];
1185 fTrackSlope0[0] = fTrackSlope[0];
1186 fTrackSlope0[1] = fTrackSlope[1];
1190 // calculate measurements
1191 fCosPhi = TMath::Cos(fPhi);
1192 fSinPhi = TMath::Sin(fPhi);
1196 fMeas[0] = fTrackPos[0] - fClustPos[0];
1197 fMeas[1] = fTrackPos[1] - fClustPos[1];
1201 fMeas[0] = - fClustPos[0];
1202 fMeas[1] = - fClustPos[1];
1207 AliDebug(1,Form("cluster: %i", iCluster));
1208 AliDebug(1,Form("x: %f\t y: %f\t z: %f\t DetElemID: %i\t ", fClustPos[0], fClustPos[1], fClustPos[2], fDetElemId));
1209 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));
1211 AliDebug(1,Form("Track Parameter: %i", iCluster));
1212 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]));
1213 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]));
1215 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]));
1217 // Set local equations
1218 LocalEquationX( doAlignment );
1219 LocalEquationY( doAlignment );
1223 return &fTrackRecord;
1227 //______________________________________________________________________________
1228 void AliMUONAlignment::ProcessTrack( AliMUONAlignmentTrackRecord* track, Bool_t doAlignment )
1231 if( !( track && doAlignment ) ) return;
1233 // loop over clusters
1234 for( Int_t index = 0; index < track->GetNRecords(); ++index )
1235 { if( AliMUONAlignmentClusterRecord* record = track->GetRecord( index ) ) LocalEquation( *record ); }
1241 //_____________________________________________________
1242 void AliMUONAlignment::LocalFit(Int_t iTrack, Double_t *lTrackParam, Int_t lSingleFit)
1245 /// Call local fit for this tracks
1246 Int_t iRes = fMillepede->LocalFit(iTrack,lTrackParam,lSingleFit);
1247 if (iRes && !lSingleFit)
1248 { fMillepede->SetNLocalEquations(fMillepede->GetNLocalEquations()+1); }
1252 //_____________________________________________________
1253 void AliMUONAlignment::GlobalFit(Double_t *parameters,Double_t *errors,Double_t *pulls)
1256 /// Call global fit; Global parameters are stored in parameters
1257 fMillepede->GlobalFit(parameters,errors,pulls);
1259 AliInfo("Done fitting global parameters!");
1260 for (int iGlob=0; iGlob<fgNDetElem; iGlob++)
1261 { 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]); }
1266 //_____________________________________________________
1267 Double_t AliMUONAlignment::GetParError(Int_t iPar)
1269 /// Get error of parameter iPar
1270 Double_t lErr = fMillepede->GetParError(iPar);
1274 //_____________________________________________________
1275 void AliMUONAlignment::PrintGlobalParameters()
1277 /// Print global parameters
1278 fMillepede->PrintGlobalParameters();
1281 //_________________________________________________________________________
1282 TGeoCombiTrans AliMUONAlignment::ReAlign(const TGeoCombiTrans & transform, const double *lMisAlignment) const
1284 /// Realign given transformation by given misalignment and return the misaligned transformation
1286 Double_t cartMisAlig[3] = {0,0,0};
1287 Double_t angMisAlig[3] = {0,0,0};
1288 // const Double_t *trans = transform.GetTranslation();
1289 // TGeoRotation *rot;
1290 // // check if the rotation we obtain is not NULL
1291 // if (transform.GetRotation()) {
1292 // rot = transform.GetRotation();
1295 // rot = new TGeoRotation("rot");
1296 // } // default constructor.
1298 cartMisAlig[0] = -TMath::Sign(1.0,transform.GetRotationMatrix()[0])*lMisAlignment[0];
1299 cartMisAlig[1] = -TMath::Sign(1.0,transform.GetRotationMatrix()[4])*lMisAlignment[1];
1300 cartMisAlig[2] = -TMath::Sign(1.0,transform.GetRotationMatrix()[8])*lMisAlignment[3];
1301 angMisAlig[2] = -TMath::Sign(1.0,transform.GetRotationMatrix()[0]*transform.GetRotationMatrix()[4])*lMisAlignment[2]*180./TMath::Pi();
1303 TGeoTranslation deltaTrans(cartMisAlig[0], cartMisAlig[1], cartMisAlig[2]);
1304 TGeoRotation deltaRot;
1305 deltaRot.RotateX(angMisAlig[0]);
1306 deltaRot.RotateY(angMisAlig[1]);
1307 deltaRot.RotateZ(angMisAlig[2]);
1309 TGeoCombiTrans deltaTransf(deltaTrans,deltaRot);
1310 TGeoHMatrix newTransfMat = transform * deltaTransf;
1312 return TGeoCombiTrans(newTransfMat);
1315 //______________________________________________________________________
1316 AliMUONGeometryTransformer *
1317 AliMUONAlignment::ReAlign(const AliMUONGeometryTransformer * transformer,
1318 const double *misAlignments, Bool_t verbose)
1321 /// Returns a new AliMUONGeometryTransformer with the found misalignments
1324 // Takes the internal geometry module transformers, copies them
1325 // and gets the Detection Elements from them.
1326 // Takes misalignment parameters and applies these
1327 // to the local transform of the Detection Element
1328 // Obtains the global transform by multiplying the module transformer
1329 // transformation with the local transformation
1330 // Applies the global transform to a new detection element
1331 // Adds the new detection element to a new module transformer
1332 // Adds the new module transformer to a new geometry transformer
1333 // Returns the new geometry transformer
1335 Double_t lModuleMisAlignment[4] = {0.,0.,0.,0.};
1336 Double_t lDetElemMisAlignment[4] = {0.,0.,0.,0.};
1337 Int_t iDetElemId = 0;
1338 Int_t iDetElemNumber = 0;
1340 AliMUONGeometryTransformer *newGeometryTransformer =
1341 new AliMUONGeometryTransformer();
1342 for (Int_t iMt = 0; iMt < transformer->GetNofModuleTransformers(); iMt++) {
1343 // module transformers
1344 const AliMUONGeometryModuleTransformer *kModuleTransformer =
1345 transformer->GetModuleTransformer(iMt, true);
1347 AliMUONGeometryModuleTransformer *newModuleTransformer =
1348 new AliMUONGeometryModuleTransformer(iMt);
1349 newGeometryTransformer->AddModuleTransformer(newModuleTransformer);
1351 TGeoCombiTrans moduleTransform =
1352 TGeoCombiTrans(*kModuleTransformer->GetTransformation());
1353 // New module transformation
1354 TGeoCombiTrans newModuleTransform = ReAlign(moduleTransform,lModuleMisAlignment);
1355 newModuleTransformer->SetTransformation(newModuleTransform);
1357 // Get delta transformation:
1358 // Tdelta = Tnew * Told.inverse
1359 TGeoHMatrix deltaModuleTransform =
1360 AliMUONGeometryBuilder::Multiply(newModuleTransform,
1361 kModuleTransformer->GetTransformation()->Inverse());
1362 // Create module mis alignment matrix
1363 newGeometryTransformer
1364 ->AddMisAlignModule(kModuleTransformer->GetModuleId(), deltaModuleTransform);
1366 AliMpExMap *detElements = kModuleTransformer->GetDetElementStore();
1369 AliInfo(Form("%i DEs in old GeometryStore %i",detElements->GetSize(), iMt));
1371 TIter next(detElements->CreateIterator());
1372 AliMUONGeometryDetElement* detElement;
1374 while ( ( detElement = static_cast<AliMUONGeometryDetElement*>(next()) ) )
1377 // make a new detection element
1378 AliMUONGeometryDetElement *newDetElement =
1379 new AliMUONGeometryDetElement(detElement->GetId(),
1380 detElement->GetVolumePath());
1381 TString lDetElemName(detElement->GetDEName());
1382 lDetElemName.ReplaceAll("DE","");
1383 iDetElemId = lDetElemName.Atoi();
1384 iDetElemNumber = iDetElemId%100;
1385 for (int iCh=0; iCh<iDetElemId/100-1; iCh++){
1386 iDetElemNumber += fgNDetElemCh[iCh];
1388 for (int i=0; i<fgNParCh; i++) {
1389 lDetElemMisAlignment[i] = 0.0;
1390 if (iMt<fgNTrkMod) {
1391 AliInfo(Form("iMt %i, iCh %i, iDe %i, iDeId %i, iDeNb %i, iPar %i",iMt, iDetElemId/100, iDe, iDetElemId, iDetElemNumber, iDetElemNumber*fgNParCh+i));
1392 lDetElemMisAlignment[i] = misAlignments[iDetElemNumber*fgNParCh+i];
1395 // local transformation of this detection element.
1396 TGeoCombiTrans localTransform
1397 = TGeoCombiTrans(*detElement->GetLocalTransformation());
1398 TGeoCombiTrans newLocalTransform = ReAlign(localTransform,lDetElemMisAlignment);
1399 newDetElement->SetLocalTransformation(newLocalTransform);
1401 // global transformation
1402 TGeoHMatrix newGlobalTransform =
1403 AliMUONGeometryBuilder::Multiply(newModuleTransform,
1405 newDetElement->SetGlobalTransformation(newGlobalTransform);
1407 // add this det element to module
1408 newModuleTransformer->GetDetElementStore()->Add(newDetElement->GetId(),
1411 // In the Alice Alignment Framework misalignment objects store
1412 // global delta transformation
1413 // Get detection "intermediate" global transformation
1414 TGeoHMatrix newOldGlobalTransform = newModuleTransform * localTransform;
1415 // Get detection element global delta transformation:
1416 // Tdelta = Tnew * Told.inverse
1417 TGeoHMatrix deltaGlobalTransform
1418 = AliMUONGeometryBuilder::Multiply(newGlobalTransform,
1419 newOldGlobalTransform.Inverse());
1421 // Create mis alignment matrix
1422 newGeometryTransformer
1423 ->AddMisAlignDetElement(detElement->GetId(), deltaGlobalTransform);
1427 AliInfo(Form("Added module transformer %i to the transformer", iMt));
1428 newGeometryTransformer->AddModuleTransformer(newModuleTransformer);
1430 return newGeometryTransformer;
1433 //______________________________________________________________________
1434 void AliMUONAlignment::SetAlignmentResolution(const TClonesArray* misAlignArray, Int_t rChId, Double_t rChResX, Double_t rChResY, Double_t rDeResX, Double_t rDeResY)
1437 /// Set alignment resolution to misalign objects to be stored in CDB
1438 Int_t chIdMin = (rChId<0)? 0 : rChId;
1439 Int_t chIdMax = (rChId<0)? 9 : rChId;
1440 Double_t chResX = rChResX;
1441 Double_t chResY = rChResY;
1442 Double_t deResX = rDeResX;
1443 Double_t deResY = rDeResY;
1445 TMatrixDSym mChCorrMatrix(6);
1446 mChCorrMatrix[0][0]=chResX*chResX;
1447 mChCorrMatrix[1][1]=chResY*chResY;
1448 // mChCorrMatrix.Print();
1450 TMatrixDSym mDECorrMatrix(6);
1451 mDECorrMatrix[0][0]=deResX*deResX;
1452 mDECorrMatrix[1][1]=deResY*deResY;
1453 // mDECorrMatrix.Print();
1455 AliAlignObjMatrix *alignMat = 0x0;
1457 for(Int_t chId=chIdMin; chId<=chIdMax; chId++) {
1461 chName1 = Form("GM%d",chId);
1462 chName2 = Form("GM%d",chId);
1464 chName1 = Form("GM%d",4+(chId-4)*2);
1465 chName2 = Form("GM%d",4+(chId-4)*2+1);
1468 for (int i=0; i<misAlignArray->GetEntries(); i++) {
1469 alignMat = (AliAlignObjMatrix*)misAlignArray->At(i);
1470 TString volName(alignMat->GetSymName());
1471 if((volName.Contains(chName1)&&
1472 ((volName.Last('/')==volName.Index(chName1)+chName1.Length())||
1473 (volName.Length()==volName.Index(chName1)+chName1.Length())))||
1474 (volName.Contains(chName2)&&
1475 ((volName.Last('/')==volName.Index(chName2)+chName2.Length())||
1476 (volName.Length()==volName.Index(chName2)+chName2.Length())))){
1477 volName.Remove(0,volName.Last('/')+1);
1478 if (volName.Contains("GM")) {
1479 // alignMat->Print("NULL");
1480 alignMat->SetCorrMatrix(mChCorrMatrix);
1481 } else if (volName.Contains("DE")) {
1482 // alignMat->Print("NULL");
1483 alignMat->SetCorrMatrix(mDECorrMatrix);