]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TRD/AliTRDalignment.cxx
Replace AliTRDCalibra
[u/mrichter/AliRoot.git] / TRD / AliTRDalignment.cxx
CommitLineData
8775e4e8 1/**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 * *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
6 * *
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 **************************************************************************/
15
16/* $Id$ */
17
18//////////////////////////////////////////////////////////////////////////////////
19// //
20// An AliTRDalignment object contains the alignment data (3 shifts and 3 tilts) //
21// for all the alignable volumes of the TRD, i.e. for 18 supermodules and 540 //
22// chambers. The class provides simple tools for reading and writing these data //
23// in different formats, and for generating fake data that can be used to //
24// simulate misalignment. //
25// The six alignment variables have the following meaning: //
26// shift in rphi //
27// shift in z //
28// shift in r //
29// tilt around rphi //
30// tilt around z //
31// tilt around r //
32// The shifts are in cm and the tilts are in degrees. //
33// The currently supported formats are: //
34// - ascii //
35// - root file containing a TClonesArray of alignment objects //
36// - offline conditions database //
37// - OCDB-like root file //
38// - geometry file (like misaligned_geometry.root) //
39// //
40// D.Miskowiec, November 2006 //
41// //
42//////////////////////////////////////////////////////////////////////////////////
43
44#include <fstream>
45#include <string>
46
47#include "TMath.h"
48#include "TFile.h"
49#include "TGeoManager.h"
50#include "TGeoPhysicalNode.h"
51#include "TClonesArray.h"
52
53#include "AliLog.h"
54#include "AliAlignObj.h"
55#include "AliAlignObjAngles.h"
56#include "AliCDBManager.h"
57#include "AliCDBStorage.h"
58#include "AliCDBMetaData.h"
59#include "AliCDBEntry.h"
60#include "AliCDBId.h"
61
62#include "AliTRDalignment.h"
63
64ClassImp(AliTRDalignment)
65
66//_____________________________________________________________________________
67AliTRDalignment::AliTRDalignment()
68 :TObject()
69 ,fRan(0)
70{
71 //
72 // constructor
73 //
74
75 SetZero();
76
77}
78
79//_____________________________________________________________________________
80AliTRDalignment::AliTRDalignment(const AliTRDalignment& source)
81 :TObject(source)
82 ,fRan(source.fRan)
83{
84 //
85 // copy constructor
86 //
87
88 for (int i = 0; i < 18; i++) {
89 SetSm(i,source.fSm[i]);
90 }
91 for (int i = 0; i < 540; i++) {
92 SetCh(i,source.fCh[i]);
93 }
94
95}
96
97//_____________________________________________________________________________
98AliTRDalignment& AliTRDalignment::operator=(const AliTRDalignment &source)
99{
100 //
101 // assignment operator
102 //
103
104 if (this != &source) {
105 for (int i = 0; i < 18; i++) SetSm(i,source.fSm[i]);
106 for (int i = 0; i < 540; i++) SetCh(i,source.fCh[i]);
107 }
108
109 return *this;
110
111}
112
113//_____________________________________________________________________________
114AliTRDalignment& AliTRDalignment::operator+=(const AliTRDalignment &source)
115{
116 //
117 // addition operator
118 //
119
120 for (int i = 0; i < 18; i++) {
121 for (int j = 0; j < 6; j++) {
122 this->fSm[i][j] =+ source.fSm[i][j];
123 }
124 }
125 for (int i = 0; i < 540; i++) {
126 for (int j = 0; j < 6; j++) {
127 this->fCh[i][j] =+ source.fCh[i][j];
128 }
129 }
130
131 return *this;
132
133}
134
135//_____________________________________________________________________________
136AliTRDalignment& AliTRDalignment::operator-=(const AliTRDalignment &source)
137{
138 //
139 // subtraction operator
140 //
141
142 for (int i = 0; i < 18; i++) {
143 for (int j = 0; j < 6; j++) {
144 fSm[i][j] -= source.fSm[i][j];
145 }
146 }
147 for (int i = 0; i < 540; i++) {
148 for (int j = 0; j < 6; j++) {
149 fCh[i][j] -= source.fCh[i][j];
150 }
151 }
152
153 return *this;
154
155}
156
157//_____________________________________________________________________________
158Bool_t AliTRDalignment::operator==(const AliTRDalignment &source) const
159{
160 //
161 // comparison operator
162 //
163
164 Bool_t areEqual = 1;
165
166 for (int i = 0; i < 18; i++) {
167 for (int j = 0; j < 6; j++) {
168 areEqual &= (fSm[i][j] == source.fSm[i][j]);
169 }
170 }
171 for (int i = 0; i < 540; i++) {
172 for (int j = 0; j < 6; j++) {
173 areEqual &= (fCh[i][j] == source.fCh[i][j]);
174 }
175 }
176
177 return areEqual;
178
179}
180
181//_____________________________________________________________________________
182void AliTRDalignment::SetSmZero()
183{
184 //
185 // reset to zero supermodule data
186 //
187
188 memset(&fSm[0][0],0,sizeof(fSm));
189
190}
191
192//_____________________________________________________________________________
193void AliTRDalignment::SetChZero()
194{
195 //
196 // reset to zero chamber data
197 //
198
199 memset(&fCh[0][0],0,sizeof(fCh));
200
201}
202
203//_____________________________________________________________________________
204void AliTRDalignment::SetSmRandom(Double_t a[6])
205{
206 //
207 // generate random gaussian supermodule data with sigmas a
208 //
209
210 double x[6];
211
212 for (Int_t i = 0; i < 18; i++) {
213 fRan.Rannor(x[0],x[1]);
214 fRan.Rannor(x[2],x[3]);
215 fRan.Rannor(x[4],x[5]);
216 for (Int_t j = 0; j < 6; j++) {
217 x[j] *= a[j];
218 }
219 SetSm(i,x);
220 //PrintSm(i);
221 }
222
223}
224
225//_____________________________________________________________________________
226void AliTRDalignment::SetChRandom(Double_t a[6])
227{
228 //
229 // generate random gaussian chamber data with sigmas a
230 //
231
232 double x[6];
233
234 for (Int_t i = 0; i < 540; i++) {
235 fRan.Rannor(x[0],x[1]);
236 fRan.Rannor(x[2],x[3]);
237 fRan.Rannor(x[4],x[5]);
238 for (Int_t j = 0; j < 6; j++) {
239 x[j] *= a[j];
240 }
241 SetCh(i,x);
242 //PrintCh(i);
243 }
244
245}
246
247//_____________________________________________________________________________
248void AliTRDalignment::SetSmFull()
249{
250 //
251 // generate random gaussian supermodule data similar to the misalignment
252 // expected from the mechanical precision
253 //
254
255 Double_t a[6];
256
257 a[0] = 0.3; // phi
258 a[1] = 0.3; // z
259 a[2] = 0.3; // r
260 a[3] = 0.4/1000.0 / TMath::Pi()*180.0; // phi
261 a[4] = 2.0/1000.0 / TMath::Pi()*180.0; // z
262 a[5] = 0.4/1000.0 / TMath::Pi()*180.0; // r
263
264 SetSmRandom(a);
265
266}
267
268//_____________________________________________________________________________
269void AliTRDalignment::SetChFull()
270{
271 //
272 // generate random gaussian chamber data similar to the misalignment
273 // expected from the mechanical precision
274 //
275
276 Double_t a[6];
277
278 a[0] = 0.1; // phi
279 a[1] = 0.1; // z
280 a[2] = 0.1; // r
281 a[3] = 1.0/1000.0 / TMath::Pi()*180.0; // phi
282 a[4] = 1.0/1000.0 / TMath::Pi()*180.0; // z
283 a[5] = 0.7/1000.0 / TMath::Pi()*180.0; // r
284
285 SetChRandom(a);
286
287}
288
289//_____________________________________________________________________________
290void AliTRDalignment::SetSmResidual()
291{
292 //
293 // generate random gaussian supermodule data similar to the misalignment
294 // remaining after full calibration
295 // I assume that it will be negligible
296 //
297
298 SetSmZero();
299
300}
301
302//_____________________________________________________________________________
303void AliTRDalignment::SetChResidual()
304{
305 //
306 // generate random gaussian chamber data similar to the misalignment
307 // remaining after full calibration
308 //
309
310 Double_t a[6];
311
312 a[0] = 0.002; // phi
313 a[1] = 0.003; // z
314 a[2] = 0.007; // r
315 a[3] = 0.3/1000.0 / TMath::Pi()*180.0; // phi
316 a[4] = 0.3/1000.0 / TMath::Pi()*180.0; // z
317 a[5] = 0.1/1000.0 / TMath::Pi()*180.0; // r
318
319 SetChRandom(a);
320
321}
322
323//_____________________________________________________________________________
324void AliTRDalignment::PrintSm(Int_t i, FILE *fp) const
325{
326 //
327 // print the supermodule data
328 //
329
330 fprintf(fp,"%4d %11.4f %11.4f %11.4f %11.5f %11.5f %11.5f %6d %s\n"
331 ,i,fSm[i][0],fSm[i][1],fSm[i][2],fSm[i][3],fSm[i][4],fSm[i][5]
332 ,0,GetSmName(i));
333
334}
335
336//_____________________________________________________________________________
337void AliTRDalignment::PrintCh(Int_t i, FILE *fp) const
338{
339 //
340 // print the chamber data
341 //
342
343 fprintf(fp,"%4d %11.4f %11.4f %11.4f %11.5f %11.5f %11.5f %6d %s\n"
344 ,i,fCh[i][0],fCh[i][1],fCh[i][2],fCh[i][3],fCh[i][4],fCh[i][5]
345 ,GetVoi(i),GetChName(i));
346
347}
348
349//_____________________________________________________________________________
350void AliTRDalignment::ReadAscii(char *filename)
351{
352 //
353 // read the alignment data from ascii file
354 //
355
356 double x[6]; // alignment data
357 int volid; // volume id
358 std::string syna; // symbolic name
359 int j; // dummy index
360
361 fstream fi(filename,fstream::in);
362 if (!fi) {
363 AliFatal(Form("cannot open input file %s",filename));
364 }
365
366 // supermodules
367
368 for (int i = 0; i < 18; i++) {
369 fi>>j>>x[0]>>x[1]>>x[2]>>x[3]>>x[4]>>x[5]>>volid>>syna;
370 if (j != i) {
371 AliError(Form("sm %d expected, %d found",i,j));
372 }
373 if (volid != 0) {
374 AliError(Form("sm %d volume id %d expected, %d found",i,0,volid));
375 }
376 std::string symnam = GetSmName(i);
377 if (syna != symnam) {
378 AliError(Form("sm %d name %s expected, %s found",i,symnam.data(),syna.data()));
379 }
380 SetSm(i,x);
381 }
382
383 // chambers
384
385 for (int i = 0; i < 540; i++) {
386 fi>>j>>x[0]>>x[1]>>x[2]>>x[3]>>x[4]>>x[5]>>volid>>syna;
387 if (j != i) {
388 AliError(Form("ch %d expected, %d found",i,j));
389 }
390 if (volid != GetVoi(i)) {
391 AliError(Form("ch %d volume id %d expected, %d found",i,GetVoi(i),volid));
392 }
393 std::string symnam = GetChName(i);
394 if (syna != symnam) {
395 AliError(Form("ch %d name %s expected, %s found",i,symnam.data(),syna.data()));
396 }
397 SetCh(i,x);
398 }
399
400 fi.close();
401
402}
403
404//_____________________________________________________________________________
405void AliTRDalignment::ReadRoot(char *filename)
406{
407 //
408 // read the alignment data from root file
409 // here I expect a fixed order and number of elements
410 // it would be much better to identify the alignment objects
411 // one by one and set the parameters of the corresponding sm or ch
412 //
413
414 TFile fi(filename,"READ");
415
416 if (fi.IsOpen()) {
417 TClonesArray *ar = (TClonesArray*) fi.Get("TRDAlignObjs");
418 ArToNumbers(ar);
419 fi.Close();
420 }
421 else {
422 AliError(Form("cannot open input file %s",filename));
423 }
424
425 return;
426
427}
428
429//_____________________________________________________________________________
430void AliTRDalignment::ReadDB(char *filename)
431{
432 //
433 // read the alignment data from database file
434 //
435
436 TFile fi(filename,"READ");
437
438 if (fi.IsOpen()) {
439 AliCDBEntry *e = (AliCDBEntry *) fi.Get("AliCDBEntry");
440 e->PrintMetaData();
441 TClonesArray *ar = (TClonesArray *) e->GetObject();
442 ArToNumbers(ar);
443 fi.Close();
444 }
445 else {
446 AliError(Form("cannot open input file %s",filename));
447 }
448
449 return;
450
451}
452
453//_____________________________________________________________________________
454void AliTRDalignment::ReadDB(char *db, char *path, Int_t run
455 , Int_t version, Int_t subversion)
456{
457 //
458 // read the alignment data from database
459 //
460
461 AliCDBManager *cdb = AliCDBManager::Instance();
462 AliCDBStorage *storLoc = cdb->GetStorage(db);
463 AliCDBEntry *e = storLoc->Get(path,run,version,subversion);
464 e->PrintMetaData();
465 TClonesArray *ar = (TClonesArray *) e->GetObject();
466 ArToNumbers(ar);
467
468}
469
470//_____________________________________________________________________________
471void AliTRDalignment::ReadGeo(char *misaligned)
472{
473 //
474 // determine misalignment by comparing original and misaligned matrix
475 // of the last node on the misaligned_geometry file
476 // an alternative longer way is in attic.C
477 //
478
479 TGeoHMatrix *ideSm[18]; // ideal
480 TGeoHMatrix *ideCh[540];
481 TGeoHMatrix *misSm[18]; // misaligned
482 TGeoHMatrix *misCh[540];
483
484 // read misaligned and original matrices
485
486 TGeoManager::Import(misaligned);
487 for (int i = 0; i < 18; i++) {
488 TGeoPNEntry *pne = gGeoManager->GetAlignableEntry(GetSmName(i));
489 if (!pne) {
490 AliError(Form("no such physical node entry: %s",GetSmName(i)));
491 return;
492 }
493 TGeoPhysicalNode *node = pne->GetPhysicalNode();
494 if (!node) {
495 AliError(Form("physical node entry %s has no physical node",GetSmName(i)));
496 return;
497 }
498 misSm[i] = new TGeoHMatrix(*node->GetNode(node->GetLevel())->GetMatrix());
499 ideSm[i] = new TGeoHMatrix(*node->GetOriginalMatrix());
500 }
501 for (int i = 0; i < 540; i++) {
502 TGeoPNEntry *pne = gGeoManager->GetAlignableEntry(GetChName(i));
503 if (!pne) {
504 AliError(Form("no such physical node entry: %s",GetChName(i)));
505 return;
506 }
507 TGeoPhysicalNode *node = pne->GetPhysicalNode();
508 if (!node) {
509 AliError(Form("physical node entry %s has no physical node",GetChName(i)));
510 return;
511 }
512 misCh[i] = new TGeoHMatrix(*node->GetNode(node->GetLevel())->GetMatrix());
513 ideCh[i] = new TGeoHMatrix(*node->GetOriginalMatrix());
514 }
515
516 // calculate the local misalignment matrices as inverse misaligned times ideal
517
518 for (int i = 0; i < 18; i++) {
519 TGeoHMatrix mat(ideSm[i]->Inverse());
520 mat.Multiply(misSm[i]);
521 double *tra = mat.GetTranslation();
522 double *rot = mat.GetRotationMatrix();
523 double pars[6];
524 pars[0] = tra[0];
525 pars[1] = tra[1];
526 pars[2] = tra[2];
527 if (TMath::Abs(rot[0])<1e-7 || TMath::Abs(rot[8])<1e-7) {
528 AliError("Failed to extract roll-pitch-yall angles!");
529 return;
530 }
531 double raddeg = TMath::RadToDeg();
532 pars[3] = raddeg * TMath::ATan2(-rot[5],rot[8]);
533 pars[4] = raddeg * TMath::ASin(rot[2]);
534 pars[5] = raddeg * TMath::ATan2(-rot[1],rot[0]);
535 SetSm(i,pars);
536 }
537
538 for (int i = 0; i < 540; i++) {
539 TGeoHMatrix mat(ideCh[i]->Inverse());
540 mat.Multiply(misCh[i]);
541 double *tra = mat.GetTranslation();
542 double *rot = mat.GetRotationMatrix();
543 double pars[6];
544 pars[0] = tra[0];
545 pars[1] = tra[1];
546 pars[2] = tra[2];
547 if(TMath::Abs(rot[0])<1e-7 || TMath::Abs(rot[8])<1e-7) {
548 AliError("Failed to extract roll-pitch-yall angles!");
549 return;
550 }
551 double raddeg = TMath::RadToDeg();
552 pars[3] = raddeg * TMath::ATan2(-rot[5],rot[8]);
553 pars[4] = raddeg * TMath::ASin(rot[2]);
554 pars[5] = raddeg * TMath::ATan2(-rot[1],rot[0]);
555 SetCh(i,pars);
556 }
557
558 // cleanup
559 for (int i = 0; i < 18; i++) delete ideSm[i];
560 for (int i = 0; i < 18; i++) delete misSm[i];
561 for (int i = 0; i < 540; i++) delete ideCh[i];
562 for (int i = 0; i < 540; i++) delete misCh[i];
563
564 return;
565
566}
567
568//_____________________________________________________________________________
569void AliTRDalignment::ReadSurveyReport(char *filename)
570{
571 // read survey report and set the supermodule parameters correspondingly
572
573 fstream fi(filename,fstream::in);
574 if (!fi) {
575 AliFatal(Form("cannot open input file %s",filename));
576 }
577
578 // to be continued...
579
580}
581
582//_____________________________________________________________________________
583void AliTRDalignment::ReadAny(char *filename)
584{
585 //
586 // read the alignment data from any kind of file
587 //
588
589 TString fist(filename);
590 if (fist.EndsWith(".txt")) {
591 ReadAscii(filename);
592 }
593 if (fist.EndsWith(".dat")) {
594 ReadAscii(filename);
595 }
596 if (fist.EndsWith(".root")) {
597 if (fist.Contains("Run")) {
598 ReadDB(filename);
599 }
600 else {
601 ReadRoot(filename);
602 }
603 }
604
605}
606
607//_____________________________________________________________________________
608void AliTRDalignment::WriteAscii(char *filename) const
609{
610 //
611 // store the alignment data on ascii file
612 //
613
614 FILE *fp = fopen(filename, "w");
615 if (!fp) {
616 AliError(Form("cannot open output file %s",filename));
617 return;
618 }
619
620 PrintSm(fp);
621 PrintCh(fp);
622
623 fclose(fp);
624
625}
626
627//_____________________________________________________________________________
628void AliTRDalignment::WriteRoot(char *filename)
629{
630 //
631 // store the alignment data on root file
632 //
633
634 TClonesArray *ar = new TClonesArray("AliAlignObjAngles",10000);
635 NumbersToAr(ar);
636 TFile fo(filename,"RECREATE");
637 if (fo.IsOpen()) {
638 fo.cd();
639 fo.WriteObject(ar,"TRDAlignObjs","kSingleKey");
640 fo.Close();
641 }
642 else {
643 AliError(Form("cannot open output file %s",filename));
644 }
645
646 delete ar;
647
648}
649
650//_____________________________________________________________________________
651void AliTRDalignment::WriteDB(char *filename, char *comment, Int_t run0, Int_t run1)
652{
653 //
654 // dumping on a DB-like file
655 //
656
657 TClonesArray *ar = new TClonesArray("AliAlignObjAngles",10000);
658 NumbersToAr(ar);
659 char *path = "di1/di2/di3";
660 AliCDBId id(path,run0,run1);
661 AliCDBMetaData *md = new AliCDBMetaData();
662 md->SetResponsible("Dariusz Miskowiec");
663 md->SetComment(comment);
664 AliCDBEntry *e = new AliCDBEntry(ar, id, md);
665 TFile fi(filename,"RECREATE");
666 if (fi.IsOpen()) {
667 e->Write();
668 fi.Close();
669 }
670 else {
671 AliError(Form("cannot open input file %s",filename));
672 }
673
674 delete e;
675 delete md;
676 delete ar;
677
678 return;
679
680}
681
682//_____________________________________________________________________________
683void AliTRDalignment::WriteDB(char *db, char *path, char *comment, Int_t run0, Int_t run1)
684{
685 //
686 // store the alignment data in database
687 //
688
689 TClonesArray *ar = new TClonesArray("AliAlignObjAngles",10000);
690 NumbersToAr(ar);
691 AliCDBManager *cdb = AliCDBManager::Instance();
692 AliCDBStorage *storLoc = cdb->GetStorage(db);
693 AliCDBMetaData *md = new AliCDBMetaData();
694 md->SetResponsible("Dariusz Miskowiec");
695 md->SetComment(comment);
696 AliCDBId id(path,run0,run1);
697 storLoc->Put(ar,id,md);
698 md->Delete();
699 delete ar;
700
701}
702
703//_____________________________________________________________________________
704void AliTRDalignment::WriteGeo(char *filename)
705{
706 //
707 // apply misalignment to (currently loaded ideal) geometry and store the
708 // resulting geometry on a root file
709 //
710
711 TClonesArray *ar = new TClonesArray("AliAlignObjAngles",10000);
712 NumbersToAr(ar);
713 for (int i = 0; i < ar->GetEntriesFast(); i++) {
714 AliAlignObj *alobj = (AliAlignObj *) ar->UncheckedAt(i);
715 alobj->ApplyToGeometry();
716 }
717 delete ar;
718 gGeoManager->Export(filename);
719
720}
721
722//_____________________________________________________________________________
723Double_t AliTRDalignment::GetSmRMS(Int_t xyz) const
724{
725 //
726 // rms fSm[][xyz]
727 //
728
729 Double_t s1 = 0.0;
730 Double_t s2 = 0.0;
731 for (int i = 0; i < 18; i++) {
732 s1 += fSm[i][xyz];
733 s2 += fSm[i][xyz]*fSm[i][xyz];
734 }
735 Double_t rms2 = s2/18.0 - s1*s1/18.0/18.0;
736
737 return rms2>0 ? sqrt(rms2) : 0.0;
738
739}
740
741//_____________________________________________________________________________
742Double_t AliTRDalignment::GetChRMS(Int_t xyz) const
743{
744 //
745 // rms fCh[][xyz]
746 //
747
748 Double_t s1 =0.0;
749 Double_t s2 =0.0;
750 for (int i = 0; i < 540; i++) {
751 s1 += fCh[i][xyz];
752 s2 += fCh[i][xyz]*fCh[i][xyz];
753 }
754 Double_t rms2 = s2/540.0 - s1*s1/540.0/540.0;
755
756 return rms2>0 ? sqrt(rms2) : 0.0;
757
758}
759
760//_____________________________________________________________________________
761void AliTRDalignment::PrintSmRMS() const
762{
763 //
764 // dump rms of fSm
765 //
766
767 printf(" %11.4f %11.4f %11.4f %11.5f %11.5f %11.5f supermodule rms\n"
768 ,GetSmRMS(0),GetSmRMS(1),GetSmRMS(2),GetSmRMS(3),GetSmRMS(4),GetSmRMS(5));
769
770}
771
772//_____________________________________________________________________________
773void AliTRDalignment::PrintChRMS() const
774{
775 //
776 // dump rms of fCh
777 //
778
779 printf(" %11.4f %11.4f %11.4f %11.5f %11.5f %11.5f chamber rms\n"
780 ,GetChRMS(0),GetChRMS(1),GetChRMS(2),GetChRMS(3),GetChRMS(4),GetChRMS(5));
781
782}
783
784//_____________________________________________________________________________
785void AliTRDalignment::ArToNumbers(TClonesArray *ar)
786{
787 //
788 // read numbers from the array of AliAlignObj objects and fill fSm and fCh
789 //
790
791 LoadIdealGeometry();
792 SetZero();
793 double pa[6];
794 for (int i = 0; i < 18; i++) {
795 AliAlignObj *aao = (AliAlignObj *) ar->At(i);
796 aao->GetLocalPars(pa,pa+3);
797 SetSm(i,pa);
798 }
799 for (int i = 0; i < 540; i++) {
800 AliAlignObj *aao = (AliAlignObj *) ar->At(18+i);
801 aao->GetLocalPars(pa,pa+3);
802 SetCh(i,pa);
803 }
804
805}
806
807//_____________________________________________________________________________
808void AliTRDalignment::NumbersToAr(TClonesArray *ar)
809{
810 //
811 // build array of AliAlignObj objects based on fSm and fCh data
812 //
813
814 LoadIdealGeometry();
815 TClonesArray &alobj = *ar;
816 int nobj = 0;
817 for (int i = 0; i < 18; i++) {
818 new(alobj[nobj]) AliAlignObjAngles(GetSmName(i)
819 ,0
820 ,fSm[i][0],fSm[i][1],fSm[i][2]
821 ,fSm[i][3],fSm[i][4],fSm[i][5]
822 ,0);
823 nobj++;
824 }
825
826 for (int i = 0; i < 540; i++) {
827 new(alobj[nobj]) AliAlignObjAngles(GetChName(i)
828 ,GetVoi(i)
829 ,fCh[i][0],fCh[i][1],fCh[i][2]
830 ,fCh[i][3],fCh[i][4],fCh[i][5]
831 ,0);
832 nobj++;
833 }
834
835}
836
837//_____________________________________________________________________________
838void AliTRDalignment::LoadIdealGeometry(char *filename)
839{
840 //
841 // load ideal geometry from filename
842 // it is needed for operations on AliAlignObj objects
843 // this needs to be straightened out
844 // particularly, sequences LoadIdealGeometry("file1"); LoadIdealGeometry("file2");
845 // do not work as one would naturally expect
846 //
847
848 static int attempt = 0; // which reload attempt is it? just to avoid endless loops
849
850 if (!gGeoManager) {
851 TGeoManager::Import(filename);
852 }
853 if (!gGeoManager) {
854 AliFatal(Form("cannot open geometry file %s",filename));
855 }
856 if (gGeoManager->GetListOfPhysicalNodes()->GetEntries()) {
857 if (attempt) {
858 AliFatal(Form("geometry on file %s is not ideal",filename));
859 }
860 AliWarning("current geometry is not ideal - it contains physical nodes");
861 AliWarning(Form("reloading geometry from %s - attempt nr %d",filename,attempt));
862 gGeoManager = 0;
863 attempt++;
864 LoadIdealGeometry(filename);
865 }
866
867 attempt = 0;
868
869}