#include <Rtypes.h>
#include "TGeoMatrix.h"
#include "TMath.h"
#include "TFile.h"
#include "TRandom.h"
#include "TGeoManager.h"
#include "TGeoVolume.h"
#include "TGeoBBox.h"
#include "TGeoTrd1.h"
#include "TGeoPhysicalNode.h"
#include "TGeoNode.h"
#include "TObjString.h"
#include "AliLog.h"
#include "AliAlignObjParams.h"
#include "AliAlignObjMatrix.h"
#include "AliCDBManager.h"
#include "AliCDBMetaData.h"
#include "AliCDBId.h"
#include "AliCDBEntry.h"
#include "AliTOFAlignment.h"
#include "AliSurveyObj.h"
#include "AliSurveyPoint.h"
#include <cstdlib>
ClassImp(AliTOFAlignment)
const Double_t AliTOFAlignment::fgkRorigTOF = 384.5;
const Double_t AliTOFAlignment::fgkX1BTOF = 124.5;
const Double_t AliTOFAlignment::fgkX2BTOF = 134.7262;
const Double_t AliTOFAlignment::fgkYBTOF = 747.2;
const Double_t AliTOFAlignment::fgkZBTOF = 29.0;
const Double_t AliTOFAlignment::fgkXFM = 38.0;
const Double_t AliTOFAlignment::fgkYFM = 457.3;
const Double_t AliTOFAlignment::fgkZFM = 11.2;
AliTOFAlignment::AliTOFAlignment():
TTask("AliTOFAlignment",""),
fNTOFAlignObj(0),
fTOFmgr(0x0),
fTOFAlignObjArray(0x0)
{
for(Int_t i=0; i<18;i++)
for(Int_t j=0; j<5; j++)
fNFMforSM[i][j]=0;
for(Int_t i=0; i<72; i++)
for (Int_t j=0; j<6; j++)
fCombFMData[i][j]=0;
for(Int_t i=0; i<18;i++)
fTOFMatrixId[i]=0;
}
AliTOFAlignment::AliTOFAlignment(const AliTOFAlignment &t):
TTask(t),
fNTOFAlignObj(t.fNTOFAlignObj),
fTOFmgr(0x0),
fTOFAlignObjArray(t.fTOFAlignObjArray)
{
for(Int_t i=0; i<18;i++)
for(Int_t j=0; j<5; j++)
fNFMforSM[i][j]=t.fNFMforSM[i][j];
for(Int_t i=0; i<72; i++)
for (Int_t j=0; j<6; j++)
fCombFMData[i][j]=t.fCombFMData[i][j];
for(Int_t i=0; i<18;i++)
fTOFMatrixId[i]=t.fTOFMatrixId[i];
}
AliTOFAlignment& AliTOFAlignment::operator=(const AliTOFAlignment &t){
if (&t == this)
return *this;
TTask::operator=(t);
fNTOFAlignObj=t.fNTOFAlignObj;
fTOFmgr=t.fTOFmgr;
fTOFAlignObjArray=t.fTOFAlignObjArray;
for(Int_t i=0; i<18;i++)
fTOFMatrixId[i]=t.fTOFMatrixId[i];
return *this;
}
AliTOFAlignment::~AliTOFAlignment() {
delete fTOFAlignObjArray;
delete fTOFmgr;
}
void AliTOFAlignment::Smear(Float_t * const tr, Float_t * const rot)
{
fTOFAlignObjArray = new TObjArray(kMaxAlignObj);
Float_t dx, dy, dz;
Float_t dpsi, dtheta, dphi;
TRandom *rnd = new TRandom(1567);
Int_t nSMTOF = 18;
AliGeomManager::ELayerID iLayer = AliGeomManager::kInvalidLayer;
UShort_t iIndex=0;
UShort_t dvoluid = AliGeomManager::LayerToVolUID(iLayer,iIndex);
Int_t i;
const Int_t kSize=100;
Char_t path[kSize];
for (i = 0; i<nSMTOF ; i++) {
snprintf(path,kSize,"/ALIC_1/B077_1/BSEGMO%i_1/BTOF%i_1",i,i);
dx = (rnd->Gaus(0.,1.))*tr[0];
dy = (rnd->Gaus(0.,1.))*tr[1];
dz = (rnd->Gaus(0.,1.))*tr[2];
dpsi = rot[0];
dtheta = rot[1];
dphi = rot[2];
AliAlignObjParams *o =new AliAlignObjParams(path, dvoluid, dx, dy, dz, dpsi, dtheta, dphi, kTRUE);
fTOFAlignObjArray->Add(o);
}
fNTOFAlignObj=fTOFAlignObjArray->GetEntries();
AliInfo(Form("Number of Alignable Volumes: %d",fNTOFAlignObj));
delete rnd;
}
void AliTOFAlignment::Align(Float_t * const tr, Float_t * const rot)
{
fTOFAlignObjArray = new TObjArray(kMaxAlignObj);
Float_t dx, dy, dz;
Float_t dpsi, dtheta, dphi;
Int_t nSMTOF = 18;
AliGeomManager::ELayerID iLayer = AliGeomManager::kInvalidLayer;
UShort_t iIndex=0;
UShort_t dvoluid = AliGeomManager::LayerToVolUID(iLayer,iIndex);
const Int_t kSize=100;
Char_t path[kSize];
Int_t i;
for (i = 0; i<nSMTOF ; i++) {
snprintf(path,kSize,"/ALIC_1/B077_1/BSEGMO%i_1/BTOF%i_1",i,i);
dx = tr[0];
dy = tr[1];
dz = tr[2];
dpsi = rot[0];
dtheta = rot[1];
dphi = rot[2];
AliAlignObjParams *o =new AliAlignObjParams(path, dvoluid, dx, dy, dz, dpsi, dtheta, dphi, kTRUE);
fTOFAlignObjArray->Add(o);
}
fNTOFAlignObj=fTOFAlignObjArray->GetEntries();
AliInfo(Form("Number of Alignable Volumes: %d",fNTOFAlignObj));
}
void AliTOFAlignment::WriteParOnCDB(const Char_t *sel, Int_t minrun, Int_t maxrun)
{
AliCDBManager *man = AliCDBManager::Instance();
const Char_t *sel1 = "AlignPar" ;
const Int_t kSize=100;
Char_t out[kSize];
snprintf(out,kSize,"%s/%s",sel,sel1);
AliCDBId idTOFAlign(out,minrun,maxrun);
AliCDBMetaData *mdTOFAlign = new AliCDBMetaData();
mdTOFAlign->SetResponsible("TOF");
AliInfo(Form("Number of Alignable Volumes: %d",fNTOFAlignObj));
man->Put(fTOFAlignObjArray,idTOFAlign,mdTOFAlign);
}
void AliTOFAlignment::ReadParFromCDB(const Char_t *sel, Int_t nrun)
{
AliCDBManager *man = AliCDBManager::Instance();
const Char_t *sel1 = "AlignPar" ;
const Int_t kSize=100;
Char_t out[kSize];
snprintf(out,kSize,"%s/%s",sel,sel1);
AliCDBEntry *entry = man->Get(out,nrun);
if (!entry) {
AliError(Form("Failed to get entry: %s",out));
return;
}
fTOFAlignObjArray=(TObjArray*)entry->GetObject();
fNTOFAlignObj=fTOFAlignObjArray->GetEntries();
AliInfo(Form("Number of Alignable Volumes from CDB: %d",fNTOFAlignObj));
}
void AliTOFAlignment::WriteSimParOnCDB(const Char_t *sel, Int_t minrun, Int_t maxrun)
{
AliCDBManager *man = AliCDBManager::Instance();
const Char_t *sel1 = "AlignSimPar" ;
const Int_t kSize=100;
Char_t out[kSize];
snprintf(out,kSize,"%s/%s",sel,sel1);
AliCDBId idTOFAlign(out,minrun,maxrun);
AliCDBMetaData *mdTOFAlign = new AliCDBMetaData();
mdTOFAlign->SetResponsible("TOF");
AliInfo(Form("Number of Alignable Volumes: %d",fNTOFAlignObj));
man->Put(fTOFAlignObjArray,idTOFAlign,mdTOFAlign);
}
void AliTOFAlignment::ReadSimParFromCDB(const Char_t *sel, Int_t nrun){
AliCDBManager *man = AliCDBManager::Instance();
const Char_t *sel1 = "AlignSimPar" ;
const Int_t kSize=100;
Char_t out[kSize];
snprintf(out,kSize,"%s/%s",sel,sel1);
AliCDBEntry *entry = man->Get(out,nrun);
if (!entry) {
AliError(Form("Failed to get entry: %s",out));
return;
}
fTOFAlignObjArray=(TObjArray*)entry->GetObject();
fNTOFAlignObj=fTOFAlignObjArray->GetEntries();
AliInfo(Form("Number of Alignable Volumes from CDB: %d",fNTOFAlignObj));
}
void AliTOFAlignment::WriteOnCDBforDC()
{
AliCDBManager *man = AliCDBManager::Instance();
AliCDBId idTOFAlign("TOF/Align/Data",0,0);
AliCDBMetaData *mdTOFAlign = new AliCDBMetaData();
mdTOFAlign->SetComment("Alignment objects for ideal geometry, i.e. applying them to TGeo has to leave geometry unchanged");
mdTOFAlign->SetResponsible("TOF");
AliInfo(Form("Number of Alignable Volumes: %d",fNTOFAlignObj));
man->Put(fTOFAlignObjArray,idTOFAlign,mdTOFAlign);
}
void AliTOFAlignment::ReadFromCDBforDC()
{
AliCDBManager *man = AliCDBManager::Instance();
AliCDBEntry *entry = man->Get("TOF/Align/Data",0);
fTOFAlignObjArray=(TObjArray*)entry->GetObject();
fNTOFAlignObj=fTOFAlignObjArray->GetEntries();
AliInfo(Form("Number of Alignable Volumes from CDB: %d",fNTOFAlignObj));
}
void AliTOFAlignment::BuildGeomForSurvey()
{
fTOFmgr = new TGeoManager("Geom","survey to alignment for TOF");
TGeoMedium *medium = 0;
TGeoVolume *top = fTOFmgr->MakeBox("TOP",medium,1000,1000,1000);
fTOFmgr->SetTopVolume(top);
TGeoTrd1 *strd1 = new TGeoTrd1(fgkX1BTOF*0.5,fgkX2BTOF*0.5, fgkYBTOF*0.5,fgkZBTOF*0.5);
TGeoVolume* trd1[18];
TGeoBBox *fmbox = new TGeoBBox(1,1,1);
TGeoVolume* fm = new TGeoVolume("FM",fmbox);
fm->SetLineColor(2);
TGeoTranslation* mAtr = new TGeoTranslation("mAtr",-fgkXFM, -fgkYFM ,fgkZFM);
TGeoTranslation* mBtr = new TGeoTranslation("mBtr",fgkXFM, -fgkYFM ,fgkZFM );
TGeoTranslation* mCtr = new TGeoTranslation("mCtr",fgkXFM, fgkYFM ,fgkZFM );
TGeoTranslation* mDtr = new TGeoTranslation("mDtr",-fgkXFM, fgkYFM ,fgkZFM );
const Int_t kSize=100;
char name[kSize];
Double_t smX = 0.;
Double_t smY = 0.;
Double_t smZ = 0.;
Float_t smR = fgkRorigTOF;
for (Int_t iSM = 0; iSM < 18; iSM++) {
Int_t mod = iSM + 13;
if (mod > 17) mod -= 18;
snprintf(name,kSize, "BTOF%d",mod);
trd1[iSM] = new TGeoVolume(name,strd1);
Float_t phi = iSM * 20.;
Float_t phi2 = 270 + phi;
if (phi2 >= 360.) phi2 -= 360.;
smX = TMath::Sin(phi*TMath::Pi()/180.)*smR;
smY = -TMath::Cos(phi*TMath::Pi()/180.)*smR;
smZ = 0.;
TGeoRotation* bTOFRot = new TGeoRotation("bTOFRot",phi,90,0.);
TGeoCombiTrans trans = *(new TGeoCombiTrans(smX,smY,smZ, bTOFRot));
TGeoMatrix* id = new TGeoHMatrix();
TGeoHMatrix transMat = *id * trans;
TGeoHMatrix *smTrans = new TGeoHMatrix(transMat);
trd1[iSM]->AddNode(fm,1,mAtr);
trd1[iSM]->AddNode(fm,2,mBtr);
trd1[iSM]->AddNode(fm,3,mCtr);
trd1[iSM]->AddNode(fm,4,mDtr);
top->AddNode(trd1[iSM],1,smTrans);
trd1[iSM]->SetVisDaughters();
trd1[iSM]->SetLineColor(iSM);
}
fTOFmgr->CloseGeometry();
fTOFmgr->GetTopVolume()->Draw();
fTOFmgr->SetVisOption(0);
fTOFmgr->SetVisLevel(6);
for (Int_t iSM = 0; iSM < 18; iSM++) {
snprintf(name,kSize, "TOP_1/BTOF%d_1", iSM);
printf("\n\n***************** TOF SuperModule: %s ****************** \n",name);
TGeoPhysicalNode* pn3 = fTOFmgr->MakePhysicalNode(name);
fTOFMatrixId[iSM] = pn3->GetMatrix();
printf("\n\n*************** The Ideal Matrix in GRS *****************\n");
fTOFMatrixId[iSM]->Print();
}
}
void AliTOFAlignment::InsertMisAlignment(Float_t * const mis)
{
Double_t lA[3]={-fgkXFM, -fgkYFM ,fgkZFM};
Double_t lB[3]={fgkXFM, -fgkYFM ,fgkZFM};
Double_t lC[3]={fgkXFM, fgkYFM ,fgkZFM};
Double_t lD[3]={-fgkXFM, fgkYFM ,fgkZFM};
const Int_t kSize=16;
char name[kSize];
for(Int_t iSM=0;iSM<18;iSM++){
snprintf(name,kSize, "TOP_1/BTOF%d_1", iSM);
fTOFmgr->cd(name);
printf("\n\n******Misaligning TOF SuperModule ************** %s \n",name);
TGeoHMatrix g3 = *fTOFmgr->GetCurrentMatrix();
AliInfo(Form("This is the ideal global trasformation of SM %i",iSM));
g3.Print();
TGeoNode* n3 = fTOFmgr->GetCurrentNode();
TGeoMatrix* l3 = n3->GetMatrix();
Double_t gA[3], gB[3], gC[3], gD[3];
g3.LocalToMaster(lA,gA);
g3.LocalToMaster(lB,gB);
g3.LocalToMaster(lC,gC);
g3.LocalToMaster(lD,gD);
Double_t dx = mis[0];
Double_t dy = mis[1];
Double_t dz = mis[2];
Double_t dphi = mis[3];
Double_t dtheta = mis[4];
Double_t dpsi = mis[5];
TGeoRotation* rrot = new TGeoRotation("rot",dphi,dtheta,dpsi);
TGeoCombiTrans localdelta = *(new TGeoCombiTrans(dx,dy,dz, rrot));
AliInfo(Form("This is the local delta trasformation for SM %i \n",iSM));
localdelta.Print();
TGeoHMatrix nlocal = *l3 * localdelta;
TGeoHMatrix* nl3 = new TGeoHMatrix(nlocal);
TGeoPhysicalNode* pn3 = fTOFmgr->MakePhysicalNode(name);
pn3->Align(nl3);
TGeoHMatrix* ng3 = pn3->GetMatrix();
printf("\n\n************* The Misaligned Matrix in GRS **************\n");
ng3->Print();
Double_t ngA[3], ngB[3], ngC[3], ngD[3];
ng3->LocalToMaster(lA,ngA);
ng3->LocalToMaster(lB,ngB);
ng3->LocalToMaster(lC,ngC);
ng3->LocalToMaster(lD,ngD);
for(Int_t coord=0;coord<3;coord++){
fCombFMData[iSM*4][2*coord]=ngA[coord];
fCombFMData[iSM*4][2*coord+1]=1;
fCombFMData[iSM*4+1][2*coord]=ngB[coord];
fCombFMData[iSM*4+1][2*coord+1]=1;
fCombFMData[iSM*4+2][2*coord]=ngC[coord];
fCombFMData[iSM*4+2][2*coord+1]=1;
fCombFMData[iSM*4+3][2*coord]=ngD[coord];
fCombFMData[iSM*4+3][2*coord+1]=1;
}
}
}
void AliTOFAlignment::WriteCombData(const Char_t *nomefile, Int_t option)
{
FILE *data;
if( (data = fopen( nomefile, "w+t" )) != NULL ){
if (option==1){
fprintf( data, "simulated data\n" );} else {
fprintf( data, "survey data\n" );}
if (option==1){
fprintf( data, "data from InsertMisAlignmentBTOF method\n");}
else {fprintf( data, "real survey data from text file (coordinate in global RS)\n");}
fprintf( data, "Point Name,XPH,YPH,ZPH,PrecisionX(mm),PrecisionY(mm),PrecisionZ(mm)\n");
fprintf( data, "> Data:\n");
for(Int_t i=0;i<72;i++){
if (fCombFMData[i][0]!=0){
fprintf( data, "SM%02iFM%i %f %f %f M Y %f %f %f\n", (i-i%4)/4, i%4, fCombFMData[i][0],fCombFMData[i][2],fCombFMData[i][4],fCombFMData[i][1]*10,fCombFMData[i][3]*10,fCombFMData[i][5]*10);
}
}
fclose( data );
}
else{
printf( "Problem opening the file\n" );
}
return;
}
void AliTOFAlignment::WriteSimSurveyData(const Char_t *nomefile)
{
FILE *data;
if( (data = fopen( nomefile, "w+t" )) != NULL )
{
fprintf( data, "> Title:\n" );
fprintf( data, "simulated data\n" );
fprintf( data, "> Date:\n" );
fprintf( data, "24.09.2007\n" );
fprintf( data, "> Subdetector:\n" );
fprintf( data, "TOF\n" );
fprintf( data, "> Report URL:\n" );
fprintf( data, "https://edms.cern.ch/document/835615\n" );
fprintf( data, "> Version:\n" );
fprintf( data, "1\n");
fprintf( data, "> General Observations:\n");
fprintf( data, "data from InsertMisAlignmentBTOF method\n");
fprintf( data, "> Coordinate System:\n");
fprintf( data, "\\ALICEPH\n");
fprintf( data, "> Units:\n");
fprintf( data, "cm\n");
fprintf( data, "> Nr Columns:\n");
fprintf( data, "9\n");
fprintf( data, "> Column Names:\n");
fprintf( data, "Point Name,XPH,YPH,ZPH,Point Type,Target Used,PrecisionX(mm),PrecisionY(mm),PrecisionZ(mm)\n");
fprintf( data, "> Data:\n");
for(Int_t i=0;i<72;i++)
if (fCombFMData[i][0]!=0)
fprintf( data, "SM%02iFM%i %f %f %f M Y %f %f %f\n", (i-i%4)/4, i%4, fCombFMData[i][0],fCombFMData[i][2],fCombFMData[i][4],fCombFMData[i][1],fCombFMData[i][3],fCombFMData[i][5]);
fclose( data );
}
else
printf( "Problem opening the file\n" );
}
void AliTOFAlignment::MakeDefData(const Int_t nf,TString namefiles[])
{
Float_t data[72][6][100];
for (Int_t i=0;i<72;i++)
for (Int_t j=0; j<6; j++)
for(Int_t k=0; k<100; k++)
data[i][j][k]=0;
Int_t nfm=0;
Int_t nsm=0;
Long64_t totdata[72]={0};
for (Int_t ii=0;ii<nf; ii++)
{
AliSurveyObj *so = new AliSurveyObj();
const Char_t *nome=namefiles[ii];
so->FillFromLocalFile(nome);
TObjArray *points = so->GetData();
Int_t nSurveyPoint=points->GetEntries();
for(Int_t jj=0;jj<nSurveyPoint;jj++){
const char* pointName= ((AliSurveyPoint *) points->At(jj))->GetPointName().Data();
nfm=atoi(&pointName[6]);
nsm=atoi(&pointName[2]);
data[nsm*4+nfm][0][totdata[nsm*4+nfm]]=((AliSurveyPoint *) points->At(jj))->GetX();
data[nsm*4+nfm][2][totdata[nsm*4+nfm]]=((AliSurveyPoint *) points->At(jj))->GetY();
data[nsm*4+nfm][4][totdata[nsm*4+nfm]]=((AliSurveyPoint *) points->At(jj))->GetZ();
data[nsm*4+nfm][1][totdata[nsm*4+nfm]]=((AliSurveyPoint *) points->At(jj))->GetPrecisionX();
data[nsm*4+nfm][3][totdata[nsm*4+nfm]]=((AliSurveyPoint *) points->At(jj))->GetPrecisionY();
data[nsm*4+nfm][5][totdata[nsm*4+nfm]]=((AliSurveyPoint *) points->At(jj))->GetPrecisionZ();
totdata[nsm*4+nfm]=totdata[nsm*4+nfm]+1;
}
delete so;
}
for(Int_t i=0; i<72 ;i++){
Float_t numx=0, numy=0,numz=0, comodox=0, comodoy=0, comodoz=0,denx=0, deny=0, denz=0;
if(totdata[i]!=0){
for(Int_t j=0; j<totdata[i]; j++){
comodox=1/(data[i][1][j]/10*data[i][1][j]/10);
numx=numx+data[i][0][j]*comodox;
denx=denx+comodox;
comodoy=1/(data[i][3][j]/10*data[i][3][j]/10);
numy=numy+data[i][2][j]*comodoy;
deny=deny+comodoy;
comodoz=1/(data[i][5][j]/10*data[i][5][j]/10);
numz=numz+data[i][4][j]*comodoz;
denz=denz+comodoz;
}
fCombFMData[i][1]=TMath::Sqrt(1/denx);
fCombFMData[i][3]=TMath::Sqrt(1/deny);
fCombFMData[i][5]=TMath::Sqrt(1/denz);
fCombFMData[i][0]=numx/denx;
fCombFMData[i][2]=numy/deny;
fCombFMData[i][4]=numz/denz;
} else continue;
}
for(Int_t i=0;i<72;i++)
if (fCombFMData[i][0]!=0){
fNFMforSM[(i-i%4)/4][i%4]=1;
fNFMforSM[(i-i%4)/4][4]=fNFMforSM[(i-i%4)/4][4]+1;
}
}
void AliTOFAlignment::ReadSurveyDataAndAlign(){
fTOFAlignObjArray = new TObjArray(kMaxAlignObj);
Float_t deltaFM0=0, deltaFM1=0, deltaFM2=0, deltaFM3=0;
for(Int_t i=0; i<18; i++){
switch(fNFMforSM[i][4]){
case 0:
printf("we don't know the position of any FM of SM %i\n",i);
break;
case 1:
printf("we know the position of only one FM for SM %i\n",i);
break;
case 2:
printf("we know the position of only 2 FM for SM %i\n",i);
break;
case 3:
if (fNFMforSM[i][0]==1 && fNFMforSM[i][1]==1 && fNFMforSM[i][2]==1){
printf("we know the position of FM A B C for SM %i\n",i);
AliTOFAlignment::AlignFromSurveyABC(i);};
if (fNFMforSM[i][0]==1 && fNFMforSM[i][1]==1 && fNFMforSM[i][3]==1){
printf("we know the position of FM A B D for SM %i\n",i);
AliTOFAlignment::AlignFromSurveyABD(i);};
if (fNFMforSM[i][0]==1 && fNFMforSM[i][2]==1 && fNFMforSM[i][3]==1){
printf("we know the position of FM A C D for SM %i\n",i);
AliTOFAlignment::AlignFromSurveyACD(i);};
if (fNFMforSM[i][1]==1 && fNFMforSM[i][2]==1 && fNFMforSM[i][3]==1){
printf("we know the position of FM B C D for SM %i\n",i);
AliTOFAlignment::AlignFromSurveyBCD(i);};
break;
case 4:
printf("we know the position of all the 4 FM for SM %i\n",i);
deltaFM0=fCombFMData[i*4][1]/TMath::Abs(fCombFMData[i*4][0])+fCombFMData[i*4][3]/TMath::Abs(fCombFMData[i*4][2])+fCombFMData[i*4][5]/TMath::Abs(fCombFMData[i*4][4]);
deltaFM1=fCombFMData[i*4+1][1]/TMath::Abs(fCombFMData[i*4+1][0])+fCombFMData[i*4+1][3]/TMath::Abs(fCombFMData[i*4+1][2])+fCombFMData[i*4+1][5]/TMath::Abs(fCombFMData[i*4+1][4]);
deltaFM2=fCombFMData[i*4+2][1]/TMath::Abs(fCombFMData[i*4+2][0])+fCombFMData[i*4+2][3]/TMath::Abs(fCombFMData[i*4+2][2])+fCombFMData[i*4+2][5]/TMath::Abs(fCombFMData[i*4+2][4]);
deltaFM3=fCombFMData[i*4+3][1]/TMath::Abs(fCombFMData[i*4+3][0])+fCombFMData[i*4+3][3]/TMath::Abs(fCombFMData[i*4+3][2])+fCombFMData[i*4+3][5]/TMath::Abs(fCombFMData[i*4+3][4]);
if(deltaFM0>=deltaFM1 && deltaFM0>=deltaFM2 && deltaFM0>=deltaFM3){
printf("to Align we use FM B,C,D");
AliTOFAlignment::AlignFromSurveyBCD(i);} else
if(deltaFM1>=deltaFM0 && deltaFM1>=deltaFM2 && deltaFM1>=deltaFM3){
printf("to Align we use FM A,C,D");
AliTOFAlignment::AlignFromSurveyACD(i);} else
if(deltaFM2>=deltaFM0 && deltaFM2>=deltaFM1 && deltaFM2>=deltaFM3){
printf("to Align we use FM A,B,D");
AliTOFAlignment::AlignFromSurveyABD(i);} else{
printf("to Align we use FM A,B,C");
AliTOFAlignment::AlignFromSurveyABC(i);}
break;
}
}
fNTOFAlignObj=fTOFAlignObjArray->GetEntries();
AliInfo(Form("Number of Alignable Volumes: %d",fNTOFAlignObj));
TFile f("TOFAlignFromSurvey.root","RECREATE");
f.cd();
f.WriteObject(fTOFAlignObjArray,"TOFAlignObjs","kSingleKey");
f.Close();
}
void AliTOFAlignment::AlignFromSurveyABC(Int_t iSM)
{
Double_t ngA[3], ngB[3], ngC[3];
for(Int_t coord=0;coord<3;coord++){
ngA[coord]= fCombFMData[iSM*4][coord*2];
ngB[coord]= fCombFMData[iSM*4+1][coord*2];
ngC[coord]= fCombFMData[iSM*4+2][coord*2];
}
printf("\n\n******Survey analysis for TOF SuperModule ************** %i \n",iSM);
Double_t ab[3], bc[3], n[3];
Double_t plane[4], s=1.;
for(Int_t i=0;i<3;i++){
ab[i] = (ngB[i] - ngA[i]);
}
for(Int_t i=0;i<3;i++){
bc[i] = (ngC[i] - ngB[i]);
}
n[0] = (ab[1] * bc[2] - ab[2] * bc[1]);
n[1] = (ab[2] * bc[0] - ab[0] * bc[2]);
n[2] = (ab[0] * bc[1] - ab[1] * bc[0]);
Double_t sizen = TMath::Sqrt( n[0]*n[0] + n[1]*n[1] + n[2]*n[2] );
if(sizen>1.e-8){
s = Double_t(1.)/sizen ;
}else{
AliInfo("Problem in normalizing the vector");
}
for(Int_t i=0;i<3;i++){
plane[i] = n[i] * s;
}
plane[3] = ( plane[0] * ngA[0] + plane[1] * ngA[1] + plane[2] * ngA[2] );
Double_t orig[3], md[3];
for(Int_t i=0;i<3;i++){
md[i] = (ngA[i] + ngC[i]) * 0.5;
}
for(Int_t i=0;i<3;i++){
orig[i] = md[i] - plane[i]*fgkZFM;
}
Double_t sx = TMath::Sqrt(ab[0]*ab[0] + ab[1]*ab[1] + ab[2]*ab[2]);
if(sx>1.e-8){
for(Int_t i=0;i<3;i++){
ab[i] /= sx;
}
}
Double_t sy = TMath::Sqrt(bc[0]*bc[0] + bc[1]*bc[1] + bc[2]*bc[2]);
if(sy>1.e-8){
for(Int_t i=0;i<3;i++){
bc[i] /= sy;
}
}
Double_t rot[9] = {ab[0],bc[0],plane[0],ab[1],bc[1],plane[1],ab[2],bc[2],plane[2]};
TGeoHMatrix ng;
ng.SetTranslation(orig);
ng.SetRotation(rot);
printf("\n\n**** The Misaligned Matrix in GRS, as from Survey data ***\n");
ng.Print();
printf("\n\n**** The ideal matrix ***\n");
fTOFMatrixId[iSM]->Print();
TGeoHMatrix gdelta =fTOFMatrixId[iSM]->Inverse();
printf("\n\n**** The inverse of the ideal matrix ***\n");
gdelta.Print();
gdelta.MultiplyLeft(&ng);
printf("\n\n**** The Delta Matrix in GRS, as from Survey data ***\n");
gdelta.Print();
Int_t index=0;
AliGeomManager::ELayerID layer = AliGeomManager::kInvalidLayer;
UShort_t dvoluid = AliGeomManager::LayerToVolUID(layer,index);
TString symname(Form("TOF/sm%02d",iSM));
AliAlignObjMatrix* o = new AliAlignObjMatrix(symname.Data(),dvoluid,gdelta,kTRUE);
fTOFAlignObjArray->Add(o);
}
void AliTOFAlignment::AlignFromSurveyABD(Int_t iSM)
{
Double_t ngA[3], ngB[3], ngD[3];
for(Int_t coord=0;coord<3;coord++){
ngA[coord]= fCombFMData[iSM*4][coord*2];
ngB[coord]= fCombFMData[iSM*4+1][coord*2];
ngD[coord]= fCombFMData[iSM*4+3][coord*2];
}
printf("\n\n******Survey analysis for TOF SuperModule ************** %i \n",iSM);
Double_t ab[3], ad[3], n[3];
Double_t plane[4], s=1.;
for(Int_t i=0;i<3;i++){
ab[i] = (ngB[i] - ngA[i]);
}
for(Int_t i=0;i<3;i++){
ad[i] = (ngD[i] - ngA[i]);
}
n[0] = (ab[1] * ad[2] - ab[2] * ad[1]);
n[1] = (ab[2] * ad[0] - ab[0] * ad[2]);
n[2] = (ab[0] * ad[1] - ab[1] * ad[0]);
Double_t sizen = TMath::Sqrt( n[0]*n[0] + n[1]*n[1] + n[2]*n[2] );
if(sizen>1.e-8){
s = Double_t(1.)/sizen ;
}else{
AliInfo("Problem in normalizing the vector");
}
for(Int_t i=0;i<3;i++){
plane[i] = n[i] * s;
}
plane[3] = ( plane[0] * ngA[0] + plane[1] * ngA[1] + plane[2] * ngA[2] );
Double_t orig[3], md[3];
for(Int_t i=0;i<3;i++){
md[i] = (ngB[i] + ngD[i]) * 0.5;
}
for(Int_t i=0;i<3;i++){
orig[i] = md[i] - plane[i]*fgkZFM;
}
Double_t sx = TMath::Sqrt(ab[0]*ab[0] + ab[1]*ab[1] + ab[2]*ab[2]);
if(sx>1.e-8){
for(Int_t i=0;i<3;i++){
ab[i] /= sx;
}
}
Double_t sy = TMath::Sqrt(ad[0]*ad[0] + ad[1]*ad[1] + ad[2]*ad[2]);
if(sy>1.e-8){
for(Int_t i=0;i<3;i++){
ad[i] /= sy;
}
}
Double_t rot[9] = {ab[0],ad[0],plane[0],ab[1],ad[1],plane[1],ab[2],ad[2],plane[2]};
TGeoHMatrix ng;
ng.SetTranslation(orig);
ng.SetRotation(rot);
printf("\n\n**** The Misaligned Matrix in GRS, as from Survey data ***\n");
ng.Print();
printf("\n\n**** The ideal matrix ***\n");
fTOFMatrixId[iSM]->Print();
TGeoHMatrix gdelta =fTOFMatrixId[iSM]->Inverse();
printf("\n\n**** The inverse of the ideal matrix ***\n");
gdelta.Print();
gdelta.MultiplyLeft(&ng);
printf("\n\n**** The Delta Matrix in GRS, as from Survey data ***\n");
gdelta.Print();
Int_t index=0;
AliGeomManager::ELayerID layer = AliGeomManager::kInvalidLayer;
UShort_t dvoluid = AliGeomManager::LayerToVolUID(layer,index);
TString symname(Form("TOF/sm%02d",iSM));
AliAlignObjMatrix* o = new AliAlignObjMatrix(symname.Data(),dvoluid,gdelta,kTRUE);
fTOFAlignObjArray->Add(o);
}
void AliTOFAlignment::AlignFromSurveyACD(Int_t iSM)
{
Double_t ngA[3], ngC[3], ngD[3];
for(Int_t coord=0;coord<3;coord++){
ngA[coord]= fCombFMData[iSM*4][coord*2];
ngC[coord]= fCombFMData[iSM*4+2][coord*2];
ngD[coord]= fCombFMData[iSM*4+3][coord*2];
}
printf("\n\n******Survey analysis for TOF SuperModule ************** %i \n",iSM);
Double_t cd[3], ad[3], n[3];
Double_t plane[4], s=1.;
for(Int_t i=0;i<3;i++){
cd[i] = (ngC[i] - ngD[i]);
}
for(Int_t i=0;i<3;i++){
ad[i] = (ngD[i] - ngA[i]);
}
n[0] = (ad[1] * cd[2] - ad[2] * cd[1]);
n[1] = (ad[2] * cd[0] - ad[0] * cd[2]);
n[2] = (ad[0] * cd[1] - ad[1] * cd[0]);
Double_t sizen = TMath::Sqrt( n[0]*n[0] + n[1]*n[1] + n[2]*n[2] );
if(sizen>1.e-8){
s = Double_t(1.)/sizen ;
}else{
AliInfo("Problem in normalizing the vector");
}
for(Int_t i=0;i<3;i++){
plane[i] = n[i] * s;
}
plane[3] = ( plane[0] * ngA[0] + plane[1] * ngA[1] + plane[2] * ngA[2] );
Double_t orig[3], md[3];
for(Int_t i=0;i<3;i++){
md[i] = (ngA[i] + ngC[i]) * 0.5;
}
for(Int_t i=0;i<3;i++){
orig[i] = md[i] + plane[i]*fgkZFM;
}
Double_t sx = TMath::Sqrt(ad[0]*ad[0] + ad[1]*ad[1] + ad[2]*ad[2]);
if(sx>1.e-8){
for(Int_t i=0;i<3;i++){
ad[i] /= sx;
}
}
Double_t sy = TMath::Sqrt(cd[0]*cd[0] + cd[1]*cd[1] + cd[2]*cd[2]);
if(sy>1.e-8){
for(Int_t i=0;i<3;i++){
cd[i] /= sy;
}
}
Double_t rot[9] = {cd[0],ad[0],-plane[0],cd[1],ad[1],-plane[1],cd[2],ad[2],-plane[2]};
TGeoHMatrix ng;
ng.SetTranslation(orig);
ng.SetRotation(rot);
printf("\n\n**** The Misaligned Matrix in GRS, as from Survey data ***\n");
ng.Print();
printf("\n\n**** The ideal matrix ***\n");
fTOFMatrixId[iSM]->Print();
TGeoHMatrix gdelta =fTOFMatrixId[iSM]->Inverse();
printf("\n\n**** The inverse of the ideal matrix ***\n");
gdelta.Print();
gdelta.MultiplyLeft(&ng);
printf("\n\n**** The Delta Matrix in GRS, as from Survey data ***\n");
gdelta.Print();
Int_t index=0;
AliGeomManager::ELayerID layer = AliGeomManager::kInvalidLayer;
UShort_t dvoluid = AliGeomManager::LayerToVolUID(layer,index);
TString symname(Form("TOF/sm%02d",iSM));
AliAlignObjMatrix* o = new AliAlignObjMatrix(symname.Data(),dvoluid,gdelta,kTRUE);
fTOFAlignObjArray->Add(o);
}
void AliTOFAlignment::AlignFromSurveyBCD(Int_t iSM)
{
Double_t ngB[3], ngC[3], ngD[3];
for(Int_t coord=0;coord<3;coord++){
ngB[coord]= fCombFMData[iSM*4+1][coord*2];
ngC[coord]= fCombFMData[iSM*4+2][coord*2];
ngD[coord]= fCombFMData[iSM*4+3][coord*2];
}
printf("\n\n******Survey analysis for TOF SuperModule ************** %i \n",iSM);
Double_t cd[3], bc[3], n[3];
Double_t plane[4], s=1.;
for(Int_t i=0;i<3;i++){
cd[i] = (ngC[i] - ngD[i]);
}
for(Int_t i=0;i<3;i++){
bc[i] = (ngC[i] - ngB[i]);
}
n[0] = (bc[1] * cd[2] - bc[2] * cd[1]);
n[1] = (bc[2] * cd[0] - bc[0] * cd[2]);
n[2] = (bc[0] * cd[1] - bc[1] * cd[0]);
Double_t sizen = TMath::Sqrt( n[0]*n[0] + n[1]*n[1] + n[2]*n[2] );
if(sizen>1.e-8){
s = Double_t(1.)/sizen ;
}else{
AliInfo("Problem in normalizing the vector");
}
for(Int_t i=0;i<3;i++){
plane[i] = n[i] * s;
}
plane[3] = ( plane[0] * ngB[0] + plane[1] * ngB[1] + plane[2] * ngB[2] );
Double_t orig[3], md[3];
for(Int_t i=0;i<3;i++){
md[i] = (ngB[i] + ngD[i]) * 0.5;
}
for(Int_t i=0;i<3;i++){
orig[i] = md[i] + plane[i]*fgkZFM;
}
Double_t sx = TMath::Sqrt(cd[0]*cd[0] + cd[1]*cd[1] + cd[2]*cd[2]);
if(sx>1.e-8){
for(Int_t i=0;i<3;i++){
cd[i] /= sx;
}
}
Double_t sy = TMath::Sqrt(bc[0]*bc[0] + bc[1]*bc[1] + bc[2]*bc[2]);
if(sy>1.e-8){
for(Int_t i=0;i<3;i++){
bc[i] /= sy;
}
}
Double_t rot[9] = {cd[0],bc[0],-plane[0],cd[1],bc[1],-plane[1],cd[2],bc[2],-plane[2]};
TGeoHMatrix ng;
ng.SetTranslation(orig);
ng.SetRotation(rot);
printf("\n\n**** The Misaligned Matrix in GRS, as from Survey data ***\n");
ng.Print();
printf("\n\n**** The ideal matrix ***\n");
fTOFMatrixId[iSM]->Print();
TGeoHMatrix gdelta =fTOFMatrixId[iSM]->Inverse();
printf("\n\n**** The inverse of the ideal matrix ***\n");
gdelta.Print();
gdelta.MultiplyLeft(&ng);
printf("\n\n**** The Delta Matrix in GRS, as from Survey data ***\n");
gdelta.Print();
Int_t index=0;
AliGeomManager::ELayerID layer = AliGeomManager::kInvalidLayer;
UShort_t dvoluid = AliGeomManager::LayerToVolUID(layer,index);
TString symname(Form("TOF/sm%02d",iSM));
AliAlignObjMatrix* o = new AliAlignObjMatrix(symname.Data(),dvoluid,gdelta,kTRUE);
fTOFAlignObjArray->Add(o);
}