1 /**************************************************************************
2 * Copyright(c) 2008-2010, 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 to convert survey tables in alignment objects
21 // origin: Marco Van Leeuwen (m.vanleeuwen1@uu.nl)
22 // Panos.Christakoglou (Panos.Christakoglou@cern.ch)
23 // Martin Poghosyan (Martin.Poghosyan@to.infn.it)
24 //////////////////////////////////////////////////////////////////////////
26 #include "Riostream.h"
29 #include "TClonesArray.h"
30 #include "TGeoManager.h"
31 #include "TGeoMatrix.h"
32 #include "TGeoPhysicalNode.h"
36 #include "AliITSSurveyToAlign.h"
37 #include "AliSurveyObj.h"
38 #include "AliSurveyPoint.h"
39 #include "AliAlignObjParams.h"
40 #include "AliGeomManager.h"
44 #include "AliCDBManager.h"
45 #include "AliCDBEntry.h"
46 #include "AliCDBStorage.h"
48 #include "AliITSgeomTGeo.h"
50 ClassImp(AliITSSurveyToAlign)
52 const Double_t AliITSSurveyToAlign::fgkLocR[6][3]={{ 3.24,0.21905,-2.4},
57 {-3.24,0.21905,-2.4}};
59 const Double_t AliITSSurveyToAlign::fgkLocL[6][3]={{-3.24,0.21905, 2.4},
64 { 3.24,0.21905, 2.4}};
67 //________________________________________________________________________
68 AliITSSurveyToAlign::AliITSSurveyToAlign(Int_t run, Int_t repSDD, Int_t repVerSDD, Int_t repModSSD, Int_t repModVerSSD, Int_t repLaddSSD, Int_t repLaddVerSSD) :
70 AliSurveyToAlignObjs(),
72 fSDDrepNumber(repSDD),
73 fSDDrepVersion(repVerSDD),
74 fSSDModuleRepNumber(repModSSD),
75 fSSDModuleRepVersion(repModVerSSD),
76 fSSDLadderRepNumber(repLaddSSD),
77 fSSDLadderRepVersion(repLaddVerSSD)
80 // default constructor
81 // Arguments are report numbers for survey data.
82 // The defaults point to reports from detector construction
86 //_________________________________________________________________________
87 AliITSSurveyToAlign::AliITSSurveyToAlign(const AliITSSurveyToAlign &align) :
89 AliSurveyToAlignObjs(align),
91 fSDDrepNumber(align.fSDDrepNumber),
92 fSDDrepVersion(align.fSDDrepVersion),
93 fSSDModuleRepNumber(align.fSSDModuleRepNumber),
94 fSSDModuleRepVersion(align.fSSDModuleRepVersion),
95 fSSDLadderRepNumber(align.fSSDLadderRepNumber),
96 fSSDLadderRepVersion(align.fSSDLadderRepVersion)
103 //__________________________________________________________________________
104 AliITSSurveyToAlign & AliITSSurveyToAlign::operator =(const AliITSSurveyToAlign & align) {
106 // assignment operator - dummy
112 //__________________________________________________________________________
113 AliITSSurveyToAlign::~AliITSSurveyToAlign() {
119 //______________________________________________________________________
120 void AliITSSurveyToAlign::Run() {
122 // Runs the full chain
123 // User should call StoreAlignObj() or StoreAlignObjCDB() afterwards to
127 // Load ideal geometry from the OCDB
128 AliCDBManager *cdb = AliCDBManager::Instance();
129 cdb->SetDefaultStorage("local://$ALICE_ROOT/OCDB");
131 AliGeomManager::LoadGeometry();
134 LoadSurveyFromAlienFile("ITS", fSDDrepNumber, fSDDrepVersion);
136 LoadSurveyFromAlienFile("ITS", fSSDModuleRepNumber, fSSDModuleRepVersion);
137 CreateAlignObjSSDModules();
138 LoadSurveyFromAlienFile("ITS", fSSDLadderRepNumber, fSSDLadderRepVersion);
139 CreateAlignObjSSDLadders();
142 void AliITSSurveyToAlign::CreateAlignObjSPD(){
144 // Create alignObjs for SPD
145 // For the moment, uses 0,0,0,0,0,0
147 for(Int_t imod = 0; imod < 240; imod++) {
148 Int_t ilayer = (imod < 80) ? AliGeomManager::kSPD1 : AliGeomManager::kSPD2;
149 Int_t imodule = (imod < 80) ? imod : imod - 80;
151 Int_t uid = AliGeomManager::LayerToVolUID(ilayer,imodule);
152 const Char_t *symname = AliGeomManager::SymName(uid);
154 new((*fAlignObjArray)[imod]) AliAlignObjParams(symname, uid, 0., 0., 0., 0., 0., 0., kTRUE);
159 void AliITSSurveyToAlign::CreateAlignObjSDD(){
161 // Create alignment objects for SDD
165 const char* symname = 0;
166 AliSurveyPoint* pt = 0;
168 Int_t iModuleIndex=240;
174 if (fSurveyPoints == 0 || fSurveyPoints->GetEntries() == 0) {
175 AliWarning("SDD survey data are not available, using zero values");
176 CreateAlignObjDummySDD();
180 for(Int_t imod = 1; imod < fSurveyPoints->GetEntries(); imod++) {
181 pt = (AliSurveyPoint*) fSurveyPoints->At(imod);
184 Int_t iLayer, iLadder, iModule, iPoint;
185 ReadPointNameSDD(pt->GetName(),iLayer, iLadder, iModule, iPoint);
187 if(iModule==iModule0)
189 fSDDmeP[iPoint][0]=pt->GetX();
190 fSDDmeP[iPoint][1]=pt->GetY();
191 fSDDmeP[iPoint][2]=pt->GetZ();
192 fSDDmeP[iPoint][3]=pt->GetPrecisionX();
193 fSDDmeP[iPoint][4]=pt->GetPrecisionY();
194 fSDDmeP[iPoint][5]=pt->GetPrecisionZ();
195 fSDDisMe[iPoint]=kTRUE;
197 if(iLayer==3) uid = AliGeomManager::LayerToVolUID(iLayer0,iModuleIndex-240);
198 if(iLayer==4) uid = AliGeomManager::LayerToVolUID(iLayer0,iModuleIndex-324);
199 symname = AliGeomManager::SymName(uid);
200 GetIdPosSDD(uid,iLayer0, iModule0, iPoint);
203 // cout << "Points red module " << imod << endl;
204 if((iModule!=iModule0)||(imod==(fSurveyPoints->GetEntries()-1)))
206 ConvertToRSofModulesAndRotSDD(iLayer0, iModule0);
215 if(nModules==2) CalcShiftSDD(x0,y0,z0);
216 if(nModules>2) CalcShiftRotSDD(tet, psi, phi, x0, y0, z0);
217 tet*=(180/TMath::Pi());
218 psi*=(180/TMath::Pi());
219 phi*=(180/TMath::Pi());
221 // printf("%s %d %f %f %f %f %f %f\n",symname, uid, x0/10., y0/10., z0/10., psi, tet, phi);
222 // cout << "Allocate alignobjparams " << imod << endl;
223 new((*fAlignObjArray)[iModuleIndex]) AliAlignObjParams(symname, uid, x0/10., y0/10., z0/10., psi, tet, phi, kFALSE);
229 iModuleIndex = AliITSgeomTGeo::GetModuleIndex(iLayer,iLadder+1,iModule+1);
230 for(Int_t i=0; i<6;i++) fSDDisMe[i]=kFALSE;
231 if(imod!=(fSurveyPoints->GetEntries()-1)) imod--;
236 void AliITSSurveyToAlign::CreateAlignObjDummySDD(){
238 // Create empty alignment objects
239 // Used when fSurveySDD == 0
241 for(Int_t imod = 0; imod < 260; imod++) {
243 Int_t ilayer = (imod < 84) ? AliGeomManager::kSDD1 : AliGeomManager::kSDD2;
244 Int_t imodule = (imod < 84) ? imod : imod - 84;
246 Int_t uid = AliGeomManager::LayerToVolUID(ilayer,imodule);
247 const Char_t *symname = AliGeomManager::SymName(uid);
249 new((*fAlignObjArray)[imod+240]) AliAlignObjParams(symname, uid, 0., 0., 0., 0., 0., 0., kTRUE);
253 void AliITSSurveyToAlign::CreateAlignObjSSDModules(){
255 // Create alignment objects for SDD
257 // Two sets of objects are made: 1 per module + 1 per ladder
260 const Float_t kMu2Cm = 1e-4;
261 const Float_t kSensLength = 7.464;
262 const Int_t kSSDMODULES = 1698;
264 if (fSurveyPoints == 0 || fSurveyPoints->GetEntries() == 0) {
265 AliWarning("SSD module survey data not available; using dummy values");
266 CreateAlignObjDummySSDModules();
270 // First do module-by-module
272 for(Int_t imod = 500; imod < kSSDMODULES + 500; imod++) {
273 Int_t iLayer, iLadder, iLaddMod;
274 AliITSgeomTGeo::GetModuleId(imod,iLayer,iLadder,iLaddMod); // returns 1-based numbers
276 TString pname="ITS/SSD";
282 AliSurveyPoint *pt1 = (AliSurveyPoint*) fSurveyPoints->FindObject(pname+"/Point0");
283 AliSurveyPoint *pt2 = (AliSurveyPoint*) fSurveyPoints->FindObject(pname+"/Point1");
285 AliWarning(Form("No Survey points for iladd %d imod %d",iLadder,iLaddMod));
289 sx = 0.5*(pt1->GetX() + pt2->GetX()) * kMu2Cm;
290 sz = 0.5*(pt1->GetZ() + pt2->GetZ()) * kMu2Cm;
292 // Minus sign to change local coordinate convention
293 Float_t theta = -(pt2->GetZ() - pt1->GetZ())*kMu2Cm/kSensLength;
295 theta *= 180./TMath::Pi();
296 Int_t iLayMod = imod - 500;
299 Int_t uid = AliGeomManager::LayerToVolUID(iLayer,iLayMod);
301 const Char_t *symname = AliGeomManager::SymName(uid);
302 if (pname.CompareTo(symname) != 0)
303 AliWarning(Form("Mapping mismatch survey point %s volume name %s",pname.Data(),symname));
305 if (imod >= 676 && imod <= 697) {
306 cout << "ilayer " << iLayer << " imod " << imod
307 << " uid " << uid << " name " << symname
308 << " survey shift " << sx << " " << 0 << " " << sz << endl
309 << " theta " << theta << endl;
312 new((*fAlignObjArray)[imod]) AliAlignObjParams(symname, uid, sx, 0, sz, 0., theta, 0., kFALSE);
316 void AliITSSurveyToAlign::CreateAlignObjSSDLadders(){
318 // Now do ladders (Torino data)
320 const Float_t kLaddLen5 = 90.27; // Layer 5: distance between mouting points
321 const Float_t kLaddLen6 = 102.0; // Layer 6: distance between mouting points
322 const Float_t Z0 = 2.927; // Distance between V mounting point and Zloc = 0
323 // = half ladder length - nom z-position of ladder from gGeoManager
324 const Float_t kMu2Cm = 1e-4;
326 TString ssdName = "ITS/SSD";
328 TObjArray *ladd_points = fSurveyPoints;
329 if (ladd_points == 0 || ladd_points->GetEntries() == 0) {
330 AliWarning("No SSD Ladder alignment points found. Skipping");
333 if (ladd_points->GetEntries()!= 2*(34+38)) {
334 AliWarning(Form("Unexpected number of survey points %d, should be 144",ladd_points->GetEntries()));
337 for (Int_t ilayer = 4; ilayer <= 5; ilayer ++) {
338 Int_t nLadder = 34; // layer 5
340 nLadder = 38; // layer 6
342 for (Int_t iLadder = 0; iLadder < nLadder; iLadder++) {
343 TString ladName = ssdName;
345 ladName += "/Ladder";
348 AliSurveyPoint *v_point = (AliSurveyPoint*) ladd_points->FindObject(ladName+"/V");
349 AliSurveyPoint *q_point = (AliSurveyPoint*) ladd_points->FindObject(ladName+"/Q");
351 AliWarning(Form("Cannot find V side point for ladder %s",ladName.Data()));
355 AliWarning(Form("Cannot find Q side point for ladder %s",ladName.Data()));
360 tmp_str.Insert(0,v_point->GetName(),3);
361 Int_t ladder = tmp_str.Atoi();
363 tmp_str.Insert(0,q_point->GetName(),3);
364 if (tmp_str.Atoi() != ladder)
365 AliError(Form("Survey data file error. Expect pairs of V,Q points. Got ladders %d %d",ladder,tmp_str.Atoi()));
367 // Note: file gives meas-nom in local offline coordinates,
368 // ie. local z = - global z and local x = - global x (for ladder 508, i.e. top ladder)
369 Double_t dx_loc = v_point->GetX() * kMu2Cm;
370 Double_t dy_loc = v_point->GetY() * kMu2Cm;
371 Double_t dz_loc = v_point->GetZ() * kMu2Cm;
374 Double_t phi = 0; // Not measured
379 // Note: local psi = -global psi, psi = atan(-(y(z1) - y(z0)) / (z1-z0))
380 // local theta = global theta = atan(dx/dz)
381 // V side is A side is large global z
382 // Q side is C side is large local z
385 theta = TMath::ATan((q_point->GetX() - v_point->GetX())*kMu2Cm/kLaddLen6);
386 psi = TMath::ATan((v_point->GetY() - q_point->GetY())*kMu2Cm/kLaddLen6);
389 theta = TMath::ATan((q_point->GetX() - v_point->GetX())*kMu2Cm/kLaddLen5);
390 psi = TMath::ATan((v_point->GetY() - q_point->GetY())*kMu2Cm/kLaddLen5);
393 // Move along ladder to local Z = 0 point
397 // Convert to degrees
398 theta *= 180./TMath::Pi();
399 psi *= 180./TMath::Pi();
400 AliDebug(1,Form("ladname %f %f %f %f %f %f ",dx_loc,dy_loc,dz_loc,psi,theta,phi));
402 new((*fAlignObjArray)[500+1698+iLadd]) AliAlignObjParams(ladName,0,dx_loc,dy_loc,dz_loc,psi,theta,phi,kFALSE);
410 void AliITSSurveyToAlign::CreateAlignObjDummySSDModules(){
412 // Create empty alignment objects
413 // Used when fSurveySSD == 0
415 for(Int_t imod = 0; imod < 1698; imod++) {
416 Int_t ilayer = (imod < 748) ? AliGeomManager::kSSD1 : AliGeomManager::kSSD2;
417 Int_t imodule = (imod < 748) ? imod : imod - 748;
419 Int_t uid = AliGeomManager::LayerToVolUID(ilayer,imodule);
420 const Char_t *symname = AliGeomManager::SymName(uid);
422 new((*fAlignObjArray)[500+imod]) AliAlignObjParams(symname, uid, 0., 0., 0., 0., 0., 0., kTRUE);
427 void AliITSSurveyToAlign::GetIdPosSDD(Int_t uid, Int_t layer, Int_t module, Int_t iPoint)
430 // Utility function used by CreateAlignObjSDD
432 TGeoHMatrix gMod = *AliGeomManager::GetMatrix(uid); //global matrix of sensor
433 TGeoPNEntry* pne = gGeoManager->GetAlignableEntryByUID(uid);
434 // TString ladderPath = AliGeomManager::SymName(uid);
435 TString ladderPath(pne->GetTitle());
436 if(ladderPath.EndsWith("/")) ladderPath.Remove(TString::kTrailing,'/');
437 ladderPath.Remove(ladderPath.Last('/'));
438 ladderPath.Remove(ladderPath.Last('/'));
439 gGeoManager->cd(ladderPath.Data());
440 TGeoHMatrix gLad = *gGeoManager->GetCurrentMatrix(); // global matrix of ladder
441 TGeoHMatrix rel = gMod; // to equal relative matrix ladder to sensor.
442 TGeoHMatrix invgLad = gLad.Inverse();
443 rel.MultiplyLeft(&invgLad);
444 TGeoRotation* rr = new TGeoRotation("rr",90,90,0,0,90,180);
445 TGeoCombiTrans* ct = 0;
446 if(layer==3) ct= new TGeoCombiTrans(25.,0.,0.,rr);
447 if(layer==4) ct= new TGeoCombiTrans(25.+7.5,0.,0.,rr);
449 rel.MultiplyLeft(ct);
451 if((layer==3)&&(module<3)) rel.LocalToMaster(fgkLocR[iPoint],fSDDidP[iPoint]);
452 if((layer==3)&&(module>2)) rel.LocalToMaster(fgkLocL[iPoint],fSDDidP[iPoint]);
453 if((layer==4)&&(module<4)) rel.LocalToMaster(fgkLocR[iPoint],fSDDidP[iPoint]);
454 if((layer==4)&&(module>3)) rel.LocalToMaster(fgkLocL[iPoint],fSDDidP[iPoint]);
456 for(Int_t i=0; i<3; i++) fSDDidP[iPoint][i]*=10;
460 void AliITSSurveyToAlign::ReadPointNameSDD(const char str[], Int_t &iLayer, Int_t &iLader, Int_t &iModul, Int_t &iPoint)
463 // Utility function used by CreateAlignObjSDD
470 if(str[7]=='2') iLayer=3;
471 if(str[7]=='3') iLayer=4;
473 if(str[15]=='0') iLader=0;
474 if(str[15]=='1') iLader=1;
475 if(str[15]=='2') iLader=2;
476 if(str[15]=='3') iLader=3;
477 if(str[15]=='4') iLader=4;
478 if(str[15]=='5') iLader=5;
479 if(str[15]=='6') iLader=6;
480 if(str[15]=='7') iLader=7;
481 if(str[15]=='8') iLader=8;
482 if(str[15]=='9') iLader=9;
485 if(str[16]=='0') {iLader=10*iLader+0; ord=1;}
486 if(str[16]=='1') {iLader=10*iLader+1; ord=1;}
487 if(str[16]=='2') {iLader=10*iLader+2; ord=1;}
488 if(str[16]=='3') {iLader=10*iLader+3; ord=1;}
489 if(str[16]=='4') {iLader=10*iLader+4; ord=1;}
490 if(str[16]=='5') {iLader=10*iLader+5; ord=1;}
491 if(str[16]=='6') {iLader=10*iLader+6; ord=1;}
492 if(str[16]=='7') {iLader=10*iLader+7; ord=1;}
493 if(str[16]=='8') {iLader=10*iLader+8; ord=1;}
494 if(str[16]=='9') {iLader=10*iLader+9; ord=1;}
496 if(str[23+ord]=='0') iModul=0;
497 if(str[23+ord]=='1') iModul=1;
498 if(str[23+ord]=='2') iModul=2;
499 if(str[23+ord]=='3') iModul=3;
500 if(str[23+ord]=='4') iModul=4;
501 if(str[23+ord]=='5') iModul=5;
502 if(str[23+ord]=='6') iModul=6;
503 if(str[23+ord]=='7') iModul=7;
504 if(str[23+ord]=='8') iModul=8;
505 if(str[23+ord]=='9') iModul=9;
507 if((str[25+ord]=='R')&&(str[26+ord]=='D')) iPoint=0;
508 if((str[25+ord]=='R')&&(str[26+ord]=='C')) iPoint=1;
509 if((str[25+ord]=='R')&&(str[26+ord]=='U')) iPoint=2;
510 if((str[25+ord]=='L')&&(str[26+ord]=='U')) iPoint=3;
511 if((str[25+ord]=='L')&&(str[26+ord]=='C')) iPoint=4;
512 if((str[25+ord]=='L')&&(str[26+ord]=='D')) iPoint=5;
517 void AliITSSurveyToAlign::ConvertToRSofModulesAndRotSDD(Int_t Layer, Int_t Module)
520 // Utility function used by CreateAlignObjSDD
531 Double_t x0=fSDDidP[1][0];
532 Double_t z0=fSDDidP[1][2]-0.52;
533 for(Int_t i=0; i<6; i++)
537 if(!fSDDisMe[i]) continue;
547 fSDDidP[i][2]=fSDDidP[i][0];
557 fSDDmeP[i][2]=fSDDmeP[i][0];
560 fSDDmeP[i][5]=fSDDmeP[i][3];
565 if(((Layer==3)&&(Module>2))||((Layer==4)&&(Module>3)))
576 void AliITSSurveyToAlign::CalcShiftSDD(Double_t &x0,Double_t &y0,Double_t &z0)
578 Double_t Xid, Yid, Zid;
579 Double_t Xme, Yme, Zme, sX2, sY2, sZ2;
580 Double_t aX=0., bX=0.;
581 Double_t aY=0., bY=0.;
582 Double_t aZ=0., bZ=0.;
583 for(Int_t iP1=0; iP1<6; iP1++)
585 if(!fSDDisMe[iP1]) continue;
592 sX2 =fSDDmeP[iP1][3]*fSDDmeP[iP1][3];
593 sY2 =fSDDmeP[iP1][4]*fSDDmeP[iP1][4];
594 sZ2 =fSDDmeP[iP1][5]*fSDDmeP[iP1][5];
612 void AliITSSurveyToAlign::CalcShiftRotSDD(Double_t &tet,Double_t &psi,Double_t &phi,Double_t &x0,Double_t &y0,Double_t &z0)
671 Double_t Xid, Yid, Zid;
672 Double_t Xme, Yme, Zme, sX2, sY2, sZ2;
674 for(Int_t iP1=0; iP1<=6; iP1++)
676 if(!fSDDisMe[iP1]) continue;
678 Xid= fSDDidP[iP1][0];
679 Yid= fSDDidP[iP1][1];
680 Zid= fSDDidP[iP1][2];
682 Xme= fSDDmeP[iP1][0];
683 Yme= fSDDmeP[iP1][1];
684 Zme= fSDDmeP[iP1][2];
686 sX2 = fSDDmeP[iP1][3]* fSDDmeP[iP1][3];
687 sY2 = fSDDmeP[iP1][4]* fSDDmeP[iP1][4];
688 sZ2 = fSDDmeP[iP1][5]* fSDDmeP[iP1][5];
690 a11+=(Zid*Zid/sX2+Xid*Xid/sZ2);
696 c1+=(Xid*(Zme-Zid)/sZ2-Zid*(Xme-Xid)/sX2);
699 a22+=(Xid*Xid/sY2+Yid*Yid/sX2);
704 c2+=(Yid*(Xme-Xid)/sX2-Xid*(Yme-Yid)/sY2);
708 a33+=(Zid*Zid/sY2+Yid*Yid/sZ2);
712 c3+=(Zid*(Yme-Yid)/sY2-Yid*(Zme-Zid)/sZ2);
739 ///////////////////////////////////////////////////////////////
741 Double_t dataC[]={a11, a12, a13, a14, a15, a16, a21, a22, a23, a24, a25, a26,
742 a31, a32, a33, a34, a35, a36, a41, a42, a43, a44, a45, a46,
743 a51, a52, a53, a54, a55, a56, a61, a62, a63, a64, a65, a66};
744 Double_t data1[]={c1, a12, a13, a14, a15, a16, c2, a22, a23, a24, a25, a26,
745 c3, a32, a33, a34, a35, a36, c4, a42, a43, a44, a45, a46,
746 c5, a52, a53, a54, a55, a56, c6, a62, a63, a64, a65, a66};
747 Double_t data2[]={a11, c1, a13, a14, a15, a16, a21, c2, a23, a24, a25, a26,
748 a31, c3, a33, a34, a35, a36, a41, c4, a43, a44, a45, a46,
749 a51, c5, a53, a54, a55, a56, a61, c6, a63, a64, a65, a66};
750 Double_t data3[]={a11, a12, c1, a14, a15, a16, a21, a22, c2, a24, a25, a26,
751 a31, a32, c3, a34, a35, a36, a41, a42, c4, a44, a45, a46,
752 a51, a52, c5, a54, a55, a56, a61, a62, c6, a64, a65, a66};
753 Double_t data4[]={a11, a12, a13, c1, a15, a16, a21, a22, a23, c2, a25, a26,
754 a31, a32, a33, c3, a35, a36, a41, a42, a43, c4, a45, a46,
755 a51, a52, a53, c5, a55, a56, a61, a62, a63, c6, a65, a66};
756 Double_t data5[]={a11, a12, a13, a14, c1, a16, a21, a22, a23, a24, c2, a26,
757 a31, a32, a33, a34, c3, a36, a41, a42, a43, a44, c4, a46,
758 a51, a52, a53, a54, c5, a56, a61, a62, a63, a64, c6, a66};
759 Double_t data6[]={a11, a12, a13, a14, a15, c1, a21, a22, a23, a24, a25, c2,
760 a31, a32, a33, a34, a35, c3, a41, a42, a43, a44, a45, c4,
761 a51, a52, a53, a54, a55, c5, a61, a62, a63, a64, a65, c6};
763 p1.SetMatrixArray(data1);
764 p2.SetMatrixArray(data2);
765 p3.SetMatrixArray(data3);
766 p4.SetMatrixArray(data4);
767 p5.SetMatrixArray(data5);
768 p6.SetMatrixArray(data6);
769 pC.SetMatrixArray(dataC);
771 // cout << "calculating determinants" << endl;
772 Double_t det0=pC.Determinant();
773 Double_t x1 = p1.Determinant()/det0;
774 Double_t x2 = p2.Determinant()/det0;
775 Double_t x3 = p3.Determinant()/det0;
776 Double_t x4 = p4.Determinant()/det0;
777 Double_t x5 = p5.Determinant()/det0;
778 Double_t x6 = p6.Determinant()/det0;
779 //cout << "calculating determinants done" << endl;
781 AliInfo("p1 singular ");
785 AliInfo("p2 singular ");
789 AliInfo("p3 singular ");
793 AliInfo("p4 singular ");
797 AliInfo("p5 singular ");
801 AliInfo("p6 singular ");