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 **************************************************************************/
16 //-----------------------------------------------------------------------------
17 /// \class AliMUONSurveyObj
18 /// Base class for the survey processing of the ALICE DiMuon spectrometer
20 /// This base object provides methods to process the survey+photogrammetry
21 /// data of the Chambers (frames) and Detection Elements of the DiMuon
22 /// Spectrometer and calculate their misalignments.
24 /// \author Javier Castillo
25 //-----------------------------------------------------------------------------
31 #include "TGeoMatrix.h"
37 #include "TGraph2DErrors.h"
41 #include "AliSurveyPoint.h"
43 #include "AliMUONSurveyObj.h"
44 #include "AliMUONSurveyUtil.h"
46 void SurveyFcn(int &npar, double *g, double &f, double *par, int iflag);
49 ClassImp(AliMUONSurveyObj)
52 AliMUONSurveyObj::AliMUONSurveyObj()
60 , fOwnerLocalTrf(kFALSE)
61 , fOwnerAlignTrf(kTRUE)
62 , fOwnerBaseTrf(kFALSE)
73 /// Default constructor
75 fSTargets = new TObjArray();
76 fSTargets->SetOwner(kFALSE);
77 fGBTargets = new TObjArray();
78 fGBTargets->SetOwner(kFALSE);
79 fLBTargets = new TObjArray();
80 fLBTargets->SetOwner(kFALSE);
82 fAlignTrf = new TGeoCombiTrans();
84 fFitter = new TFitter(100);
87 AliMUONSurveyObj::~AliMUONSurveyObj() {
105 if(fOwnerLocalTrf && fLocalTrf) {
109 if(fOwnerAlignTrf && fAlignTrf) {
113 if(fOwnerBaseTrf && fBaseTrf) {
123 void AliMUONSurveyObj::AddStickerTarget(AliSurveyPoint *stPoint){
124 /// Add sticker target
126 fSTargets->Add(ConvertPointUnits(stPoint,0.1));
128 fSTargets->Add(stPoint);
132 void AliMUONSurveyObj::AddGButtonTarget(AliSurveyPoint *btPoint){
133 /// Add global button target
135 fGBTargets->Add(ConvertPointUnits(btPoint,0.1));
137 fGBTargets->Add(btPoint);
141 void AliMUONSurveyObj::AddLButtonTarget(AliSurveyPoint *btPoint){
142 /// Add local button target target; AliSurveyPoint
144 fLBTargets->Add(ConvertPointUnits(btPoint,0.1));
146 fLBTargets->Add(btPoint);
150 void AliMUONSurveyObj::AddLButtonTarget(TVector3 *btVector){
151 /// Add local button target target; TVector3
152 fLBTargets->Add(btVector);
155 Int_t AliMUONSurveyObj::AddStickerTargets(TObjArray *pArray, TString stBaseName, Int_t lTargetMax){
156 /// Add a maximum of lTargetMax sticker targets with stBaseName from pArray
158 AliError(Form("Survey points array is empty %p!",pArray));
161 if (stBaseName.IsNull()) {
162 AliError(Form("Need base name for sticker targets %s!",stBaseName.Data()));
167 AliSurveyPoint *pointSST = 0x0;
171 for (int iPoint=0; iPoint<lTargetMax; iPoint++) {
172 TString stFullName(stBaseName);
173 stNumber = Form("%d",iPoint+1);
174 if(lTargetMax>9&&iPoint+1<10) {
177 stFullName+=stNumber;
179 pointSST = (AliSurveyPoint *)pArray->FindObject(stFullName.Data());
182 AddStickerTarget(pointSST);
183 AliInfo(Form("Added survey sticker target %s at index %d",pointSST->GetName(),stIndex));
188 AliInfo(Form("Found %d sticker targets with base name %s",fSTargets->GetEntries(),stBaseName.Data()));
192 Int_t AliMUONSurveyObj::AddGButtonTargets(TObjArray *pArray, TString btBaseName, Int_t lTargetMax){
193 /// Add a maximum of lTargetMax global button targets with stBaseName from pArray
194 printf("%s \n",btBaseName.Data());
196 AliError(Form("Survey points array is empty %p!",pArray));
199 if (btBaseName.IsNull()) {
200 AliError(Form("Need base name for button targets %s!",btBaseName.Data()));
205 AliSurveyPoint *pointSBT = 0x0;
209 for (int iPoint=0; iPoint<lTargetMax; iPoint++) {
210 TString btFullName(btBaseName);
211 btNumber = Form("%d",iPoint+1);
212 if(lTargetMax>9&&iPoint+1<10) {
215 btFullName+=btNumber;
216 printf("%s \n",btFullName.Data());
217 pointSBT = (AliSurveyPoint *)pArray->FindObject(btFullName.Data());
220 AddGButtonTarget(pointSBT);
221 AliInfo(Form("Added survey button target %s at index %d",pointSBT->GetName(),btIndex));
226 AliInfo(Form("Found %d button targets with base name %s",fGBTargets->GetEntries(),btBaseName.Data()));
230 Int_t AliMUONSurveyObj::AddLButtonTargets(TObjArray *pArray, TString btBaseName, Int_t lTargetMax){
231 /// Add a maximum of lTargetMax local button targets with stBaseName from pArray
232 printf("%s \n",btBaseName.Data());
234 AliError(Form("Local points array is empty %p!",pArray));
237 if (btBaseName.IsNull()) {
238 AliError(Form("Need base name for button targets %s!",btBaseName.Data()));
243 AliSurveyPoint *pointSBT = 0x0;
247 for (int iPoint=0; iPoint<lTargetMax; iPoint++) {
248 TString btFullName(btBaseName);
249 btNumber = Form("%d",iPoint+1);
250 if(lTargetMax>9&&iPoint+1<10) {
253 btFullName+=btNumber;
254 printf("%s \n",btFullName.Data());
255 pointSBT = (AliSurveyPoint *)pArray->FindObject(btFullName.Data());
258 AddLButtonTarget(pointSBT);
259 AliInfo(Form("Added local button target %s at index %d",pointSBT->GetName(),btIndex));
264 AliInfo(Form("Found %d local button targets with base name %s",fLBTargets->GetEntries(),btBaseName.Data()));
268 Int_t AliMUONSurveyObj::GetNStickerTargets() {
269 /// return number of sticker targets
270 return fSTargets->GetEntriesFast();
273 AliSurveyPoint* AliMUONSurveyObj::GetStickerTarget(Int_t stIndex){
274 /// return sticker target at stIndex
275 if (stIndex<0||stIndex>=fSTargets->GetEntriesFast()) {
276 AliError(Form("No sticker target at index %d",stIndex));
280 return (AliSurveyPoint*)fSTargets->At(stIndex);
284 Int_t AliMUONSurveyObj::GetNGButtonTargets() {
285 /// return number of global button targets
286 return fGBTargets->GetEntriesFast();
289 AliSurveyPoint* AliMUONSurveyObj::GetGButtonTarget(Int_t btIndex){
290 /// return global button target at btIndex
291 if (btIndex<0||btIndex>=fGBTargets->GetEntriesFast()) {
292 AliError(Form("No surveyed button target at index %d",btIndex));
296 return (AliSurveyPoint*)fGBTargets->At(btIndex);
300 Int_t AliMUONSurveyObj::GetNLButtonTargets() {
301 /// return number of local button targets
302 return fGBTargets->GetEntriesFast();
305 AliSurveyPoint* AliMUONSurveyObj::GetLButtonTarget(Int_t btIndex){
306 /// return local button target at btIndex
307 if (btIndex<0||btIndex>=fLBTargets->GetEntriesFast()) {
308 AliError(Form("No surveyed button target at index %d",btIndex));
312 if(fLBTargets->At(btIndex)->IsA()==TVector3::Class()){
313 TVector3 *lBT = (TVector3*)fLBTargets->At(btIndex);
315 return (new AliSurveyPoint(TString("local"),(float)lBT->X(),(float)lBT->Y(),(float)lBT->Z(),(float)0.,(float)0.,(float)0.,'B',kTRUE));
316 } else if(fLBTargets->At(btIndex)->IsA()==AliSurveyPoint::Class()) {
317 AliSurveyPoint *lBT = (AliSurveyPoint*)fLBTargets->At(btIndex);
320 AliError(Form("Unexpected class %s ! Valid classes are TVector3 or AliSurveyPoint",fLBTargets->At(btIndex)->ClassName()));
326 void AliMUONSurveyObj::SetPlane(TString pName, Double_t xMin, Double_t xMax, Double_t yMin, Double_t yMax){
327 /// Set the plane function for the plane fitting
332 fPlane = new TF2(pName,this,&AliMUONSurveyObj::EqPlane,xMin,xMax,yMin,yMax,3,"AliMUONSurveyObj","EqPlane");
335 void AliMUONSurveyObj::SetPlaneParameters(Double_t p0, Double_t p1, Double_t p2) {
336 /// Set the parameters of plane function for the plane fitting
338 AliError("Must use SetPlane before SetPlaneParameters!!!");
341 fPlane->SetParameter(0,p0);
342 fPlane->SetParameter(1,p1);
343 fPlane->SetParameter(2,p2);
347 void AliMUONSurveyObj::DrawSTargets() {
348 /// Draw a graph of the sticker targets
349 TGraph2DErrors *gST = new TGraph2DErrors(3);
350 AliSurveyPoint *pST = 0x0;
351 for (Int_t iPoint=0; iPoint<GetNStickerTargets(); iPoint++) {
352 pST = GetStickerTarget(iPoint);
353 // pST->PrintPoint();
354 gST->SetPoint(iPoint,pST->GetX(),pST->GetY(),pST->GetZ());
355 gST->SetPointError(iPoint,pST->GetPrecisionX(),pST->GetPrecisionY(),pST->GetPrecisionZ());
357 gST->DrawClone("P0");
362 Double_t AliMUONSurveyObj::FitPlane() {
363 /// Fit plane to sticker targets
365 AliError("Must use SetPlane before FitPlane!!!");
368 if (fSTargets->GetEntriesFast()<3) {
369 AliError("Not enough sticker targets (%d) for plane fitting!!!");
373 Double_t pl[3] = {0};
374 Double_t pg[3] = {0};
376 TGraph2DErrors *gST = new TGraph2DErrors(3);
377 AliSurveyPoint *pST = 0x0;
378 for (Int_t iPoint=0; iPoint<GetNStickerTargets(); iPoint++) {
379 pST = GetStickerTarget(iPoint);
380 // pST->PrintPoint();
384 fBaseTrf->MasterToLocal(pg,pl);
385 gST->SetPoint(iPoint,pl[0],pl[1],pl[2]);
386 printf("%d %f %f %f\n",iPoint,pl[0],pl[1],pl[2]);
387 gST->SetPointError(iPoint,pST->GetPrecisionX(),pST->GetPrecisionY(),pST->GetPrecisionZ());
393 return fPlane->GetChisquare();
396 Double_t AliMUONSurveyObj::SurveyChi2(Double_t *par){
397 /// Returns the chisquare between local2global transform of local button targets and their surveyed position
398 TGeoTranslation transTemp;
399 TGeoRotation rotTemp;
400 TGeoCombiTrans trfTemp;
404 trfTemp.SetTranslation(transTemp);
405 trfTemp.SetRotation(rotTemp);
407 trfTemp.RotateZ(TMath::RadToDeg()*par[5]);
408 trfTemp.RotateY(TMath::RadToDeg()*par[4]);
409 trfTemp.RotateX(TMath::RadToDeg()*par[3]);
410 trfTemp.SetTranslation(par[0],par[1],par[2]);
412 TGeoHMatrix matGlo = (*fBaseTrf)*trfTemp;
413 TGeoCombiTrans trfGlo(matGlo);
415 Double_t pl[3] = {0};
416 Double_t pg[3] = {0};
418 for(Int_t iPoint=0; iPoint<fGBTargets->GetEntries(); iPoint++){
419 AliSurveyPoint *gBT = (AliSurveyPoint*)fGBTargets->At(iPoint);
420 if(fLBTargets->At(iPoint)->IsA()==TVector3::Class()){
421 TVector3 *lBT = (TVector3*)fLBTargets->At(iPoint);
425 } else if(fLBTargets->At(iPoint)->IsA()==AliSurveyPoint::Class()) {
426 AliSurveyPoint *lBT = (AliSurveyPoint*)fLBTargets->At(iPoint);
431 AliError(Form("Unexpected class %s ! Valid classes are TVector3 or AliSurveyPoint",fLBTargets->At(iPoint)->ClassName()));
435 trfGlo.LocalToMaster(pl, pg);
436 // printf("%d %f %f %f\n",iPoint,pg[0],pg[1],pg[2]);
437 if(fLBTargets->At(iPoint)->IsA()==TVector3::Class()){
438 lChi2 += (pg[0]-gBT->GetX())*(pg[0]-gBT->GetX())/(gBT->GetPrecisionX()*gBT->GetPrecisionX());
439 lChi2 += (pg[1]-gBT->GetY())*(pg[1]-gBT->GetY())/(gBT->GetPrecisionY()*gBT->GetPrecisionY());
440 lChi2 += (pg[2]-gBT->GetZ())*(pg[2]-gBT->GetZ())/(gBT->GetPrecisionZ()*gBT->GetPrecisionZ());
441 } else if(fLBTargets->At(iPoint)->IsA()==AliSurveyPoint::Class()) {
442 AliSurveyPoint *lBT = (AliSurveyPoint*)fLBTargets->At(iPoint);
443 lChi2 += (pg[0]-gBT->GetX())*(pg[0]-gBT->GetX())/(gBT->GetPrecisionX()*gBT->GetPrecisionX()+lBT->GetPrecisionX()*lBT->GetPrecisionX());
444 lChi2 += (pg[1]-gBT->GetY())*(pg[1]-gBT->GetY())/(gBT->GetPrecisionY()*gBT->GetPrecisionY()+lBT->GetPrecisionY()*lBT->GetPrecisionY());
445 lChi2 += (pg[2]-gBT->GetZ())*(pg[2]-gBT->GetZ())/(gBT->GetPrecisionZ()*gBT->GetPrecisionZ()+lBT->GetPrecisionZ()*lBT->GetPrecisionZ());
447 AliError(Form("Unexpected class %s ! Valid classes are TVector3 or AliSurveyPoint",fLBTargets->At(iPoint)->ClassName()));
455 //_____________________________________________________________________________
456 void SurveyFcn(int &npar, double *g, double &f, double *par, int iflag) {
458 /// Standard function as needed by Minuit-like minimization procedures.
459 /// For the set of parameters par calculates and returns chi-squared.
462 // smuggle a C++ object into a C function
463 AliMUONSurveyObj *aSurveyObj = (AliMUONSurveyObj*) gMinuit->GetObjectFit();
465 f = aSurveyObj->SurveyChi2(par);
468 if (g) {} // no warnings about unused stuff...
472 //_____________________________________________________________________________
473 Int_t AliMUONSurveyObj::SurveyToAlign(TGeoCombiTrans &quadTransf, Double_t *parErr, Double_t psi, Double_t tht, Double_t epsi, Double_t etht) {
474 /// Main function to obtain the misalignments from the surveyed position of the button targets;
475 if (fGBTargets->GetEntries()!=fLBTargets->GetEntries()){
476 AliError(Form("Different number of button targets: %d survey points and %d local coord!",
477 fGBTargets->GetEntries(),fLBTargets->GetEntries()));
482 gMinuit->SetObjectFit(this);
483 fitter.SetFCN(SurveyFcn);
484 fitter.SetParameter(0,"dx",0,0.1,-20,20);
485 fitter.SetParameter(1,"dy",0,0.1,-20,20);
486 fitter.SetParameter(2,"dz",0,0.1,0,0);
487 fitter.SetParameter(3,"rx",psi,0.0001,-0.05,0.05);
488 fitter.SetParameter(4,"ry",tht,0.0001,-0.05,0.05);
489 // fitter.SetParameter(3,"rx",psi,0.0001,psi-5*epsi,psi+5*epsi);
490 // fitter.SetParameter(4,"ry",tht,0.0001,tht-5*etht,tht+5*etht);
491 fitter.SetParameter(5,"rz",0,0.0001,0,0);
493 if(psi) fitter.FixParameter(3);
494 if(tht) fitter.FixParameter(4);
498 fitter.ExecuteCommand("SET PRINT", arglist, 1);
499 fitter.ExecuteCommand("SET ERR", arglist, 1);
501 //fitter.ExecuteCommand("SIMPLEX", arglist, 1);
502 // fitter.ExecuteCommand("MINIMIZE", arglist, 1);
503 fitter.ExecuteCommand("MIGRAD", arglist, 1);
504 fitter.ExecuteCommand("IMPROVE", arglist, 1);
505 // fitter.ExecuteCommand("MINOS", arglist, 1);
506 // fitter.ExecuteCommand("CALL 3", arglist,0);
508 for (int j=0; j<3; j++) printf("%10.3f ",fitter.GetParameter(j));
509 for (int j=3; j<6; j++) printf("%10.3f ",1000*fitter.GetParameter(j));
511 for (int j=0; j<3; j++) printf("%10.3f ",fitter.GetParError(j));
512 for (int j=3; j<6; j++) printf("%10.3f ",1000*fitter.GetParError(j));
516 quadTransf.RotateZ(TMath::RadToDeg()*fitter.GetParameter(5));
517 quadTransf.RotateY(TMath::RadToDeg()*fitter.GetParameter(4));
518 quadTransf.RotateX(TMath::RadToDeg()*fitter.GetParameter(3));
519 quadTransf.SetTranslation(fitter.GetParameter(0),fitter.GetParameter(1),fitter.GetParameter(2));
521 for(Int_t iPar=0; iPar<6; iPar++){
522 parErr[iPar] = fitter.GetParError(iPar);
524 if(epsi) parErr[3] = epsi;
525 if(etht) parErr[4] = etht;
531 //_____________________________________________________________________________
532 Int_t AliMUONSurveyObj::SurveyToAlign(Double_t psi, Double_t tht, Double_t epsi, Double_t etht) {
533 /// Main function to obtain the misalignments from the surveyed position of the button targets;
534 if (fGBTargets->GetEntries()!=fLBTargets->GetEntries()){
535 AliError(Form("Different number of button targets: %d survey points and %d local coord!",
536 fGBTargets->GetEntries(),fLBTargets->GetEntries()));
540 // TFitter fitter(100);
541 gMinuit->SetObjectFit(this);
542 fFitter->SetFCN(SurveyFcn);
543 fFitter->SetParameter(0,"dx",0,0.1,-20,20);
544 fFitter->SetParameter(1,"dy",0,0.1,-20,20);
545 fFitter->SetParameter(2,"dz",0,0.1,0,0);
547 fFitter->SetParameter(3,"rx",psi,epsi,psi-5*epsi,psi+5*epsi);
549 fFitter->SetParameter(3,"rx",psi,0.0001,-0.05,0.05);
551 fFitter->SetParameter(4,"ry",tht,etht,tht-5*etht,tht+5*etht);
553 fFitter->SetParameter(4,"ry",tht,0.0001,-0.05,0.05);
554 // fFitter->SetParameter(3,"rx",psi,0.0001,psi-5*epsi,psi+5*epsi);
555 // fFitter->SetParameter(4,"ry",tht,0.0001,tht-5*etht,tht+5*etht);
556 fFitter->SetParameter(5,"rz",0,0.0001,0,0);
558 if(psi) fFitter->FixParameter(3);
559 if(tht) fFitter->FixParameter(4);
563 fFitter->ExecuteCommand("SET PRINT", arglist, 1);
564 fFitter->ExecuteCommand("SET ERR", arglist, 1);
566 //fFitter->ExecuteCommand("SIMPLEX", arglist, 1);
567 // fFitter->ExecuteCommand("MINIMIZE", arglist, 1);
568 fFitter->ExecuteCommand("MIGRAD", arglist, 1);
569 fFitter->ExecuteCommand("IMPROVE", arglist, 1);
570 // fFitter->ExecuteCommand("MINOS", arglist, 1);
571 // fFitter->ExecuteCommand("CALL 3", arglist,0);
573 for (int j=0; j<3; j++) printf("%10.3f ",fFitter->GetParameter(j));
574 for (int j=3; j<6; j++) printf("%10.3f ",1000*fFitter->GetParameter(j));
576 for (int j=0; j<3; j++) printf("%10.3f ",fFitter->GetParError(j));
577 for (int j=3; j<6; j++) printf("%10.3f ",1000*fFitter->GetParError(j));
581 fAlignTrf->RotateZ(TMath::RadToDeg()*fFitter->GetParameter(5));
582 fAlignTrf->RotateY(TMath::RadToDeg()*fFitter->GetParameter(4));
583 fAlignTrf->RotateX(TMath::RadToDeg()*fFitter->GetParameter(3));
584 fAlignTrf->SetTranslation(fFitter->GetParameter(0),fFitter->GetParameter(1),fFitter->GetParameter(2));
586 if(epsi) fFitter->ReleaseParameter(3); // To get error
587 if(etht) fFitter->ReleaseParameter(4); // To get error
589 TGeoCombiTrans lGlobalTrf = TGeoCombiTrans((*fBaseTrf)*(*fAlignTrf));
590 AliSurveyPoint *pointGBT;
591 AliSurveyPoint *pointLBT;
592 Double_t pl[3] = {0};
593 Double_t pg[3] = {0};
594 for (int iPoint=0; iPoint<GetNGButtonTargets(); iPoint++){
595 pointGBT=GetGButtonTarget(iPoint);
596 pointLBT=GetLButtonTarget(iPoint);
597 pl[0] = pointLBT->GetX();
598 pl[1] = pointLBT->GetY();
599 pl[2] = pointLBT->GetZ();
600 lGlobalTrf.LocalToMaster(pl,pg);
601 printf("Point %d local: %.3f %.3f %.3f\n",iPoint,pl[0],pl[1],pl[2]);
602 printf("Point %d global: %.3f %.3f %.3f\n",iPoint,pg[0],pg[1],pg[2]);
603 printf("Point %d survey: %.3f %.3f %.3f\n",iPoint,pointGBT->GetX(),pointGBT->GetY(),pointGBT->GetZ());
610 Double_t AliMUONSurveyObj::EvalFunction(const TF2 *lFunction, Int_t iP1, Int_t iP2, const Char_t *lCoord) {
611 /// Evaluate the given function at the given points for the given coordinate
613 AliError("No function given!!!");
616 AliSurveyPoint *gP1 = GetGButtonTarget(iP1);
617 AliSurveyPoint *gP2 = GetGButtonTarget(iP2);
620 AliError("Missing global button target!!!");
624 // AliInfo(Form("Function %s parameters %f %f %f %f %f %f",lFunction->GetName(),lFunction->GetParameter(0),lFunction->GetParameter(1),lFunction->GetParameter(2),lFunction->GetParameter(3),lFunction->GetParameter(4),lFunction->GetParameter(5)));
625 Double_t pl1[3] = {0};
626 Double_t pl2[3] = {0};
627 Double_t pg1[3] = {0};
628 Double_t pg2[3] = {0};
630 pg1[0] = gP1->GetX();
631 pg1[1] = gP1->GetY();
632 pg1[2] = gP1->GetZ();
633 pg2[0] = gP2->GetX();
634 pg2[1] = gP2->GetY();
635 pg2[2] = gP2->GetZ();
637 fBaseTrf->MasterToLocal(pg1,pl1);
638 fBaseTrf->MasterToLocal(pg2,pl2);
644 lVal = lFunction->Eval(pl1[0],pl2[0]);
645 // lVal = lFunction->Eval(gP1->GetX(),gP2->GetX());
646 // AliInfo(Form("case X, lVal = %f",lVal));
651 lVal = lFunction->Eval(pl1[1],pl2[1]);
652 // lVal = lFunction->Eval(gP1->GetY(),gP2->GetY());
653 // AliInfo(Form("case Y, lVal = %f",lVal));
658 lVal = lFunction->Eval(pl1[2],pl2[2]);
659 // lVal = lFunction->Eval(gP1->GetZ(),gP2->GetZ());
660 // AliInfo(Form("case Z, lVal = %f",lVal));
665 AliError(Form("Coordinate %s is not valid, options are X Y Z",lCoord));
671 void AliMUONSurveyObj::CalculateTranslation(TF2 *xFunc, TF2 *yFunc, TF2 *zFunc, Int_t iP1, Int_t iP2, Double_t *lCenTemp) {
672 /// Calculate the center translation using analytic functions
673 lCenTemp[0] = EvalFunction(xFunc,iP1,iP2,"X");
674 lCenTemp[1] = EvalFunction(yFunc,iP1,iP2,"Y");
675 lCenTemp[2] = EvalFunction(zFunc,iP1,iP2,"Z");
679 Double_t AliMUONSurveyObj::CalculateGlobalDiff(TGeoCombiTrans &lTransf, Int_t nPoints, TArrayD &lDiff){
680 /// By hand computation of distance between local2global transform of target position and its surveyed position
681 if (nPoints > GetNGButtonTargets()) {
682 nPoints = GetNGButtonTargets();
685 for(Int_t iVal=0; iVal<nPoints*(3+1)+1; iVal++){
689 Double_t pl[3] = {0};
690 Double_t pg[3] = {0};
691 Double_t pml[3] = {0};
692 Double_t pmg[3] = {0};
693 AliSurveyPoint *gBT = 0x0;
694 for(Int_t iPoint=0; iPoint<nPoints; iPoint++){
695 gBT = GetGButtonTarget(iPoint);
696 if(!gBT||!fLBTargets->At(iPoint)){
697 AliError(Form("The local or global target %d is missing!",iPoint));
698 lDiff[nPoints*(3+1)] = 1.e7;
699 return lDiff[nPoints*(3+1)];
701 if(fLBTargets->At(iPoint)->IsA()==TVector3::Class()){
702 TVector3 *lBT = (TVector3*)fLBTargets->At(iPoint);
706 } else if(fLBTargets->At(iPoint)->IsA()==AliSurveyPoint::Class()) {
707 AliSurveyPoint *lBT = (AliSurveyPoint*)fLBTargets->At(iPoint);
712 AliError(Form("Unexpected class %s ! Valid classes are TVector3 or AliSurveyPoint",fLBTargets->At(iPoint)->ClassName()));
716 lTransf.LocalToMaster(pl,pg);
717 pmg[0] = gBT->GetX();
718 pmg[1] = gBT->GetY();
719 pmg[2] = gBT->GetZ();
720 fBaseTrf->MasterToLocal(pmg,pml);
721 // printf("l %d %f %f %f\n",iPoint,pl[0],pl[1],pl[2]);
722 // printf("g %d %f %f %f\n",iPoint,pg[0],pg[1],pg[2]);
723 // printf("ml %d %f %f %f\n",iPoint,pml[0],pml[1],pml[2]);
724 // printf("mg %d %f %f %f\n",iPoint,gBT->GetX(),gBT->GetY(),gBT->GetZ());
725 lDiff[iPoint*(3+1)+0] = (pml[0]-pg[0]);
726 lDiff[iPoint*(3+1)+1] = (pml[1]-pg[1]);
727 lDiff[iPoint*(3+1)+2] = (pml[2]-pg[2]);
729 lDiff[iPoint*(3+1)+3] = TMath::Sqrt(lDiff[iPoint*(3+1)+0]*lDiff[iPoint*(3+1)+0]+
730 lDiff[iPoint*(3+1)+1]*lDiff[iPoint*(3+1)+1]+
731 lDiff[iPoint*(3+1)+2]*lDiff[iPoint*(3+1)+2]);
733 lDiff[nPoints*(3+1)] += lDiff[iPoint*(3+1)+3]*lDiff[iPoint*(3+1)+3];
736 lDiff[nPoints*(3+1)] = TMath::Sqrt(lDiff[nPoints*(3+1)]);
737 return lDiff[nPoints*(3+1)];
740 Int_t AliMUONSurveyObj::CalculateBestTransf(Int_t iP1, Int_t iP2, Double_t *lXYZ, Double_t *lPTP) {
741 /// By hand calculation of the best local to global transform using 2 button targets
742 Double_t lPsi = lPTP[0];
743 Double_t lTht = lPTP[1];
745 Double_t pl1[3] = {0};
746 Double_t pl2[3] = {0};
748 if(!fLBTargets->At(iP1)||!fLBTargets->At(iP2)){
749 AliError(Form("Local target %d or %d is missing!",iP1,iP2));
753 if(fLBTargets->At(iP1)->IsA()==TVector3::Class()){
754 TVector3 *lBT1 = (TVector3*)fLBTargets->At(iP1);
758 } else if(fLBTargets->At(iP1)->IsA()==AliSurveyPoint::Class()) {
759 AliSurveyPoint *lBT1 = (AliSurveyPoint*)fLBTargets->At(iP1);
764 AliError(Form("Unexpected class %s ! Valid classes are TVector3 or AliSurveyPoint",fLBTargets->At(iP1)->ClassName()));
767 if(fLBTargets->At(iP2)->IsA()==TVector3::Class()){
768 TVector3 *lBT2 = (TVector3*)fLBTargets->At(iP2);
772 } else if(fLBTargets->At(iP2)->IsA()==AliSurveyPoint::Class()) {
773 AliSurveyPoint *lBT2 = (AliSurveyPoint*)fLBTargets->At(iP2);
778 AliError(Form("Unexpected class %s ! Valid classes are TVector3 or AliSurveyPoint",fLBTargets->At(iP2)->ClassName()));
783 AliMUONSurveyUtil *surveyUtil = AliMUONSurveyUtil::Instance();
786 const char *fxcName = "fXcn00";
787 TF2 **fXc = new TF2*[2];
789 fXc[0] = new TF2(fxcName,surveyUtil,&AliMUONSurveyUtil::XnCenter,fXMin,fXMax,fYMin,fYMax,7,"AliMUONSurveyUtil","XnCenter");
791 fXc[1] = new TF2(fxcName,surveyUtil,&AliMUONSurveyUtil::XpCenter,fXMin,fXMax,fYMin,fYMax,7,"AliMUONSurveyUtil","XpCenter");
794 const char *fycName = "fYcn00";
795 TF2 **fYc = new TF2*[2];
797 fYc[0] = new TF2(fycName,surveyUtil,&AliMUONSurveyUtil::YnCenter,fYMin,fYMax,fYMin,fYMax,8,"AliMUONSurveyUtil","YnCenter");
799 fYc[1] = new TF2(fycName,surveyUtil,&AliMUONSurveyUtil::YpCenter,fYMin,fYMax,fYMin,fYMax,8,"AliMUONSurveyUtil","YpCenter");
802 const char *fzcName = "fZcn00";
803 TF2 **fZc = new TF2*[2];
805 fZc[0] = new TF2(fzcName,surveyUtil,&AliMUONSurveyUtil::ZnCenter,fZMin,fZMax,fZMin,fZMax,8,"AliMUONSurveyUtil","ZnCenter");
807 fZc[1] = new TF2(fzcName,surveyUtil,&AliMUONSurveyUtil::ZpCenter,fZMin,fZMax,fZMin,fZMax,8,"AliMUONSurveyUtil","ZpCenter");
809 // Phi rotation using xglobal coords functions
810 const char *fphixName = "fPhiXnn00";
811 TF2 ***fPhiX = new TF2**[2];
812 for (Int_t iX =0; iX<2; iX++) {
813 fPhiX[iX] = new TF2*[2];
815 fphixName = "fPhiXnn";
816 fPhiX[0][0] = new TF2(fphixName,surveyUtil,&AliMUONSurveyUtil::PhiXnn,fXMin,fXMax,fXMin,fXMax,7,"AliMUONSurveyUtil","PhiXnn");
817 fphixName = "fPhiXnp";
818 fPhiX[0][1] = new TF2(fphixName,surveyUtil,&AliMUONSurveyUtil::PhiXnp,fXMin,fXMax,fXMin,fXMax,7,"AliMUONSurveyUtil","PhiXnp");
819 fphixName = "fPhiXpn";
820 fPhiX[1][0] = new TF2(fphixName,surveyUtil,&AliMUONSurveyUtil::PhiXpn,fXMin,fXMax,fXMin,fXMax,7,"AliMUONSurveyUtil","PhiXpn");
821 fphixName = "fPhiXpp";
822 fPhiX[1][1] = new TF2(fphixName,surveyUtil,&AliMUONSurveyUtil::PhiXpp,fXMin,fXMax,fXMin,fXMax,7,"AliMUONSurveyUtil","PhiXpp");
824 // Phi rotation using yglobal coords functions
825 const char *fphiyName = "fPhiYnn00";
826 TF2 ***fPhiY = new TF2**[2];
827 for (Int_t iY =0; iY<2; iY++) {
828 fPhiY[iY] = new TF2*[2];
830 fphiyName = "fPhiYnn";
831 fPhiY[0][0] = new TF2(fphiyName,surveyUtil,&AliMUONSurveyUtil::PhiYnn,fYMin,fYMax,fYMin,fYMax,8,"AliMUONSurveyUtil","PhiYnn");
832 fphiyName = "fPhiYnp";
833 fPhiY[0][1] = new TF2(fphiyName,surveyUtil,&AliMUONSurveyUtil::PhiYnp,fYMin,fYMax,fYMin,fYMax,8,"AliMUONSurveyUtil","PhiYnp");
834 fphiyName = "fPhiYpn";
835 fPhiY[1][0] = new TF2(fphiyName,surveyUtil,&AliMUONSurveyUtil::PhiYpn,fYMin,fYMax,fYMin,fYMax,8,"AliMUONSurveyUtil","PhiYpn");
836 fphiyName = "fPhiYpp";
837 fPhiY[1][1] = new TF2(fphiyName,surveyUtil,&AliMUONSurveyUtil::PhiYpp,fYMin,fYMax,fYMin,fYMax,8,"AliMUONSurveyUtil","PhiYpp");
840 // Set Parameters of functions
841 for(Int_t iS=0; iS<2; iS++){
842 fXc[iS]->SetParameters(pl1[0],pl1[1],pl1[2],pl2[0],pl2[1],pl2[2],lTht);
843 fYc[iS]->SetParameters(pl1[0],pl1[1],pl1[2],pl2[0],pl2[1],pl2[2],lPsi,lTht);
844 fZc[iS]->SetParameters(pl1[0],pl1[1],pl1[2],pl2[0],pl2[1],pl2[2],lPsi,lTht);
845 // fXc[iS]->SetParameters(lBT1->X(),lBT1->Y(),lBT1->Z(),lBT2->X(),lBT2->Y(),lBT2->Z(),lTht);
846 // fYc[iS]->SetParameters(lBT1->X(),lBT1->Y(),lBT1->Z(),lBT2->X(),lBT2->Y(),lBT2->Z(),lPsi,lTht);
847 // fZc[iS]->SetParameters(lBT1->X(),lBT1->Y(),lBT1->Z(),lBT2->X(),lBT2->Y(),lBT2->Z(),lPsi,lTht);
848 for(Int_t jS=0; jS<2; jS++){
849 fPhiX[iS][jS]->SetParameters(pl1[0],pl1[1],pl1[2],pl2[0],pl2[1],pl2[2],lTht);
850 fPhiY[iS][jS]->SetParameters(pl1[0],pl1[1],pl1[2],pl2[0],pl2[1],pl2[2],lPsi,lTht);
851 // fPhiX[iS][jS]->SetParameters(lBT1->X(),lBT1->Y(),lBT1->Z(),lBT2->X(),lBT2->Y(),lBT2->Z(),lTht);
852 // fPhiY[iS][jS]->SetParameters(lBT1->X(),lBT1->Y(),lBT1->Z(),lBT2->X(),lBT2->Y(),lBT2->Z(),lPsi,lTht);
856 Double_t lCenTemp[3];
857 Double_t lRotTemp[3];
859 TGeoCombiTrans trfTemp;
861 Int_t nPoints = GetNGButtonTargets();
863 TArrayD lDiffTemp(nPoints*(3+1)+1);
864 TArrayD lDiffMin(nPoints*(3+1)+1);
866 for(Int_t i=0; i<nPoints*(3+1)+1; i++){
867 lDiffMin[i]=1000000.;
872 // Calculate Detection Element Center from button targets
875 // Trying 2x*2y*2z*(2phi+2phi) possibilities
876 for(Int_t iX=0; iX<2; iX++){
877 for(Int_t iY=0; iY<2; iY++){
878 for(Int_t iZ=0; iZ<2; iZ++){
879 CalculateTranslation(fXc[iX],fYc[iY],fZc[iZ],iP1,iP2,lCenTemp);
883 for(Int_t iP=0; iP<2; iP++){
884 lRotTemp[2] = EvalFunction(fPhiX[iX][iP],iP1,iP2,"X");
887 trfTemp.RotateZ(TMath::RadToDeg()*lRotTemp[2]);
888 trfTemp.RotateY(TMath::RadToDeg()*lRotTemp[1]);
889 trfTemp.RotateX(TMath::RadToDeg()*lRotTemp[0]);
890 trfTemp.SetTranslation(lCenTemp);
892 if(CalculateGlobalDiff(trfTemp,nPoints,lDiffTemp)<lDiffMin[nPoints*(3+1)]){
894 for(Int_t i=0; i<nPoints*(3+1)+1; i++){
895 printf(" %f",lDiffTemp[i]);
898 printf(" : mycenX%dY%dZ%d(%f,%f,%f); rotx%d(%f,%f,%f)\n",iX,iY,iZ,lCenTemp[0],lCenTemp[1],lCenTemp[2],iP,lRotTemp[0],lRotTemp[1],lRotTemp[2]);
899 printf("Transformation improved ...\n");
900 for (int i=0; i<3; i++) {
901 lXYZ[i] = lCenTemp[i];
903 lPTP[2] = lRotTemp[2];
904 for(Int_t i=0; i<nPoints*(3+1)+1; i++){
905 lDiffMin[i]=lDiffTemp[i];
909 for(Int_t iP=0; iP<2; iP++){
910 lRotTemp[2] = EvalFunction(fPhiY[iY][iP],iP1,iP2,"Y");
913 trfTemp.RotateZ(TMath::RadToDeg()*lRotTemp[2]);
914 trfTemp.RotateY(TMath::RadToDeg()*lRotTemp[1]);
915 trfTemp.RotateX(TMath::RadToDeg()*lRotTemp[0]);
916 trfTemp.SetTranslation(lCenTemp);
918 if(CalculateGlobalDiff(trfTemp,nPoints,lDiffTemp)<lDiffMin[nPoints*(3+1)]){
920 for(Int_t i=0; i<nPoints*(3+1)+1; i++){
921 printf(" %f",lDiffTemp[i]);
924 printf(" : mycenX%dY%dZ%d(%f,%f,%f); roty%d(%f,%f,%f)\n",iX,iY,iZ,lCenTemp[0],lCenTemp[1],lCenTemp[2],iP,lRotTemp[0],lRotTemp[1],lRotTemp[2]);
925 printf("Transformation improved ...\n");
926 for (int i=0; i<3; i++) {
927 lXYZ[i] = lCenTemp[i];
929 lPTP[2] = lRotTemp[2];
930 for(Int_t i=0; i<nPoints*(3+1)+1; i++){
931 lDiffMin[i]=lDiffTemp[i];
939 for (Int_t i=0; i<2; i++) {
943 for (Int_t j=0; j<2; j++) {
956 if (lDiffMin[nPoints*(3+1)]>20) return 0;
961 void AliMUONSurveyObj::CalculateMeanTransf(Double_t *lXYZ, Double_t *lPTP) {
962 /// By hand calculation of the mean (for nPairs of targets) of the best local to global transform using 2 button targets
969 Int_t nPoints = GetNGButtonTargets();
970 // Loop over all possible pairs of button tragets
971 for(Int_t iP1=0; iP1<nPoints; iP1++){
972 for(Int_t iP2=iP1+1; iP2<nPoints; iP2++){
973 printf("%d and %d\n",iP1,iP2);
975 if(CalculateBestTransf(iP1,iP2,lXYZ,lPTP)) {
994 void AliMUONSurveyObj::PrintLocalTrf() {
995 /// Print the local transformation
996 Double_t lRotTemp[3];
997 AliMUONSurveyUtil::MatrixToAngles(fLocalTrf->GetRotationMatrix(),lRotTemp);
998 printf("(%.3f %.3f %.3f), (%.6f %.6f %.6f)\n",fLocalTrf->GetTranslation()[0],fLocalTrf->GetTranslation()[1],fLocalTrf->GetTranslation()[2],lRotTemp[0],lRotTemp[1],lRotTemp[2]);
1001 void AliMUONSurveyObj::PrintAlignTrf() {
1002 /// Print the alignment transformation
1003 Double_t lRotTemp[3];
1004 AliMUONSurveyUtil::MatrixToAngles(fAlignTrf->GetRotationMatrix(),lRotTemp);
1005 printf("(%.3f %.3f %.3f), (%.6f %.6f %.6f)\n",fAlignTrf->GetTranslation()[0],fAlignTrf->GetTranslation()[1],fAlignTrf->GetTranslation()[2],lRotTemp[0],lRotTemp[1],lRotTemp[2]);
1008 void AliMUONSurveyObj::FillSTHistograms(TString baseNameC, TH2 *hSTc, TString baseNameA, TH2 *hSTa) {
1009 /// Fill sticker target histograms for monitoring
1010 if(baseNameC.IsNull()||!hSTc){
1011 AliError("Need base name for points on side C and/or a histogram for them!");
1014 AliSurveyPoint *pointST = 0x0;
1015 for (Int_t iPoint=0; iPoint<GetNStickerTargets(); iPoint++) {
1016 pointST = GetStickerTarget(iPoint);
1017 if (!pointST) continue;
1018 if (pointST->GetPointName().Contains(baseNameC)){
1019 hSTc->Fill(pointST->GetX(),pointST->GetY(),-pointST->GetZ());
1020 } else if ((!baseNameA.IsNull()) &&
1021 (pointST->GetPointName().Contains(baseNameA))) {
1023 AliError("Base name for points on side A provided but no histogram for them!");
1026 hSTa->Fill(pointST->GetX(),pointST->GetY(),-pointST->GetZ());
1031 Double_t AliMUONSurveyObj::GetAlignResX() {
1032 /// Returns the uncertainty of the x translation parameter
1034 AliError("There is no fitter for this object! X resolution will be 0.");
1037 return fFitter->GetParError(0);
1040 Double_t AliMUONSurveyObj::GetAlignResY() {
1041 /// Returns the uncertainty of the y translation parameter
1043 AliError("There is no fitter for this object! Y resolution will be 0.");
1046 return fFitter->GetParError(1);
1049 AliSurveyPoint* AliMUONSurveyObj::ConvertPointUnits(AliSurveyPoint *stPoint, Float_t lFactor) {
1050 /// Return the AliSurveyPoint with new units. Default is from mm -> cm
1051 return new AliSurveyPoint(stPoint->GetPointName(),
1052 lFactor*stPoint->GetX(),lFactor*stPoint->GetY(),lFactor*stPoint->GetZ(),
1053 lFactor*stPoint->GetPrecisionX(),lFactor*stPoint->GetPrecisionY(),lFactor*stPoint->GetPrecisionZ(),
1054 stPoint->GetType(), stPoint->GetTarget());