]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/AliTRDalignment.cxx
Correct character constantness
[u/mrichter/AliRoot.git] / TRD / AliTRDalignment.cxx
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 // An AliTRDalignment object contains the alignment data (3 shifts and 3     //
20 // tilts) for all the alignable volumes of the TRD, i.e. for 18 supermodules //
21 // and 540 chambers. The class provides simple tools for reading and writing //
22 // these data in different formats, and for generating fake data that can be //
23 // used to simulate misalignment.                                            //
24 // The six alignment variables have the following meaning:                   //
25 // shift in rphi                                                             //
26 // shift in z                                                                //
27 // shift in r                                                                //
28 // tilt around rphi                                                          //
29 // tilt around z                                                             //
30 // tilt around r                                                             //
31 // The shifts are in cm and the tilts are in degrees.                        //
32 // The currently supported formats are:                                      //
33 // - ascii                                                                   //
34 // - root file containing a TClonesArray of alignment objects                //
35 // - offline conditions database                                             //
36 // - OCDB-like root file                                                     //
37 // - geometry file (like misaligned_geometry.root)                           //
38 //                                                                           //
39 // Some examples of usage (in an aliroot session):                           //
40 // AliTRDalignment a,b,c,d,e;                                                //
41 // double xsm[]={0,0,0,-70,0,0};                                             //
42 // double xch[]={0,0,-50,0,0,0};                                             //
43 // a.SetSm(4,xsm);                                                           // 
44 // a.SetCh(120,xch);                                                         //
45 // a.WriteAscii("kuku.dat");                                                 //
46 // TGeoManager::Import("geometry.root"); a.WriteRoot("kuku.root");           //
47 // TGeoManager::Import("geometry.root"); a.WriteDB("kukudb.root",0,0);       //
48 // TGeoManager::Import("geometry.root");                                     //
49 // a.WriteDB("local://$ALICE_ROOT", "TRD/Align/Data", 0,0);                  //
50 // TGeoManager::Import("geometry.root"); a.WriteGeo("kukugeometry.root");    //
51 //                                                                           //
52 // b.ReadAscii("kuku.dat");                                                  //
53 // TGeoManager::Import("geometry.root"); c.ReadRoot("kuku.root");            //
54 // TGeoManager::Import("geometry.root"); d.ReadDB("kukudb.root");            //
55 // TGeoManager::Import("kukugeometry.root"); e.ReadCurrentGeo();             //
56 //                                                                           //
57 // e.PrintSm(4);                                                             //
58 // e.PrintCh(120);                                                           // 
59 // a.PrintRMS();                                                             //
60 // b.PrintRMS();                                                             //
61 // e.PrintRMS();                                                             //
62 //                                                                           //
63 //                                                                           //
64 // D.Miskowiec, November 2006                                                //
65 //                                                                           //
66 ///////////////////////////////////////////////////////////////////////////////
67
68 #include <iostream>
69 #include <fstream>
70 //#include <string>
71
72 #include "TMath.h"
73 #include "TFile.h"
74 #include "TGeoManager.h"
75 #include "TGeoPhysicalNode.h"
76 #include "TClonesArray.h"
77 #include "TString.h"
78 #include "TFitter.h"
79 #include "TMinuit.h"
80
81 #include "AliLog.h"
82 #include "AliAlignObj.h"
83 #include "AliAlignObjParams.h"
84 #include "AliCDBManager.h"
85 #include "AliCDBStorage.h"
86 #include "AliCDBMetaData.h"
87 #include "AliCDBEntry.h"
88 #include "AliCDBId.h"
89 #include "AliSurveyObj.h"
90 #include "AliSurveyPoint.h"
91
92 #include "AliTRDalignment.h"
93
94 void trdAlignmentFcn(Int_t &npar, Double_t *gin, Double_t &f, Double_t *x, Int_t iflag);
95
96 ClassImp(AliTRDalignment)
97
98 //_____________________________________________________________________________
99 AliTRDalignment::AliTRDalignment() 
100   :TObject()
101   ,fComment()
102   ,fRan(0)
103 {
104   //
105   // constructor
106   //
107
108   SetZero();
109
110   for (int i=0; i<18; i++) for (int j=0; j<2; j++) for (int k=0; k<2; k++) for (int l=0; l<2; l++) {
111     fSurveyX[i][j][k][l] = 0.0;
112     fSurveyY[i][j][k][l] = 0.0;
113     fSurveyZ[i][j][k][l] = 0.0;
114     fSurveyEX[i][j][k][l] = 0.0;
115     fSurveyEY[i][j][k][l] = 0.0;
116     fSurveyEZ[i][j][k][l] = 0.0;
117   }
118
119   // Initialize the nominal positions of the survey points 
120   // in the local frame of supermodule (where y is the long side, 
121   // z corresponds to the radius in lab, and x to the phi in lab).
122   // Four survey marks are on each z-side of the supermodule. 
123   //               A           B
124   //           ----o-----------o----        x |
125   //           \                   /          |
126   //            \                 /           |
127   //             \               /            |
128   //              \             /             |
129   //               ---o-----o---              -------------->
130   //                  C     D                              y
131   // 
132   // For the purpose of this explanation lets define the origin such that 
133   // the supermodule occupies 0 < x < 77.9 cm. Then the coordinates (x,y) 
134   // are (in cm) 
135   // A (76.2,-30.25)
136   // B (76.2,+30.25)
137   // C ( 2.2,-22.5 )
138   // D ( 2.2,+22.5 )
139   // 
140
141   double x[2] = {22.5,30.25};                   // lab phi, or tracking-y
142   double y[2] = {353.0, -353.0};                // lab z; inc. 2 cm survey target offset
143   double z[2] = {-(77.9/2.0-2.0),77.9/2.0-1.5}; // lab r, or better tracking-x
144
145   for (int j=0; j<2; j++) for (int k=0; k<2; k++) for (int l=0; l<2; l++) {
146     fSurveyX0[j][k][l] = -TMath::Power(-1,l) * x[k];
147     fSurveyY0[j][k][l] = y[j];
148     fSurveyZ0[j][k][l] = z[k];
149   }
150
151 }
152
153 //_____________________________________________________________________________
154 AliTRDalignment::AliTRDalignment(const AliTRDalignment& source) 
155   :TObject(source)
156   ,fComment(source.fComment)
157   ,fRan(source.fRan)
158 {
159   //
160   // copy constructor
161   //
162
163   for (int i=0; i<18; i++) SetSm(i,source.fSm[i]);
164   for (int i=0; i<540; i++) SetCh(i,source.fCh[i]);
165   for (int i=0; i<18; i++) for (int j=0; j<2; j++) for (int k=0; k<2; k++) for (int l=0; l<2; l++) {
166     fSurveyX[i][j][k][l] = source.fSurveyX[i][j][k][l];
167     fSurveyY[i][j][k][l] = source.fSurveyY[i][j][k][l];
168     fSurveyZ[i][j][k][l] = source.fSurveyZ[i][j][k][l];
169     fSurveyEX[i][j][k][l] = source.fSurveyEX[i][j][k][l];
170     fSurveyEY[i][j][k][l] = source.fSurveyEY[i][j][k][l];
171     fSurveyEZ[i][j][k][l] = source.fSurveyEZ[i][j][k][l];
172   }
173   for (int j=0; j<2; j++) for (int k=0; k<2; k++) for (int l=0; l<2; l++) {
174     fSurveyX0[j][k][l] = source.fSurveyX0[j][k][l];
175     fSurveyY0[j][k][l] = source.fSurveyY0[j][k][l];
176     fSurveyZ0[j][k][l] = source.fSurveyZ0[j][k][l];
177   }
178
179 }
180
181 //_____________________________________________________________________________
182 AliTRDalignment& AliTRDalignment::operator=(const AliTRDalignment &source) 
183 {
184   //
185   // assignment operator
186   //
187
188   if (this != &source) {
189     for (int i = 0; i <  18; i++) SetSm(i,source.fSm[i]);
190     for (int i = 0; i < 540; i++) SetCh(i,source.fCh[i]);
191     for (int i=0; i<18; i++) for (int j=0; j<2; j++) for (int k=0; k<2; k++) for (int l=0; l<2; l++) {
192       fSurveyX[i][j][k][l] = source.fSurveyX[i][j][k][l];
193       fSurveyY[i][j][k][l] = source.fSurveyY[i][j][k][l];
194       fSurveyZ[i][j][k][l] = source.fSurveyZ[i][j][k][l];
195       fSurveyEX[i][j][k][l] = source.fSurveyEX[i][j][k][l];
196       fSurveyEY[i][j][k][l] = source.fSurveyEY[i][j][k][l];
197       fSurveyEZ[i][j][k][l] = source.fSurveyEZ[i][j][k][l];
198     }
199     for (int j=0; j<2; j++) for (int k=0; k<2; k++) for (int l=0; l<2; l++) {
200       fSurveyX0[j][k][l] = source.fSurveyX0[j][k][l];
201       fSurveyY0[j][k][l] = source.fSurveyY0[j][k][l];
202       fSurveyZ0[j][k][l] = source.fSurveyZ0[j][k][l];
203     }
204     fComment = source.fComment;
205   }
206
207   return *this;
208
209 }
210
211 //_____________________________________________________________________________
212 AliTRDalignment& AliTRDalignment::operator*=(double fac) 
213 {
214   //
215   // multiplication operator
216   //
217
218   for (int i = 0; i <  18; i++) for (int j = 0; j < 6; j++) this->fSm[i][j] *= fac;
219   for (int i = 0; i < 540; i++) for (int j = 0; j < 6; j++) this->fCh[i][j] *= fac;
220
221   return *this;
222
223 }
224
225 //_____________________________________________________________________________
226 AliTRDalignment& AliTRDalignment::operator+=(const AliTRDalignment &source) 
227 {
228   //
229   // addition operator
230   //
231
232   for (int i = 0; i <  18; i++) for (int j = 0; j < 6; j++) this->fSm[i][j] += source.fSm[i][j];
233   for (int i = 0; i < 540; i++) for (int j = 0; j < 6; j++) this->fCh[i][j] += source.fCh[i][j];
234
235   return *this;
236
237 }
238
239 //_____________________________________________________________________________
240 AliTRDalignment& AliTRDalignment::operator-=(const AliTRDalignment &source) 
241 {
242   //
243   // subtraction operator
244   //
245
246   for (int i = 0; i <  18; i++) for (int j = 0; j < 6; j++) fSm[i][j] -= source.fSm[i][j];
247   for (int i = 0; i < 540; i++) for (int j = 0; j < 6; j++) fCh[i][j] -= source.fCh[i][j];
248
249   return *this;
250
251 }
252
253 //_____________________________________________________________________________
254 Bool_t AliTRDalignment::operator==(const AliTRDalignment &source) const
255 {
256   //
257   // comparison operator
258   //
259
260   Bool_t areEqual = 1;
261
262   for (int i = 0; i <  18; i++) for (int j = 0; j < 6; j++) areEqual &= (fSm[i][j] == source.fSm[i][j]);
263   for (int i = 0; i < 540; i++) for (int j = 0; j < 6; j++) areEqual &= (fCh[i][j] == source.fCh[i][j]);
264
265   return areEqual;
266
267 }
268
269 //_____________________________________________________________________________
270 void AliTRDalignment::SetSmZero() 
271 {
272   //
273   // reset to zero supermodule data
274   //
275
276   memset(&fSm[0][0],0,sizeof(fSm));
277
278 }
279
280 //_____________________________________________________________________________
281 void AliTRDalignment::SetChZero() 
282 {
283   //
284   // reset to zero chamber data
285   //
286
287   memset(&fCh[0][0],0,sizeof(fCh));
288
289 }
290
291 //_____________________________________________________________________________
292 void AliTRDalignment::SetSmRandom(double a[6]) 
293 {
294   //
295   // generate random gaussian supermodule data with sigmas a
296   //
297
298   double x[6];
299
300   for (int i = 0; i < 18; i++) {
301     fRan.Rannor(x[0],x[1]);
302     fRan.Rannor(x[2],x[3]);
303     fRan.Rannor(x[4],x[5]);
304     for (int j = 0; j < 6; j++) x[j] *= a[j];
305     SetSm(i,x);
306     //PrintSm(i);
307   }
308
309 }
310
311 //_____________________________________________________________________________
312 void AliTRDalignment::SetChRandom(double a[6]) 
313 {
314   //
315   // generate random gaussian chamber data with sigmas a
316   //
317
318   double x[6];
319
320   for (int i = 0; i < 540; i++) {
321     fRan.Rannor(x[0],x[1]);
322     fRan.Rannor(x[2],x[3]);
323     fRan.Rannor(x[4],x[5]);
324     for (int j = 0; j < 6; j++) x[j] *= a[j];
325     SetCh(i,x);
326     //PrintCh(i);
327   }
328
329 }
330
331 //_____________________________________________________________________________
332 void AliTRDalignment::SetSmFull() 
333 {
334   //
335   // generate random gaussian supermodule data similar to the misalignment 
336   // expected from the mechanical precision 
337   //
338
339   double a[6];
340
341   a[0] = 0.3; // phi
342   a[1] = 0.3; // z
343   a[2] = 0.3; // r
344   a[3] = 0.4/1000.0 / TMath::Pi()*180.0; // phi
345   a[4] = 2.0/1000.0 / TMath::Pi()*180.0; // z
346   a[5] = 0.4/1000.0 / TMath::Pi()*180.0; // r
347
348   SetSmRandom(a);
349
350 }
351
352 //_____________________________________________________________________________
353 void AliTRDalignment::SetChFull() 
354 {
355   //
356   // generate random gaussian chamber data similar to the misalignment 
357   // expected from the mechanical precision 
358   //
359
360   double a[6];
361
362   a[0] = 0.1; // phi
363   a[1] = 0.1; // z
364   a[2] = 0.1; // r
365   a[3] = 1.0/1000.0 / TMath::Pi()*180.0; // phi
366   a[4] = 1.0/1000.0 / TMath::Pi()*180.0; // z
367   a[5] = 0.7/1000.0 / TMath::Pi()*180.0; // r
368
369   SetChRandom(a);
370
371 }
372
373 //_____________________________________________________________________________
374 void AliTRDalignment::SetSmResidual() 
375 {
376   //
377   // generate random gaussian supermodule data similar to the misalignment 
378   // remaining after full calibration
379   // I assume that it will be negligible
380   //
381
382   SetSmZero();
383
384 }
385
386 //_____________________________________________________________________________
387 void AliTRDalignment::SetChResidual() 
388 {
389   //
390   // generate random gaussian chamber data similar to the misalignment 
391   // remaining after full calibration
392   //
393
394   double a[6];
395
396   a[0] = 0.002; // phi
397   a[1] = 0.003; // z
398   a[2] = 0.007; // r
399   a[3] = 0.3/1000.0 / TMath::Pi()*180.0; // phi
400   a[4] = 0.3/1000.0 / TMath::Pi()*180.0; // z
401   a[5] = 0.1/1000.0 / TMath::Pi()*180.0; // r
402
403   SetChRandom(a);
404
405 }
406
407 //_____________________________________________________________________________
408 void AliTRDalignment::PrintSm(int i, FILE *fp) const 
409 {
410   //
411   // print the supermodule data
412   //
413
414   fprintf(fp,"%4d   %11.4f %11.4f  %11.4f      %11.5f  %11.5f  %11.5f   %6d  %s\n"
415          ,i,fSm[i][0],fSm[i][1],fSm[i][2],fSm[i][3],fSm[i][4],fSm[i][5]
416          ,0,GetSmName(i));
417
418 }
419
420 //_____________________________________________________________________________
421 void AliTRDalignment::PrintCh(int i, FILE *fp) const 
422 {
423   //
424   // print the chamber data
425   //
426
427   fprintf(fp,"%4d   %11.4f %11.4f  %11.4f      %11.5f  %11.5f  %11.5f   %6d  %s\n"
428          ,i,fCh[i][0],fCh[i][1],fCh[i][2],fCh[i][3],fCh[i][4],fCh[i][5]
429          ,GetVoi(i),GetChName(i));
430
431 }
432
433 //_____________________________________________________________________________
434 void AliTRDalignment::ReadAscii(char *filename) 
435 {
436   //
437   // read the alignment data from ascii file
438   //
439
440   double x[6];      // alignment data
441   int volid;        // volume id
442   std::string syna; // symbolic name
443   int j;            // dummy index
444
445   fstream fi(filename,fstream::in);
446   if (!fi) {
447     AliError(Form("cannot open input file %s",filename));
448     return;
449   }
450
451   // supermodules
452
453   for (int i = 0; i < 18; i++) {
454     fi>>j>>x[0]>>x[1]>>x[2]>>x[3]>>x[4]>>x[5]>>volid>>syna;
455     if (j != i) AliError(Form("sm %d expected, %d found",i,j));
456     if (volid != 0) AliError(Form("sm %d volume id %d expected, %d found",i,0,volid));
457     std::string symnam = GetSmName(i);
458     if (syna != symnam) AliError(Form("sm %d name %s expected, %s found",i,symnam.data(),syna.data()));
459     SetSm(i,x);
460   }
461
462   // chambers
463
464   for (int i = 0; i < 540; i++) {
465     fi>>j>>x[0]>>x[1]>>x[2]>>x[3]>>x[4]>>x[5]>>volid>>syna;
466     if (j != i) AliError(Form("ch %d expected, %d found",i,j));
467     if (volid != GetVoi(i)) AliError(Form("ch %d volume id %d expected, %d found",i,GetVoi(i),volid));
468     std::string symnam = GetChName(i);
469     if (syna != symnam) AliError(Form("ch %d name %s expected, %s found",i,symnam.data(),syna.data()));
470     SetCh(i,x);
471   }
472
473   fi.close();
474
475 }
476
477 //_____________________________________________________________________________
478 void AliTRDalignment::ReadCurrentGeo() 
479 {
480   //
481   // use currently loaded geometry to determine misalignment by comparing 
482   // original and misaligned matrix of the last node
483   // Now, original, does not mean "ideal". It is the matrix before the alignment. 
484   // So, if alignment was applied more than once, the numbers extracted will 
485   // represent just the last alignment. -- check this!
486   //
487
488   TGeoPNEntry *pne;
489   TGeoHMatrix *ideSm[18];  // ideal
490   TGeoHMatrix *misSm[18];  // misaligned
491   for (int i = 0; i < 18; i++) if ((pne = gGeoManager->GetAlignableEntry(GetSmName(i)))) {
492
493     // read misaligned and original matrices
494
495     TGeoPhysicalNode *node = pne->GetPhysicalNode();
496     if (!node) AliError(Form("physical node entry %s has no physical node",GetSmName(i)));
497     if (!node) continue;
498     misSm[i] = new TGeoHMatrix(*node->GetNode(node->GetLevel())->GetMatrix());
499     ideSm[i] = new TGeoHMatrix(*node->GetOriginalMatrix());
500
501     // calculate the local misalignment matrices as inverse misaligned times ideal
502
503     TGeoHMatrix mat(ideSm[i]->Inverse()); 
504     mat.Multiply(misSm[i]);
505     double *tra = mat.GetTranslation();
506     double *rot = mat.GetRotationMatrix();
507     double pars[6];
508     pars[0] = tra[0];
509     pars[1] = tra[1];
510     pars[2] = tra[2];
511     if (TMath::Abs(rot[0])<1e-7 || TMath::Abs(rot[8])<1e-7) AliError("Failed to extract roll-pitch-yall angles!");
512     double raddeg = TMath::RadToDeg();
513     pars[3] = raddeg * TMath::ATan2(-rot[5],rot[8]);
514     pars[4] = raddeg * TMath::ASin(rot[2]);
515     pars[5] = raddeg * TMath::ATan2(-rot[1],rot[0]);
516     SetSm(i,pars);
517
518     // cleanup
519
520     delete ideSm[i];
521     delete misSm[i];
522   }
523
524   TGeoHMatrix *ideCh[540]; // ideal
525   TGeoHMatrix *misCh[540]; // misaligned
526   for (int i = 0; i < 540; i++) if ((pne = gGeoManager->GetAlignableEntry(GetChName(i)))) {
527
528     // read misaligned and original matrices
529
530     TGeoPhysicalNode *node = pne->GetPhysicalNode();
531     if (!node) AliError(Form("physical node entry %s has no physical node",GetChName(i)));
532     if (!node) continue;
533     misCh[i] = new TGeoHMatrix(*node->GetNode(node->GetLevel())->GetMatrix());
534     ideCh[i] = new TGeoHMatrix(*node->GetOriginalMatrix());
535
536     // calculate the local misalignment matrices as inverse misaligned times ideal
537
538     TGeoHMatrix mat(ideCh[i]->Inverse()); 
539     mat.Multiply(misCh[i]);
540     double *tra = mat.GetTranslation();
541     double *rot = mat.GetRotationMatrix();
542     double pars[6];
543     pars[0] = tra[0];
544     pars[1] = tra[1];
545     pars[2] = tra[2];
546     if(TMath::Abs(rot[0])<1e-7 || TMath::Abs(rot[8])<1e-7) {
547       AliError("Failed to extract roll-pitch-yall angles!");
548       return;
549     }
550     double raddeg = TMath::RadToDeg();
551     pars[3] = raddeg * TMath::ATan2(-rot[5],rot[8]);
552     pars[4] = raddeg * TMath::ASin(rot[2]);
553     pars[5] = raddeg * TMath::ATan2(-rot[1],rot[0]);
554     SetCh(i,pars);
555
556     // cleanup
557     delete ideCh[i];
558     delete misCh[i];
559   }
560
561   return;
562
563 }
564
565 //_____________________________________________________________________________
566 void AliTRDalignment::ReadRoot(char *filename) 
567 {
568   //
569   // read the alignment data from root file
570   //
571
572   TFile fi(filename,"READ");
573
574   if (fi.IsOpen()) {
575     TClonesArray *ar = (TClonesArray*) fi.Get("TRDAlignObjs");
576     ArToNumbers(ar);
577     fi.Close();
578   } 
579   else AliError(Form("cannot open input file %s",filename));
580
581   return;
582
583 }
584
585 //_____________________________________________________________________________
586 void AliTRDalignment::ReadDB(char *filename) 
587 {
588   //
589   // read the alignment data from database file
590   //
591
592   TFile fi(filename,"READ");
593
594   if (fi.IsOpen()) {
595     AliCDBEntry  *e  = (AliCDBEntry *) fi.Get("AliCDBEntry");
596     e->PrintMetaData();
597     fComment.SetString(e->GetMetaData()->GetComment());
598     TClonesArray *ar = (TClonesArray *) e->GetObject();
599     ArToNumbers(ar);
600     fi.Close();
601   } 
602   else AliError(Form("cannot open input file %s",filename));
603
604   return;
605
606 }
607
608 //_____________________________________________________________________________
609 void AliTRDalignment::ReadDB(char *db, char *path, int run
610                            , int version, int subversion)
611 {
612   //
613   // read the alignment data from database
614   //
615
616   AliCDBManager *cdb     = AliCDBManager::Instance();
617   AliCDBStorage *storLoc = cdb->GetStorage(db);
618   AliCDBEntry   *e       = storLoc->Get(path,run,version,subversion);
619   if (e) {
620     e->PrintMetaData();
621     fComment.SetString(e->GetMetaData()->GetComment());
622     TClonesArray  *ar      =  (TClonesArray *) e->GetObject();
623     ArToNumbers(ar);
624   }
625 }
626
627 //_____________________________________________________________________________
628 Bool_t AliTRDalignment::DecodeSurveyPointName(TString pna, Int_t &sm, Int_t &iz, 
629                                               Int_t &ir, Int_t &iphi) {
630   // decode the survey point name and extract the sm, z, r and phi indices
631   
632   if (pna(0,6)!="TRD_sm") {
633     AliError(Form("unexpected point name: %s",pna.Data()));
634     return kFALSE;
635   }
636   sm = atoi(pna(6,2).Data()); // supermodule number
637   iz = -1;
638   if (pna(8) == 'a') iz=0; // anticlockwise, positive z
639   if (pna(8) == 'c') iz=1; // clockwise, negative z
640   ir = -1;
641   if (pna(9) == 'l') ir=0; // low radius
642   if (pna(9) == 'h') ir=1; // high radius
643   iphi = atoi(pna(10,0).Data()); // phi within supermodule
644   if (sm>=0 && sm<18 && iz>=0 && iz<2 && ir>=0 && ir<2 && iphi>=0 && iphi<2) return kTRUE;
645   AliError(Form("cannot decode point name: %s",pna.Data()));
646   return kFALSE;
647 }
648 //_____________________________________________________________________________
649 void AliTRDalignment::ReadSurveyReport(char *filename) 
650 {
651   //
652   // Read survey report and store the numbers in fSurveyX, fSurveyY, fSurveyZ, 
653   // and fSurveyE.  Store the survey info in the fComment.
654   // Each supermodule has 8 survey points. The point names look like 
655   // TRD_sm08ah0 and have the following meaning. 
656   //
657   // sm00..17 mean supermodule 0 through 17, following the phi.
658   // Supermodule 00 is between phi=0 and phi=20 degrees.
659   //
660   // a or c denotes the anticlockwise and clockwise end of the supermodule
661   // in z. Clockwise end is where z is negative and where the muon arm sits.
662   //
663   // l or h denote low radius and high radius holes
664   //
665   // 0 or 1 denote the hole at smaller and at larger phi, respectively.
666   //
667
668   // read the survey file
669
670   fstream in(filename,fstream::in);
671   if (!in) {
672     AliError(Form("cannot open input file %s",filename));
673     return;
674   }
675
676   // loop through the lines of the file until the beginning of data
677
678   TString title,date,subdetector,url,version,observations,system,units;
679   while (1) {
680     char pee=in.peek();
681     if (pee==EOF) break; 
682     TString line;
683     line.ReadLine(in);
684     if (line.Contains("Title:"))        title.ReadLine(in);
685     if (line.Contains("Date:"))         date.ReadLine(in);
686     if (line.Contains("Subdetector:"))  subdetector.ReadLine(in);
687     if (line.Contains("URL:"))          url.ReadLine(in);
688     if (line.Contains("Version:"))      version.ReadLine(in);
689     if (line.Contains("Observations:")) observations.ReadLine(in);
690     if (line.Contains("System:"))       system.ReadLine(in);
691     if (line.Contains("Units:"))        units.ReadLine(in);
692     if (line.Contains("Data:"))         break;
693   }
694
695   // check what we found so far (watch out, they have \r at the end)
696
697   std::cout<<"title .........."<<title<<std::endl;
698   std::cout<<"date ..........."<<date<<std::endl;
699   std::cout<<"subdetector ...."<<subdetector<<std::endl;
700   std::cout<<"url ............"<<url<<std::endl;
701   std::cout<<"version ........"<<version<<std::endl;
702   std::cout<<"observations ..."<<observations<<std::endl;
703   std::cout<<"system ........."<<system<<std::endl;
704   std::cout<<"units .........."<<units<<std::endl;
705
706   if (!subdetector.Contains("TRD")) {
707     AliWarning(Form("Not a TRD survey file, subdetector = %s",subdetector.Data()));
708     return;
709   }
710   double tocm = 0; // we want to have it in cm
711   if (units.Contains("mm"))      tocm = 0.1;
712   else if (units.Contains("cm")) tocm = 1.0;
713   else if (units.Contains("m"))  tocm = 100.0;
714   else if (units.Contains("pc")) tocm = 3.24078e-15;
715   else {
716     AliError(Form("unexpected units: %s",units.Data()));
717     return;
718   }
719   if (!system.Contains("ALICEPH")) {
720     AliError(Form("wrong system: %s, should be ALICEPH",system.Data()));
721     return;
722   }
723
724   // scan the rest of the file which should contain list of surveyed points
725   // for every point, decode the point name and store the numbers in the right 
726   // place in the arrays fSurveyX etc.
727
728   while (1) {
729     TString pna; // point name
730     char type, target;
731     double x,y,z,precision;
732     
733     in >> pna >> x >> y >> z >> type >> target >> precision;  
734     if (in.fail()) break;
735     Int_t i,j,k,l;
736     if (DecodeSurveyPointName(pna,i,j,k,l)) {
737       fSurveyX[i][j][k][l] = tocm*x;
738       fSurveyY[i][j][k][l] = tocm*y;
739       fSurveyZ[i][j][k][l] = tocm*z;
740       fSurveyEX[i][j][k][l] = precision/10; // "precision" is supposed to be in mm
741       fSurveyEY[i][j][k][l] = precision/10; // "precision" is supposed to be in mm
742       fSurveyEZ[i][j][k][l] = precision/10; // "precision" is supposed to be in mm
743       // if, at some point, separate precision numbers for x,y,z show up in the 
744       // survey reports the function will fail here
745       std::cout << "decoded "<<pna<<" "
746                 <<fSurveyX[i][j][k][l]<<" "
747                 <<fSurveyY[i][j][k][l]<<" "
748                 <<fSurveyZ[i][j][k][l]<<" "
749                 <<fSurveyEX[i][j][k][l]<<" "
750                 <<fSurveyEY[i][j][k][l]<<" "
751                 <<fSurveyEZ[i][j][k][l]<<" "<<std::endl;        
752     } else AliError(Form("cannot decode point name: %s",pna.Data()));
753   }
754   in.close();
755   TString info = "Survey "+title+" "+date+" "+url+" "+version+" "+observations;
756   info.ReplaceAll("\r","");
757   fComment.SetString(info.Data());
758                          
759 }
760
761 //_____________________________________________________________________________
762 void AliTRDalignment::ReadSurveyReport(AliSurveyObj *so) 
763 {
764   //
765   // Read survey report and store the numbers in fSurveyX, fSurveyY, fSurveyZ, 
766   // and fSurveyE.  Store the survey info in the fComment.
767   // Each supermodule has 8 survey points. The point names look like 
768   // TRD_sm08ah0 and have the following meaning. 
769   //
770   // sm00..17 mean supermodule 0 through 17, following the phi.
771   // Supermodule 00 is between phi=0 and phi=20 degrees.
772   //
773   // a or c denotes the anticlockwise and clockwise end of the supermodule
774   // in z. Clockwise end is where z is negative and where the muon arm sits.
775   //
776   // l or h denote low radius and high radius holes
777   //
778   // 0 or 1 denote the hole at smaller and at larger phi, respectively.
779   //
780
781   // read and process the data from the survey object
782
783   Int_t size = so->GetEntries();
784   printf("-> %d\n", size);
785
786   TString title        = so->GetReportTitle();
787   TString date         = so->GetReportDate();
788   TString subdetector  = so->GetDetector();
789   TString url          = so->GetURL();
790   TString report       = so->GetReportNumber();
791   TString version      = so->GetReportVersion();
792   TString observations = so->GetObservations();
793   TString system       = so->GetCoordSys();
794   TString units        = so->GetUnits();
795
796   // check what we found so far (watch out, they have \r at the end)
797
798   std::cout<<"title .........."<<title<<std::endl;
799   std::cout<<"date ..........."<<date<<std::endl;
800   std::cout<<"subdetector ...."<<subdetector<<std::endl;
801   std::cout<<"url ............"<<url<<std::endl;
802   std::cout<<"version ........"<<version<<std::endl;
803   std::cout<<"observations ..."<<observations<<std::endl;
804   std::cout<<"system ........."<<system<<std::endl;
805   std::cout<<"units .........."<<units<<std::endl;
806
807   if (!subdetector.Contains("TRD")) {
808     AliWarning(Form("Not a TRD survey file, subdetector = %s",subdetector.Data()));
809     return;
810   }
811   double tocm = 0; // we want to have it in cm
812   if (units.Contains("mm"))      tocm = 0.1;
813   else if (units.Contains("cm")) tocm = 1.0;
814   else if (units.Contains("m"))  tocm = 100.0;
815   else if (units.Contains("pc")) tocm = 3.24078e-15;
816   else {
817     AliError(Form("unexpected units: %s",units.Data()));
818     return;
819   }
820   if (!system.Contains("ALICEPH")) {
821     AliError(Form("wrong system: %s, should be ALICEPH",system.Data()));
822     return;
823   }
824
825   // for every survey point, decode the point name and store the numbers in 
826   // the right place in the arrays fSurveyX etc.
827
828   TObjArray *points = so->GetData();
829   for (int ip = 0; ip<points->GetEntries(); ++ip) {
830     AliSurveyPoint *po = (AliSurveyPoint *) points->At(ip);
831     TString pna = po->GetPointName();
832     Int_t i,j,k,l;
833     if (DecodeSurveyPointName(pna,i,j,k,l)) {
834       fSurveyX[i][j][k][l] = tocm*po->GetX();
835       fSurveyY[i][j][k][l] = tocm*po->GetY();
836       fSurveyZ[i][j][k][l] = tocm*po->GetZ();
837       fSurveyEX[i][j][k][l] = po->GetPrecisionX()/10; // "precision" is supposed to be in mm
838       fSurveyEY[i][j][k][l] = po->GetPrecisionY()/10;
839       fSurveyEZ[i][j][k][l] = po->GetPrecisionZ()/10;
840       std::cout << "decoded "<<pna<<" "
841                 <<fSurveyX[i][j][k][l]<<" "
842                 <<fSurveyY[i][j][k][l]<<" "
843                 <<fSurveyZ[i][j][k][l]<<" "
844                 <<fSurveyEX[i][j][k][l]<<" "
845                 <<fSurveyEY[i][j][k][l]<<" "
846                 <<fSurveyEZ[i][j][k][l]<<" "<<std::endl;        
847     } else AliError(Form("cannot decode point name: %s",pna.Data()));
848   }
849
850   TString info = "Survey "+title+" "+date+" "+url+" "+report+" "+version+" "+observations;
851   info.ReplaceAll("\r","");
852   fComment.SetString(info.Data());                       
853 }
854
855 //_____________________________________________________________________________
856 double AliTRDalignment::SurveyChi2(int i, double *a) {
857
858   //
859   // Compare the survey results to the ideal positions of the survey marks
860   // in the local frame of supermodule. When transforming, use the alignment 
861   // parameters a[6]. Return chi-squared.
862   //
863
864   if (!IsGeoLoaded()) return 0;
865   printf("Survey of supermodule %d\n",i);
866   AliAlignObjParams al(GetSmName(i),0,a[0],a[1],a[2],a[3],a[4],a[5],0);
867   TGeoPNEntry      *pne  = gGeoManager->GetAlignableEntry(GetSmName(i));
868   if (!pne) AliError(Form("no such physical node entry: %s",GetSmName(i)));
869   TGeoPhysicalNode *node = pne->GetPhysicalNode();
870   if (!node) AliError(Form("physical node entry %s has no physical node",GetSmName(i))); 
871
872   //  al.ApplyToGeometry();    
873   //  node = pne->GetPhysicalNode(); // changed in the meantime
874   //  TGeoHMatrix *ma = node->GetMatrix();
875
876   // a less destructive method (it does not modify geometry), gives the same result:
877
878   TGeoHMatrix *ma = new TGeoHMatrix();
879   al.GetLocalMatrix(*ma);
880   ma->MultiplyLeft(node->GetMatrix()); // global trafo, modified by a[]
881
882   double chi2=0;
883   printf("              sm   z   r  phi    x (lab phi)  y (lab z)   z (lab r)   all in cm\n");
884   for (int j=0; j<2; j++) for (int k=0; k<2; k++) for (int l=0; l<2; l++) {
885     if (fSurveyEX[i][j][k][l] == 0.0 
886         && fSurveyEY[i][j][k][l] == 0.0 
887         && fSurveyEZ[i][j][k][l] == 0.0) continue; // no data for this survey point
888     double master[3] = {fSurveyX[i][j][k][l],fSurveyY[i][j][k][l],fSurveyZ[i][j][k][l]};
889     double local[3];
890     ma->MasterToLocal(master,local);
891     double dx = local[0]-fSurveyX0[j][k][l];
892     double dy = local[1]-fSurveyY0[j][k][l];
893     double dz = local[2]-fSurveyZ0[j][k][l];
894     chi2 += dx*dx/fSurveyEX[i][j][k][l]/fSurveyEX[i][j][k][l];
895     chi2 += dy*dy/fSurveyEY[i][j][k][l]/fSurveyEY[i][j][k][l];
896     chi2 += dz*dz/fSurveyEZ[i][j][k][l]/fSurveyEZ[i][j][k][l];
897     printf("local survey %3d %3d %3d %3d %12.3f %12.3f %12.3f\n",i,j,k,l,local[0],local[1],local[2]);
898     printf("local ideal                  %12.3f %12.3f %12.3f\n",fSurveyX0[j][k][l],
899            fSurveyY0[j][k][l],fSurveyZ0[j][k][l]);
900     printf("difference                   %12.3f %12.3f %12.3f\n",dx,dy,dz);
901   }
902   printf("chi2 = %.2f\n",chi2);
903   return chi2;
904 }
905
906 //_____________________________________________________________________________
907 void trdAlignmentFcn(int &npar, double *g, double &f, double *par, int iflag) {
908
909   // 
910   // Standard function as needed by Minuit-like minimization procedures. 
911   // For the set of parameters par calculates and returns chi-squared.
912   //
913
914   // smuggle a C++ object into a C function
915   AliTRDalignment *alignment = (AliTRDalignment*) gMinuit->GetObjectFit(); 
916
917   f = alignment->SurveyChi2(par);
918   if (iflag==3) {}
919   if (npar) {} 
920   if (g) {} // no warnings about unused stuff...
921
922 }
923
924 //_____________________________________________________________________________
925 void AliTRDalignment::SurveyToAlignment(int i,char *flag) {
926
927   //
928   // Find the supermodule alignment parameters needed to make the survey 
929   // results coincide with the ideal positions of the survey marks.
930   // The string flag should look like "101000"; the six characters corresponds 
931   // to the six alignment parameters and 0/1 mean that the parameter should 
932   // be fixed/released in the fit. 
933
934   if (strlen(flag)!=6) {
935     AliError(Form("unexpected flag: %s",flag));
936     return;
937   }
938
939   printf("Finding alignment matrix for supermodule %d\n",i);
940   fIbuffer[0] = i; // store the sm number in the buffer so minuit can see it
941
942   TFitter fitter(100);
943   gMinuit->SetObjectFit(this);
944   fitter.SetFCN(trdAlignmentFcn);
945   fitter.SetParameter(0,"dx",0,0.5,0,0);
946   fitter.SetParameter(1,"dy",0,0.5,0,0);
947   fitter.SetParameter(2,"dz",0,0.5,0,0);
948   fitter.SetParameter(3,"rx",0,0.1,0,0);
949   fitter.SetParameter(4,"ry",0,0.1,0,0);
950   fitter.SetParameter(5,"rz",0,0.1,0,0);
951
952   for (int j=0; j<6; j++) if (flag[j]=='0') fitter.FixParameter(j);
953
954   double arglist[100];
955   arglist[0] = 2;
956   fitter.ExecuteCommand("SET PRINT", arglist, 1);
957   fitter.ExecuteCommand("SET ERR", arglist, 1);
958   arglist[0]=50;
959   //fitter.ExecuteCommand("SIMPLEX", arglist, 1);
960   fitter.ExecuteCommand("MINIMIZE", arglist, 1);
961   fitter.ExecuteCommand("CALL 3", arglist,0);
962   double a[6];
963   for (int j=0; j<6; j++) a[j] = fitter.GetParameter(j);
964   SetSm(i,a);
965   for (int j=0; j<6; j++) printf("%10.3f ",fitter.GetParameter(j));   
966   printf("\n");
967   for (int j=0; j<6; j++) printf("%10.3f ",fitter.GetParError(j));
968   printf("\n");
969
970 }
971
972 //_____________________________________________________________________________
973 void AliTRDalignment::ReadAny(char *filename) 
974 {
975   //
976   // read the alignment data from any kind of file
977   //
978
979   TString fist(filename);
980   if (fist.EndsWith(".txt")) ReadAscii(filename);
981   if (fist.EndsWith(".dat")) ReadAscii(filename);
982   if (fist.EndsWith(".root")) {
983     if (fist.Contains("Run")) ReadDB(filename);
984     else ReadRoot(filename);
985   }
986
987 }
988
989 //_____________________________________________________________________________
990 void AliTRDalignment::WriteAscii(char *filename) const
991 {
992   //
993   // store the alignment data on ascii file
994   //
995
996   FILE *fp = fopen(filename, "w");
997   if (!fp) {
998     AliError(Form("cannot open output file %s",filename));
999     return;
1000   }
1001
1002   PrintSm(fp);
1003   PrintCh(fp);
1004   
1005   fclose(fp);
1006
1007 }
1008
1009 //_____________________________________________________________________________
1010 void AliTRDalignment::WriteRoot(char *filename) 
1011 {
1012   //
1013   // store the alignment data on root file
1014   //
1015
1016   TClonesArray *ar = new TClonesArray("AliAlignObjParams",10000);
1017   NumbersToAr(ar);
1018   TFile fo(filename,"RECREATE");
1019   if (fo.IsOpen()) {
1020     fo.cd();
1021     fo.WriteObject(ar,"TRDAlignObjs","kSingleKey");
1022     fo.Close();
1023   } 
1024   else AliError(Form("cannot open output file %s",filename));
1025
1026   delete ar;
1027
1028 }
1029
1030 //_____________________________________________________________________________
1031 void AliTRDalignment::WriteDB(char *filename, int run0, int run1) 
1032 {
1033   //
1034   // dumping on a DB-like file
1035   //
1036
1037   TClonesArray   *ar = new TClonesArray("AliAlignObjParams",10000);
1038   NumbersToAr(ar);
1039   const Char_t *path = "TRD/Align/Data";
1040   AliCDBId id(path,run0,run1);
1041   AliCDBMetaData *md = new AliCDBMetaData();
1042   md->SetResponsible("Dariusz Miskowiec");
1043   md->SetComment(fComment.GetString().Data());
1044   AliCDBEntry    *e  = new AliCDBEntry(ar, id, md);
1045   TFile fi(filename,"RECREATE");
1046   if (fi.IsOpen()) {
1047     e->Write();
1048     fi.Close();
1049   } 
1050   else AliError(Form("cannot open input file %s",filename));
1051
1052   delete e;
1053   delete md;
1054   delete ar;
1055
1056   return;
1057
1058 }
1059
1060 //_____________________________________________________________________________
1061 void AliTRDalignment::WriteDB(char *db, char *path, int run0, int run1) 
1062 {
1063   //
1064   // store the alignment data in database
1065   //
1066
1067   TClonesArray   *ar      = new TClonesArray("AliAlignObjParams",10000);
1068   NumbersToAr(ar);
1069   AliCDBManager  *cdb     = AliCDBManager::Instance();
1070   AliCDBStorage  *storLoc = cdb->GetStorage(db);
1071   AliCDBMetaData *md      = new AliCDBMetaData();
1072   md->SetResponsible("Dariusz Miskowiec");
1073   md->SetComment(fComment.GetString().Data());
1074   AliCDBId id(path,run0,run1);
1075   storLoc->Put(ar,id,md);
1076   md->Delete();
1077   delete ar;
1078
1079 }
1080
1081 //_____________________________________________________________________________
1082 void AliTRDalignment::WriteGeo(char *filename) 
1083 {
1084   //
1085   // apply misalignment to current geometry and store the 
1086   // resulting geometry on a root file
1087   //
1088
1089   TClonesArray *ar = new TClonesArray("AliAlignObjParams",10000);
1090   NumbersToAr(ar);
1091   delete ar;
1092   gGeoManager->Export(filename);
1093
1094 }
1095
1096 //_____________________________________________________________________________
1097 double AliTRDalignment::GetSmRMS(int xyz) const 
1098 {
1099   //
1100   // rms fSm[][xyz]
1101   //
1102
1103   double s1 = 0.0;
1104   double s2 = 0.0;
1105   for (int i = 0; i < 18; i++) {
1106     s1 += fSm[i][xyz];
1107     s2 += fSm[i][xyz]*fSm[i][xyz];
1108   }
1109   double rms2 = s2/18.0 - s1*s1/18.0/18.0;
1110
1111   return rms2>0 ? sqrt(rms2) : 0.0;
1112
1113 }
1114
1115 //_____________________________________________________________________________
1116 double AliTRDalignment::GetChRMS(int xyz) const
1117 {
1118   //
1119   // rms fCh[][xyz]
1120   //
1121
1122   double s1 =0.0;
1123   double s2 =0.0;
1124   for (int i = 0; i < 540; i++) {
1125     s1 += fCh[i][xyz];
1126     s2 += fCh[i][xyz]*fCh[i][xyz];
1127   }
1128   double rms2 = s2/540.0 - s1*s1/540.0/540.0;
1129
1130   return rms2>0 ? sqrt(rms2) : 0.0;
1131
1132 }
1133
1134 //_____________________________________________________________________________
1135 void AliTRDalignment::PrintSmRMS() const
1136 {
1137   //
1138   // dump rms of fSm
1139   //
1140
1141   printf("       %11.4f %11.4f  %11.4f      %11.5f  %11.5f  %11.5f  supermodule rms\n"
1142         ,GetSmRMS(0),GetSmRMS(1),GetSmRMS(2),GetSmRMS(3),GetSmRMS(4),GetSmRMS(5));
1143
1144 }
1145
1146 //_____________________________________________________________________________
1147 void AliTRDalignment::PrintChRMS() const
1148 {
1149   //
1150   // dump rms of fCh
1151   //
1152
1153   printf("       %11.4f %11.4f  %11.4f      %11.5f  %11.5f  %11.5f  chamber rms\n"
1154         ,GetChRMS(0),GetChRMS(1),GetChRMS(2),GetChRMS(3),GetChRMS(4),GetChRMS(5));
1155
1156 }
1157
1158 //_____________________________________________________________________________
1159 void AliTRDalignment::ArToNumbers(TClonesArray *ar) 
1160 {
1161   //
1162   // for each of the alignment objects in array ar extract the six local 
1163   // alignment parameters; recognize by name to which supermodule or chamber 
1164   // the alignment object pertains; set the respective fSm or fCh
1165   //
1166
1167   ar->Sort();
1168   if (!IsGeoLoaded()) return;
1169   for (int i = 0; i < ar->GetEntries(); i++) {
1170     AliAlignObj *aao = (AliAlignObj *) ar->At(i);
1171     aao->ApplyToGeometry();
1172   }
1173   SetZero();
1174   ReadCurrentGeo();
1175
1176 }
1177
1178 //_____________________________________________________________________________
1179 void AliTRDalignment::NumbersToAr(TClonesArray *ar) 
1180 {
1181   //
1182   // build array of AliAlignObj objects based on fSm and fCh data
1183   // at the same time, apply misalignment to the currently loaded geometry
1184   // it is important to apply misalignment of supermodules before creating 
1185   // alignment objects for chambers
1186   //
1187
1188   if (!IsGeoLoaded()) return;
1189   TClonesArray &alobj = *ar;
1190   int nobj = 0;
1191   for (int i = 0; i <  18; i++) {      
1192       new(alobj[nobj]) AliAlignObjParams(GetSmName(i)
1193                                         ,0 
1194                                         ,fSm[i][0],fSm[i][1],fSm[i][2]
1195                                         ,fSm[i][3],fSm[i][4],fSm[i][5]
1196                                         ,0);
1197     ((AliAlignObj *) alobj[nobj])->ApplyToGeometry();
1198     nobj++;
1199   }
1200
1201   for (int i = 0; i < 540; i++) {
1202     if (gGeoManager->GetAlignableEntry(GetChName(i))) {
1203       new(alobj[nobj]) AliAlignObjParams(GetChName(i)
1204                                          ,GetVoi(i)
1205                                          ,fCh[i][0],fCh[i][1],fCh[i][2]
1206                                          ,fCh[i][3],fCh[i][4],fCh[i][5]
1207                                          ,0);
1208       ((AliAlignObj *) alobj[nobj])->ApplyToGeometry();
1209       nobj++;
1210     }
1211   }
1212   AliInfo("current geometry modified");
1213
1214 }
1215
1216 //_____________________________________________________________________________
1217 int AliTRDalignment::IsGeoLoaded() 
1218 {
1219   //
1220   // check whether a geometry is loaded
1221   // issue a warning if geometry is not ideal
1222   //
1223
1224   if (gGeoManager) {
1225     if (gGeoManager->GetListOfPhysicalNodes()->GetEntries()) AliWarning("current geometry is not ideal");
1226     return 1;
1227   } else {
1228     AliError("first load geometry by calling TGeoManager::Import(filename)");
1229     return 0;
1230   }
1231
1232 }
1233
1234 //_____________________________________________________________________________