]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/MUONSurveyCh4.C
Compilation on Windoiws/Cygwin
[u/mrichter/AliRoot.git] / MUON / MUONSurveyCh4.C
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 /// \ingroup macros
19 /// \file MUONSurveyCh8L.C
20 /// \brief Macro to process survey and photogrammetry data of chamber 4
21 /// 
22 /// Macro loads the survey data from .txt file using AliSurveyObj.
23 /// Macro uses AliMUONSurvey... classes.
24 /// The transformations of the detection elements or chambers can be obtained 
25 /// in three ways:
26 /// - A: Fit of local to global transformation using the fixed button targets.
27 /// - B: Fit a plane to the sticker targets -> psi, theta
28 ///      Use above psi and theta and fit remaining 4 parameters using the fixed 
29 ///      button targets
30 /// -  C: Fit a plane to the sticker targets -> psi, theta
31 ///      Use above psi and theta to calculate xc, yc, zc and phi by solving 
32 ///      the equations from a local to global transformation of the fixed button 
33 ///      targets
34 ///
35 /// Methods A and B are prefered to C, and B is better if sticker targets are 
36 /// available and lie on a plane!
37 /// For slats only methods B and C can be used.
38 /// Various histograms are filled and printed for monitoring.
39 /// MisAlignment object is then created.
40 ///
41 /// \author Javier Castillo
42
43 #if !defined(__CINT__) || defined(__MAKECINT__)
44
45 #include "AliMUONGeometryTransformer.h"
46 #include "AliMUONGeometryModuleTransformer.h"
47 #include "AliMUONGeometryDetElement.h"
48 #include "AliMUONGeometryMisAligner.h"
49 #include "AliMUONSurveyChamber.h"
50 #include "AliMUONSurveyDetElem.h"
51 #include "AliMUONSurveyUtil.h"
52
53 #include "AliSurveyObj.h"
54 #include "AliSurveyPoint.h"
55 #include "AliGeomManager.h"
56 #include "AliCDBManager.h"
57 #include "AliCDBMetaData.h"
58 #include "AliCDBId.h"
59 #include "AliAlignObjMatrix.h"
60 #include "AliAlignObj.h"
61
62 #include <TROOT.h>
63 #include <TGeoManager.h>
64 #include <TClonesArray.h>
65 #include <TObjArray.h>
66 #include <TArrayD.h>
67 #include <TMath.h>
68 #include <TString.h>
69 #include <TFitter.h>
70 #include <TH2.h>
71 #include <TF2.h>
72 #include <TGraphErrors.h>
73 #include <TCanvas.h>
74 #include <TPad.h>
75 #include <TLine.h>
76 #include <TPostScript.h>
77 #include <TPaveLabel.h>
78 #include <TStyle.h>
79 #include <TFile.h>
80 #include <TMatrixDSym.h>
81
82 #include <fstream>
83
84 #endif
85
86 void MUONSurveyCh4() {
87   
88   TString str;
89   TString title;
90   
91   Bool_t bMonitor = kTRUE;
92   Bool_t bOCDB = kTRUE;
93   Bool_t saveps = kFALSE;
94   const int cWidth = (int)(700*(29./21.));
95   const int cHeight = 700;
96   const int filetype  = 112; // landscape  
97
98   Int_t chId = 3;
99   Int_t nChs = 2;
100   Int_t nDetElems = 4;
101   Int_t nDetElemsI = 2;
102   //  Int_t nDetElemsO = 2;
103   Int_t iDetElemToDetElemId[4] = {400,401,402,403};
104   Int_t iDetElemPseudoIdToDetElem[4] = {0,1,2,3};
105   Int_t iDetElemsOfChamber[2][2] = {{0,3},
106                                     {1,2}};
107   
108   TObjArray *myChamberArray = new TObjArray(); 
109   myChamberArray->Add(new AliMUONSurveyChamber(chId));  
110   myChamberArray->Add(new AliMUONSurveyChamber(chId));  
111
112   AliMUONSurveyChamber *myChamberI = 0x0;
113   AliMUONSurveyChamber *myChamberO = 0x0;
114   AliMUONSurveyChamber *myChamber = 0x0;
115   AliMUONSurveyDetElem *myDetElem = 0x0;
116
117   myChamberI = (AliMUONSurveyChamber*)myChamberArray->At(0);
118   myChamberI->GetSurveyObj()->FillFromLocalFile("../Reports/AliceSt2_TC4_4971d.txt");
119   myChamberO = (AliMUONSurveyChamber*)myChamberArray->At(1);
120   myChamberO->GetSurveyObj()->FillFromLocalFile("../Reports/AliceSt2_TC4_4971d.txt");
121
122   myChamber = myChamberI;
123   myChamber->PrintSurveyReport();
124
125   // Chamber & DetElems button targets local coordinates
126   AliSurveyObj *lSO = new AliSurveyObj();    
127   lSO->FillFromLocalFile("$ALICE_ROOT/MUON/data/MUONTrackingTargetsLocal.txt");
128
129   // Set survey targets for chambers 
130   myChamberI->AddGButtonTargets(TString("T1_4520"),3);
131   myChamberI->AddStickerTargets(TString("T1_4620"),6);
132   myChamberI->AddLButtonTargets(lSO->GetData(),Form("%d_520",chId+1),3);
133   myChamberI->AddGButtonTargets(TString("T1_4510"),3);
134   myChamberI->AddStickerTargets(TString("T1_4610"),6);
135   myChamberI->AddLButtonTargets(lSO->GetData(),Form("%d_510",chId+1),3);
136   myChamberO->AddGButtonTargets(TString("T1_4520"),3);
137   myChamberO->AddStickerTargets(TString("T1_4620"),6);
138   myChamberO->AddLButtonTargets(lSO->GetData(),Form("%d_520",chId+1),3);
139   myChamberO->AddGButtonTargets(TString("T1_4510"),3);
140   myChamberO->AddStickerTargets(TString("T1_4610"),6);
141   myChamberO->AddLButtonTargets(lSO->GetData(),Form("%d_510",chId+1),3);
142
143   // Set survey targets for detection elements 
144   for (int iCh =0; iCh<=1; iCh++) {
145     myChamber = (AliMUONSurveyChamber*)myChamberArray->At(iCh);
146     for (int iDetElemI=0; iDetElemI<nDetElemsI; iDetElemI++){
147       Int_t iDetElem = iDetElemsOfChamber[iCh][iDetElemI];
148       myChamber->AddSurveyDetElem(iDetElemToDetElemId[iDetElem]);
149       TString baseName;
150       myDetElem =  myChamber->GetDetElem(iDetElemI);
151       if (myDetElem) {
152         baseName = Form("T1_453%d",iDetElem+1);
153         myDetElem->AddGButtonTargets(baseName,3);
154         baseName = Form("T1_463%d",iDetElem+1);
155         myDetElem->AddStickerTargets(baseName,6);
156         baseName = Form("%d_53%d",chId+1,iDetElem+1);
157         myDetElem->AddLButtonTargets(lSO->GetData(),baseName,3);
158       }
159     }
160   }
161   printf("All targets added! \n");
162
163   Double_t **lCenDetElem = new Double_t*[nDetElems+nChs];
164   Double_t **lRotDetElem = new Double_t*[nDetElems+nChs];
165   Double_t **lDiffCenDetElem0 = new Double_t*[nDetElems+nChs];
166   Double_t **lDiffRotDetElem0 = new Double_t*[nDetElems+nChs];
167   Double_t **lDiffThCenDetElem0 = new Double_t*[nDetElems+nChs];
168   Double_t **lDiffThRotDetElem0 = new Double_t*[nDetElems+nChs];
169   Double_t **lDeltaDiffCenDetElem0 = new Double_t*[nDetElems+nChs];
170   Double_t **lDeltaDiffRotDetElem0 = new Double_t*[nDetElems+nChs];
171
172   for (int iDetElem=0; iDetElem<nDetElems+nChs; iDetElem++) {
173     lCenDetElem[iDetElem] = new Double_t[3];
174     lRotDetElem[iDetElem] = new Double_t[3];
175     lDiffCenDetElem0[iDetElem] = new Double_t[3];
176     lDiffRotDetElem0[iDetElem] = new Double_t[3];
177     lDiffThCenDetElem0[iDetElem] = new Double_t[3];
178     lDiffThRotDetElem0[iDetElem] = new Double_t[3];
179     lDeltaDiffCenDetElem0[iDetElem] = new Double_t[3];
180     lDeltaDiffRotDetElem0[iDetElem] = new Double_t[3];
181   }
182
183   TGeoCombiTrans dtrfDetElem[nDetElems+nChs];
184   TGeoCombiTrans localTrfDetElem[nDetElems+nChs];
185   TGeoCombiTrans localTrfThDetElem[nDetElems+nChs];
186
187   // Import TGeo geometry 
188   char* geoFilename = "geometry.root";
189   if ( ! AliGeomManager::GetGeometry() ) {
190     AliGeomManager::LoadGeometry(geoFilename);
191     if (! AliGeomManager::GetGeometry() ) {
192       printf("MUONSurveyCh%d: getting geometry from file %s failed\n", chId+1,geoFilename);
193       return;
194     }
195     cout << "geometry imported" << endl;
196   }
197
198   AliMUONGeometryTransformer *transform = new AliMUONGeometryTransformer();
199   transform->LoadGeometryData();
200
201   TGeoCombiTrans trfThChamber;
202   TGeoCombiTrans trfThDetElem;
203
204   for (int iCh =0; iCh<=1; iCh++) {
205     myChamber = (AliMUONSurveyChamber*)myChamberArray->At(iCh);
206
207     trfThChamber = TGeoCombiTrans(*transform->GetModuleTransformerByDEId(iDetElemToDetElemId[iCh*nDetElemsI])->GetTransformation());
208     trfThChamber.SetTranslation(trfThChamber.GetTranslation()[0]*10,
209                                 trfThChamber.GetTranslation()[1]*10,
210                                 trfThChamber.GetTranslation()[2]*10);
211     trfThChamber.Print();
212     myChamber->SetLocalTransformation(new TGeoCombiTrans(trfThChamber),kTRUE);
213
214     // Set Chamber plane function
215     cout << "Setting plane for Chamber" << iCh+1 << " ..." << endl;
216     myChamber->SetPlane(Form("fChamber%d",iCh+1));
217     myChamber->SetPlaneParameters(0.,0.,0.);
218
219     // Fit a plane to sticker targets
220     Double_t lCChi2 = myChamber->FitPlane();
221     cout << "... chi2 = " << lCChi2 << " ..." << endl; 
222     
223     // Get best transformation from fitting method 
224     // (local to global transformation)
225     cout << "Trying fitting method for chamber " << iCh << endl;
226     myChamber->SurveyToAlign();    
227     //    myChamber->SurveyToAlign(myChamber->GetPlane()->GetParameter(0),myChamber->GetPlane()->GetParameter(1),myChamber->GetPlane()->GetParError(0),myChamber->GetPlane()->GetParError(1));    
228     
229     for (int iDetElemI=0; iDetElemI<nDetElemsI; iDetElemI++){
230       myDetElem =  myChamber->GetDetElem(iDetElemI);
231       Int_t iDetElem = iDetElemsOfChamber[iCh][iDetElemI];
232
233       trfThDetElem.Clear();
234       trfThDetElem = TGeoCombiTrans(*transform->GetDetElement(iDetElemToDetElemId[iDetElem])->GetLocalTransformation());
235       trfThDetElem.SetTranslation(trfThDetElem.GetTranslation()[0]*10,
236                                   trfThDetElem.GetTranslation()[1]*10,
237                                   trfThDetElem.GetTranslation()[2]*10);
238       trfThDetElem.Print();
239
240       myDetElem->SetLocalTransformation(new TGeoCombiTrans(trfThDetElem),kTRUE);
241
242       for (int iCor=0; iCor<3; iCor++){
243         lCenDetElem[iDetElem][iCor] = 0;
244         lRotDetElem[iDetElem][iCor] = 0;
245       }
246
247       if (bMonitor){
248         // MONITOR: Draw graph with sticker targets for plane fitting
249         myDetElem->DrawSTargets();
250         gPad->Update();
251       }
252
253       // Get DetElem transformation. 
254       // Psi and Theta are obtained by fitting a plane to the sticker targets.
255       // Then Xc, Yc, Zc and Phi are obtained by solving the equations to the ref.
256       // syst. transformation of the button targets
257
258       // Set DetElem plane function
259       cout << "Setting plane for DetElem" << iDetElem+1 << " ..." << endl;
260       myDetElem->SetPlane(Form("fDetElem%d",iDetElem+1));
261       myDetElem->SetPlaneParameters(0.,0.,30.)
262 ;
263       // Fit a plane to sticker targets
264       Double_t lChi2 = myDetElem->FitPlane();
265       cout << "... chi2 = " << lChi2 << " ..." << endl; 
266       
267       lRotDetElem[iDetElem][0]=TMath::ATan(myDetElem->GetPlane()->GetParameter(0));
268       lRotDetElem[iDetElem][1]=TMath::ATan(myDetElem->GetPlane()->GetParameter(1));
269       
270       // Calculate Mean transformation using previous plane fit 
271       // and pairs of button targets
272       myDetElem->CalculateMeanTransf(lCenDetElem[iDetElem],lRotDetElem[iDetElem]);
273     
274       cout << "DetElem" << iDetElem+1 << " : mycen(" << lCenDetElem[iDetElem][0] << "," << lCenDetElem[iDetElem][1] << "," << lCenDetElem[iDetElem][2] << "); rot(" << lRotDetElem[iDetElem][0] << "," << lRotDetElem[iDetElem][1] << "," << lRotDetElem[iDetElem][2] << ")"  << endl;          
275       
276      
277       // Get best transformation from fitting method 
278       // (local to global transformation)
279       cout << "Trying fitting method for DetElemId " << iDetElemToDetElemId[iDetElem] << endl;
280       //      myDetElem->SurveyToAlign();     
281       myDetElem->SurveyToAlign(lRotDetElem[iDetElem][0],lRotDetElem[iDetElem][1],myDetElem->GetPlane()->GetParError(0),myDetElem->GetPlane()->GetParError(1));    
282     }
283   }
284
285   // Print found transformation of Detection Element (plane fit + loop)        
286   for (int iCh =0; iCh<=1; iCh++) {
287     for (int iDetElemI=0; iDetElemI<nDetElemsI; iDetElemI++){
288       Int_t iDetElem = iDetElemsOfChamber[iCh][iDetElemI];
289       cout << "DetElem" << iDetElem+1 << " : mycen(" << lCenDetElem[iDetElem][0] << "," << lCenDetElem[iDetElem][1] << "," << lCenDetElem[iDetElem][2] << "); rot(" << lRotDetElem[iDetElem][0] << "," << lRotDetElem[iDetElem][1] << "," << lRotDetElem[iDetElem][2] << ")"  << endl;          
290     }
291   }
292
293   // Print Theoretical transformation of Detection Element
294   for (int iCh =0; iCh<=1; iCh++) {
295     myChamber = (AliMUONSurveyChamber*)myChamberArray->At(iCh);
296     for (int iDetElemI=0; iDetElemI<nDetElemsI; iDetElemI++){
297       myChamber->GetDetElem(iDetElemI)->PrintLocalTrf();
298     }
299   }
300
301   // Print found delta transformation of Detection Element
302   for (int iCh =0; iCh<=1; iCh++) {
303     myChamber = (AliMUONSurveyChamber*)myChamberArray->At(iCh);      
304     for (int iDetElemI=0; iDetElemI<nDetElemsI; iDetElemI++){
305       myChamber->GetDetElem(iDetElemI)->PrintAlignTrf();
306     }
307   }
308
309   //
310   // Compare transformations to expected ones 
311   //
312   Int_t iDetElemToPos[18] = {0, 1, 2, 3};
313  
314   TGraphErrors *gDeltaDiffCenXDetElem0 = new TGraphErrors(nDetElems);
315   TGraphErrors *gDeltaDiffCenYDetElem0 = new TGraphErrors(nDetElems);
316   TGraphErrors *gDeltaDiffCenZDetElem0 = new TGraphErrors(nDetElems);
317   TGraphErrors *gDeltaDiffPsiDetElem0 = new TGraphErrors(nDetElems);
318   TGraphErrors *gDeltaDiffThtDetElem0 = new TGraphErrors(nDetElems);
319   TGraphErrors *gDeltaDiffPhiDetElem0 = new TGraphErrors(nDetElems);
320
321   for (int iCh =0; iCh<=1; iCh++) {
322     myChamber = (AliMUONSurveyChamber*)myChamberArray->At(iCh);      
323     for (int iDetElemI=0; iDetElemI<nDetElemsI; iDetElemI++){
324       myChamber->GetDetElem(iDetElemI)->GetAlignTrf()->Print();
325     }
326   }
327
328   for (int iCh =0; iCh<=1; iCh++) {
329     myChamber = (AliMUONSurveyChamber*)myChamberArray->At(iCh);
330     // Store delta transformations to use for CDB entry creation
331     dtrfDetElem[nDetElems+iCh].Clear();
332     dtrfDetElem[nDetElems+iCh] = *(myChamber->GetAlignTrf());
333
334     // Those are for checks and visualizations
335     localTrfDetElem[nDetElems+iCh].Clear();
336     localTrfDetElem[nDetElems+iCh] = (*(myChamber->GetLocalTrf())*(*(myChamber->GetAlignTrf())));
337     localTrfDetElem[nDetElems+iCh].Print();
338     lDiffCenDetElem0[nDetElems+iCh] = (Double_t*)localTrfDetElem[nDetElems+iCh].GetTranslation();
339     AliMUONSurveyUtil::MatrixToAngles(localTrfDetElem[nDetElems+iCh].GetRotationMatrix(),lDiffRotDetElem0[nDetElems+iCh]);
340     
341     localTrfThDetElem[nDetElems+iCh].Clear();
342     localTrfThDetElem[nDetElems+iCh] = (*(myChamber->GetLocalTrf()));
343     localTrfThDetElem[nDetElems+iCh].Print();
344     lDiffThCenDetElem0[nDetElems+iCh] = (Double_t*)localTrfThDetElem[nDetElems+iCh].GetTranslation();
345     AliMUONSurveyUtil::MatrixToAngles(localTrfThDetElem[nDetElems+iCh].GetRotationMatrix(),lDiffThRotDetElem0[nDetElems+iCh]);
346
347     for (int iCor=0; iCor<3; iCor++){
348       lDeltaDiffCenDetElem0[nDetElems+iCh][iCor] = lDiffCenDetElem0[nDetElems+iCh][iCor]-lDiffThCenDetElem0[nDetElems+iCh][iCor];
349       lDeltaDiffRotDetElem0[nDetElems+iCh][iCor] = lDiffRotDetElem0[nDetElems+iCh][iCor]-lDiffThRotDetElem0[nDetElems+iCh][iCor];
350       if (lDeltaDiffRotDetElem0[nDetElems+iCh][iCor]>TMath::Pi()) lDeltaDiffRotDetElem0[nDetElems+iCh][iCor]-=TMath::TwoPi();
351       if (lDeltaDiffRotDetElem0[nDetElems+iCh][iCor]<-TMath::Pi()) lDeltaDiffRotDetElem0[nDetElems+iCh][iCor]+=TMath::TwoPi();
352     }      
353
354     for (int iDetElemI=0; iDetElemI<nDetElemsI; iDetElemI++){
355       Int_t iDetElem = iDetElemsOfChamber[iCh][iDetElemI];
356       myDetElem =  myChamber->GetDetElem(iDetElemI);
357       // Store delta transformations to use for CDB entry creation
358       dtrfDetElem[iDetElem].Clear();
359       dtrfDetElem[iDetElem] = *(myDetElem->GetAlignTrf());
360
361       // Those are for checks and visualizations
362       localTrfDetElem[iDetElem].Clear();
363       localTrfDetElem[iDetElem] = (*(myDetElem->GetLocalTrf())*(*(myDetElem->GetAlignTrf())));
364       //      localTrfDetElem[iDetElem] = (*(myDetElem->GetBaseTrf())*(*(myDetElem->GetAlignTrf())));
365       localTrfDetElem[iDetElem].Print();
366       lDiffCenDetElem0[iDetElem] = (Double_t*)localTrfDetElem[iDetElem].GetTranslation();
367       AliMUONSurveyUtil::MatrixToAngles(localTrfDetElem[iDetElem].GetRotationMatrix(),lDiffRotDetElem0[iDetElem]);
368 //       lDiffCenDetElem0[iDetElem] = lCenDetElem[iDetElem];
369 //       lDiffRotDetElem0[iDetElem] = lRotDetElem[iDetElem];
370
371       localTrfThDetElem[iDetElem].Clear();
372       localTrfThDetElem[iDetElem] = (*(myDetElem->GetLocalTrf()));
373       //      localTrfThDetElem[iDetElem] = (*(myDetElem->GetBaseTrf()));
374       localTrfThDetElem[iDetElem].Print();
375       lDiffThCenDetElem0[iDetElem] = (Double_t*)localTrfThDetElem[iDetElem].GetTranslation();
376       AliMUONSurveyUtil::MatrixToAngles(localTrfThDetElem[iDetElem].GetRotationMatrix(),lDiffThRotDetElem0[iDetElem]);
377       
378       for (int iCor=0; iCor<3; iCor++){
379         lDeltaDiffCenDetElem0[iDetElem][iCor] = lDiffCenDetElem0[iDetElem][iCor]-lDiffThCenDetElem0[iDetElem][iCor];
380         lDeltaDiffRotDetElem0[iDetElem][iCor] = lDiffRotDetElem0[iDetElem][iCor]-lDiffThRotDetElem0[iDetElem][iCor];
381         if (lDeltaDiffRotDetElem0[iDetElem][iCor]>TMath::Pi()) lDeltaDiffRotDetElem0[iDetElem][iCor]-=TMath::TwoPi();
382         if (lDeltaDiffRotDetElem0[iDetElem][iCor]<-TMath::Pi()) lDeltaDiffRotDetElem0[iDetElem][iCor]+=TMath::TwoPi();
383       }
384
385       gDeltaDiffCenXDetElem0->SetPoint(iDetElem,lDeltaDiffCenDetElem0[iDetElem][0],iDetElemToPos[iDetElem]+1);
386       gDeltaDiffCenYDetElem0->SetPoint(iDetElem,lDeltaDiffCenDetElem0[iDetElem][1],iDetElemToPos[iDetElem]+1);
387       gDeltaDiffCenZDetElem0->SetPoint(iDetElem,lDeltaDiffCenDetElem0[iDetElem][2],iDetElemToPos[iDetElem]+1);
388       gDeltaDiffPsiDetElem0->SetPoint(iDetElem,1e3*lDeltaDiffRotDetElem0[iDetElem][0],iDetElemToPos[iDetElem]+1);
389       gDeltaDiffThtDetElem0->SetPoint(iDetElem,1e3*lDeltaDiffRotDetElem0[iDetElem][1],iDetElemToPos[iDetElem]+1);
390       gDeltaDiffPhiDetElem0->SetPoint(iDetElem,1e3*lDeltaDiffRotDetElem0[iDetElem][2],iDetElemToPos[iDetElem]+1);
391       gDeltaDiffCenXDetElem0->SetPointError(iDetElem,myDetElem->GetFitter()->GetParError(0),0.);
392       gDeltaDiffCenYDetElem0->SetPointError(iDetElem,myDetElem->GetFitter()->GetParError(1),0.);
393       gDeltaDiffCenZDetElem0->SetPointError(iDetElem,myDetElem->GetFitter()->GetParError(2),0.);
394       gDeltaDiffPsiDetElem0->SetPointError(iDetElem,1e3*myDetElem->GetFitter()->GetParError(3),0.);
395       gDeltaDiffThtDetElem0->SetPointError(iDetElem,1e3*myDetElem->GetFitter()->GetParError(4),0.);
396       gDeltaDiffPhiDetElem0->SetPointError(iDetElem,1e3*myDetElem->GetFitter()->GetParError(5),0.);  
397     }
398   }
399
400   // Apply the found misalignments to the geometry
401   AliMUONGeometryTransformer *newTransform = AliMUONSurveyUtil::ReAlign(transform,chId,nDetElems,iDetElemPseudoIdToDetElem,dtrfDetElem,true); 
402   newTransform->WriteTransformations(Form("transform2ReAlignSurveyCh%d.dat",chId+1));
403
404   if(bOCDB){
405     // Generate realigned data in local cdb
406     const TClonesArray* array = newTransform->GetMisAlignmentData();
407     
408     // Set the alignment resolution in the align objects for this chamber
409     Double_t chResX = 0.1*myChamberI->GetAlignResX();
410     Double_t chResY = 0.1*myChamberI->GetAlignResY();
411     Double_t deResX = 0.1*(myChamberI->GetMeanDetElemAlignResX()+myChamberO->GetMeanDetElemAlignResX())/2.;
412     Double_t deResY = 0.1*(myChamberI->GetMeanDetElemAlignResY()+myChamberO->GetMeanDetElemAlignResY())/2.;
413     printf("Chamber alignment resolutions: resX=%f , resY=%f\n",chResX,chResY); 
414     printf("Detection Elements alignment resolutions: resX=%f , resY=%f\n",deResX,deResY); 
415     chResX = TMath::Sqrt(0.1*0.1+chResX*chResX);
416     chResY = TMath::Sqrt(0.1*0.1+chResY*chResY);
417     AliMUONSurveyUtil::SetAlignmentResolution(array,chId,chResX,chResY,deResX,deResY);
418     
419     // CDB manager
420     AliCDBManager* cdbManager = AliCDBManager::Instance();
421     cdbManager->SetDefaultStorage(Form("local://ReAlignSurveyCh%dCDB",chId+1));
422     
423     AliCDBMetaData* cdbData = new AliCDBMetaData();
424     cdbData->SetResponsible("Dimuon Offline project");
425     cdbData->SetComment("MUON alignment objects with survey realignment");
426     AliCDBId id("MUON/Align/Data", 0, AliCDBRunRange::Infinity()); 
427     cdbManager->Put(const_cast<TClonesArray*>(array), id, cdbData);
428   }
429
430   for(Int_t iCor=0; iCor<3; iCor++){
431     for(Int_t iDetElem=0; iDetElem<nDetElems; iDetElem++){
432       cout << lDeltaDiffCenDetElem0[iDetElem][iCor] << " " << lDiffCenDetElem0[iDetElem][iCor] << " " << lDiffThCenDetElem0[iDetElem][iCor] << endl;
433     }
434     cout << endl;
435   }
436   cout << endl;
437   for(Int_t iCor=0; iCor<3; iCor++){
438     for(Int_t iDetElem=0; iDetElem<nDetElems; iDetElem++){
439       cout << lDeltaDiffRotDetElem0[iDetElem][iCor] << " " << lDiffRotDetElem0[iDetElem][iCor] << " " << lDiffThRotDetElem0[iDetElem][iCor] << endl;
440     }
441     cout << endl;
442   }
443
444   TH1F *myDetElemDeltaDiffCenX = new TH1F("myDetElemDeltaDiffCenX","myDetElemDeltaDiffCenX",100,-10,10);
445   myDetElemDeltaDiffCenX->SetMaximum(nDetElems+1);
446   myDetElemDeltaDiffCenX->SetMinimum(0);
447   TH1F *myDetElemDeltaDiffCenY = new TH1F("myDetElemDeltaDiffCenY","myDetElemDeltaDiffCenY",100,-10,10);
448   myDetElemDeltaDiffCenY->SetMaximum(nDetElems+1);
449   myDetElemDeltaDiffCenY->SetMinimum(0);
450   TH1F *myDetElemDeltaDiffCenZ = new TH1F("myDetElemDeltaDiffCenZ","myDetElemDeltaDiffCenZ",100,-20,20);
451   myDetElemDeltaDiffCenZ->SetMaximum(nDetElems+1);
452   myDetElemDeltaDiffCenZ->SetMinimum(0);
453
454   TH1F *myDetElemDeltaDiffRotX = new TH1F("myDetElemDeltaDiffRotX","myDetElemDeltaDiffRotX",100,-15,15);
455   myDetElemDeltaDiffRotX->SetMaximum(nDetElems+1);
456   myDetElemDeltaDiffRotX->SetMinimum(0);
457   TH1F *myDetElemDeltaDiffRotY = new TH1F("myDetElemDeltaDiffRotY","myDetElemDeltaDiffRotY",100,-15,15);
458   myDetElemDeltaDiffRotY->SetMaximum(nDetElems+1);
459   myDetElemDeltaDiffRotY->SetMinimum(0);
460   TH1F *myDetElemDeltaDiffRotZ = new TH1F("myDetElemDeltaDiffRotZ","myDetElemDeltaDiffRotZ",100,-5,5);
461   myDetElemDeltaDiffRotZ->SetMaximum(nDetElems+1);
462   myDetElemDeltaDiffRotZ->SetMinimum(0);
463
464   //
465   // ******** Starting plots 
466   //
467   TCanvas *canvas;
468   TPad *pad;
469   TPaveLabel *theTitle;
470   gStyle->SetPalette(1);
471
472   TPostScript *ps = 0;
473   
474   if( saveps ){
475     ps = new TPostScript(Form("surveyChamber%d",chId+1),filetype); 
476     ps->NewPage();
477   }
478
479   // Observed misalignments
480   str = Form("Chamber %d",chId+1);
481   TCanvas *cvn2 = new TCanvas("cvn2",str,cWidth,cHeight);
482   canvas = cvn2;
483   canvas->Range(0,0,21,29);
484   
485   title = Form(" MisAlignments Chamber %d - PL2G - In Single Frame - DE2Ok - DE3Ok - DE4Ok ",chId+1);
486   TPaveLabel *theTitle2 = new TPaveLabel(3,27.0,18,28.5,title,"br");
487   theTitle = theTitle2;
488   theTitle->SetFillColor(18);
489   theTitle->SetTextFont(32);
490   theTitle->SetTextSize(0.4);
491   theTitle->SetTextColor(1);
492   theTitle->Draw();
493  
494   TPad *pad2 = new TPad("pad2","pad2",0.01,0.01,0.98,0.91,0);
495   pad = pad2;
496   pad->Draw();
497   TLine *ch0Line = new TLine(0,1,0,2);
498   TLine *ch1Line = new TLine(0,1,0,2);
499   ch1Line->SetLineStyle(2);
500   pad->Divide(3,2);
501   
502   pad->cd(1);
503   myDetElemDeltaDiffCenX->Draw();
504   myDetElemDeltaDiffCenX->SetXTitle("#Delta[xc_{i}^{m}-xc_{i}^{th}] (mm)");
505   myDetElemDeltaDiffCenX->SetYTitle("DetElem arbitrary ordering");
506   gDeltaDiffCenXDetElem0->SetMarkerStyle(20);
507   gDeltaDiffCenXDetElem0->Draw("P");
508   ch0Line->DrawLine(lDeltaDiffCenDetElem0[nDetElems+0][0],0.5,lDeltaDiffCenDetElem0[nDetElems+0][0],1.5);
509   ch0Line->DrawLine(lDeltaDiffCenDetElem0[nDetElems+0][0],3.5,lDeltaDiffCenDetElem0[nDetElems+0][0],4.5);
510   ch1Line->DrawLine(lDeltaDiffCenDetElem0[nDetElems+1][0],1.5,lDeltaDiffCenDetElem0[nDetElems+1][0],3.5);
511
512   pad->cd(2);
513   myDetElemDeltaDiffCenY->Draw();
514   myDetElemDeltaDiffCenY->SetXTitle("#Delta[yc_{i}^{m}-yc_{i}^{th}] (mm)");
515   myDetElemDeltaDiffCenY->SetYTitle("DetElem arbitrary ordering");
516   gDeltaDiffCenYDetElem0->SetMarkerStyle(20);
517   gDeltaDiffCenYDetElem0->Draw("P");
518   ch0Line->DrawLine(lDeltaDiffCenDetElem0[nDetElems+0][1],0.5,lDeltaDiffCenDetElem0[nDetElems+0][1],1.5);
519   ch0Line->DrawLine(lDeltaDiffCenDetElem0[nDetElems+0][1],3.5,lDeltaDiffCenDetElem0[nDetElems+0][1],4.5);
520   ch1Line->DrawLine(lDeltaDiffCenDetElem0[nDetElems+1][1],1.5,lDeltaDiffCenDetElem0[nDetElems+1][1],3.5);
521
522   pad->cd(3);
523   myDetElemDeltaDiffCenZ->Draw();
524   myDetElemDeltaDiffCenZ->SetXTitle("#Delta[zc_{i}^{m}-zc_{i}^{th}] (mm)");
525   myDetElemDeltaDiffCenZ->SetYTitle("DetElem arbitrary ordering");
526   gDeltaDiffCenZDetElem0->SetMarkerStyle(20);
527   gDeltaDiffCenZDetElem0->Draw("P");
528   ch0Line->DrawLine(lDeltaDiffCenDetElem0[nDetElems+0][2],0.5,lDeltaDiffCenDetElem0[nDetElems+0][2],1.5);
529   ch0Line->DrawLine(lDeltaDiffCenDetElem0[nDetElems+0][2],3.5,lDeltaDiffCenDetElem0[nDetElems+0][2],4.5);
530   ch1Line->DrawLine(lDeltaDiffCenDetElem0[nDetElems+1][2],1.5,lDeltaDiffCenDetElem0[nDetElems+1][2],3.5);
531
532   pad->cd(4);
533   myDetElemDeltaDiffRotX->Draw();
534   myDetElemDeltaDiffRotX->SetXTitle("#Delta[#psi_{i}^{m}-#psi_{i}^{th}] (mrad)");
535   myDetElemDeltaDiffRotX->SetYTitle("DetElem arbitrary ordering");
536   gDeltaDiffPsiDetElem0->SetMarkerStyle(20);
537   gDeltaDiffPsiDetElem0->Draw("P");
538   ch0Line->DrawLine(1e3*lDeltaDiffRotDetElem0[nDetElems+0][0],0.5,1e3*lDeltaDiffRotDetElem0[nDetElems+0][0],1.5);
539   ch0Line->DrawLine(1e3*lDeltaDiffRotDetElem0[nDetElems+0][0],3.5,1e3*lDeltaDiffRotDetElem0[nDetElems+0][0],4.5);
540   ch1Line->DrawLine(1e3*lDeltaDiffRotDetElem0[nDetElems+1][0],1.5,1e3*lDeltaDiffRotDetElem0[nDetElems+1][0],3.5);
541
542   pad->cd(5);
543   myDetElemDeltaDiffRotY->Draw();
544   myDetElemDeltaDiffRotY->SetXTitle("#Delta[#theta_{i}^{m}-#theta_{i}^{th}] (mrad)");
545   myDetElemDeltaDiffRotY->SetYTitle("DetElem arbitrary ordering");
546   gDeltaDiffThtDetElem0->SetMarkerStyle(20);
547   gDeltaDiffThtDetElem0->Draw("P");
548   ch0Line->DrawLine(1e3*lDeltaDiffRotDetElem0[nDetElems+0][1],0.5,1e3*lDeltaDiffRotDetElem0[nDetElems+0][1],1.5);
549   ch0Line->DrawLine(1e3*lDeltaDiffRotDetElem0[nDetElems+0][1],3.5,1e3*lDeltaDiffRotDetElem0[nDetElems+0][1],4.5);
550   ch1Line->DrawLine(1e3*lDeltaDiffRotDetElem0[nDetElems+1][1],1.5,1e3*lDeltaDiffRotDetElem0[nDetElems+1][1],3.5);
551
552   pad->cd(6);
553   myDetElemDeltaDiffRotZ->Draw();
554   myDetElemDeltaDiffRotZ->SetXTitle("#Delta[#phi_{i}^{m}-#phi_{i}^{th}] (mrad)");
555   myDetElemDeltaDiffRotZ->SetYTitle("DetElem arbitrary ordering");
556   gDeltaDiffPhiDetElem0->SetMarkerStyle(20);
557   gDeltaDiffPhiDetElem0->Draw("P");
558   ch0Line->DrawLine(1e3*lDeltaDiffRotDetElem0[nDetElems+0][2],0.5,1e3*lDeltaDiffRotDetElem0[nDetElems+0][2],1.5);
559   ch0Line->DrawLine(1e3*lDeltaDiffRotDetElem0[nDetElems+0][2],3.5,1e3*lDeltaDiffRotDetElem0[nDetElems+0][2],4.5);
560   ch1Line->DrawLine(1e3*lDeltaDiffRotDetElem0[nDetElems+1][2],1.5,1e3*lDeltaDiffRotDetElem0[nDetElems+1][2],3.5);
561
562   pad->Update();
563  
564   if(bMonitor){    
565     // MONITOR: Histograms for monitoring
566     TH2F *hCPSTa = new TH2F("hCPSTa","hCPSTa",100,-2000,2000,100,-2000,2000);
567     TH2F *hCPSTc = new TH2F("hCPSTc","hCPSTc",100,-2000,2000,100,-2000,2000);
568     TH2F *hSSTa = new TH2F("hSSTa","hSSTa",100,-2000,2000,100,-2000,2000);
569     TH2F *hSSTc = new TH2F("hSSTc","hSSTc",100,-2000,2000,100,-2000,2000);
570     TH2F *hSSTap = new TH2F("hSSTap","hSSTap",800,-2000,2000,800,-2000,2000);
571     TH2F *hSSTcp = new TH2F("hSSTcp","hSSTcp",800,-2000,2000,800,-2000,2000);
572     
573     // MONITOR: Fill histograms with chambers and slats sticker target positions
574     for (int iCh =0; iCh<=1; iCh++) {
575       myChamber = (AliMUONSurveyChamber*)myChamberArray->At(iCh);            
576       cout << "Filling histograms of sticker target for chamber" << iCh+1 << " ..." << endl;
577       //    myChamber->FillCPSTHistograms(TString("C"),hCPSTc);
578       myChamber->FillDESTHistograms(TString("T1_46"),hSSTc);
579
580       for (int iDetElemI=0; iDetElemI<nDetElemsI; iDetElemI++){
581         //      Int_t iDetElem = iDetElemsOfChamber[iCh][iDetElemI];
582         myDetElem =  myChamber->GetDetElem(iDetElemI);
583         // MONITOR: Fill slat plane for fit monitor.
584         Double_t pl[3] = {0};
585         Double_t pg[3] = {0};
586         AliSurveyPoint *pointSBT0 = myDetElem->GetLButtonTarget(0);
587         AliSurveyPoint *pointSBT1 = myDetElem->GetLButtonTarget(1);
588         if(pointSBT0&&pointSBT1) {
589           if (pointSBT0->GetX()>pointSBT1->GetX()){
590             pointSBT0=myDetElem->GetLButtonTarget(1);
591             pointSBT1=myDetElem->GetLButtonTarget(0);
592           }
593           Double_t lX = pointSBT0->GetX();
594           while(lX<pointSBT1->GetX()){
595             Double_t lY = pointSBT0->GetY()-200;
596             while(lY<pointSBT0->GetY()+200){
597               pl[0] = lX;  pl[1] = lY;  pl[2] = 0.;
598               (TGeoCombiTrans((*(myDetElem->GetBaseTrf()))*(*(myDetElem->GetAlignTrf())))).LocalToMaster(pl,pg);
599               if(myDetElem->GetGButtonTarget(0)->GetPointName().Contains("A")){
600                 hSSTap->Fill(pg[0],pg[1],-pg[2]);
601               }
602               else {
603                 hSSTcp->Fill(pg[0],pg[1],-pg[2]);
604               }
605               lY+=hSSTap->GetYaxis()->GetBinWidth(1);
606             }
607             lX+=hSSTap->GetXaxis()->GetBinWidth(1);
608           }
609         }
610       }
611     }
612
613     if( saveps ){
614       ps->NewPage();
615     }
616
617     // View from side A
618     str = Form("Chamber %d - side A",chId+1);
619     TCanvas *cvn0 = new TCanvas("cvn0",str,cWidth,cHeight);
620     canvas = cvn0;
621     canvas->Range(0,0,21,29);
622   
623     title = Form(" Deformations of chamber %d - View from side A ",chId+1);
624     TPaveLabel *theTitle0 = new TPaveLabel(3,27.0,18,28.5,title,"br");
625     theTitle = theTitle0;
626     theTitle->SetFillColor(18);
627     theTitle->SetTextFont(32);
628     theTitle->SetTextSize(0.4);
629     theTitle->SetTextColor(1);
630     theTitle->Draw();
631  
632     TPad *pad0 = new TPad("pad0","pad0",0.01,0.01,0.98,0.91,0);
633     pad = pad0;
634     pad->Draw();
635     pad->Divide(2,2);
636
637     pad->cd(1);
638     hCPSTa->SetMinimum(9500);
639     hCPSTa->SetMaximum(9750);
640     hCPSTa->Draw("lego2z");
641
642     pad->cd(2);
643     hSSTa->SetMinimum(9500);
644     hSSTa->SetMaximum(9800);
645     hSSTa->Draw("lego2z");
646
647     pad->cd(3);
648
649     pad->cd(4);
650     hSSTap->SetMinimum(9500);
651     hSSTap->SetMaximum(9800);
652     hSSTap->Draw("surf2z");
653
654     pad->Update();
655     if(saveps){
656       ps->NewPage();
657     }
658
659     // Inv Mass, Multiplicity
660     str = Form("Chamber %d - side C",chId+1);
661     TCanvas *cvn1 = new TCanvas("cvn1",str,cWidth,cHeight);
662     canvas = cvn1;
663     canvas->Range(0,0,21,29);
664
665     title = Form(" Deformations of chamber %d - View from side C ",chId+1);
666     TPaveLabel *theTitle1 = new TPaveLabel(3,27.0,18,28.5,title,"br");
667     theTitle = theTitle1;
668     theTitle->SetFillColor(18);
669     theTitle->SetTextFont(32);
670     theTitle->SetTextSize(0.4);
671     theTitle->SetTextColor(1);
672     theTitle->Draw();
673  
674     TPad *pad1 = new TPad("pad1","pad1",0.01,0.01,0.98,0.91,0);
675     pad = pad1;
676     pad->Draw();
677     pad->Divide(2,2);
678
679     pad->cd(1);
680     hCPSTc->SetMinimum(9550);
681     hCPSTc->SetMaximum(9850);
682     hCPSTc->Draw("lego2z");
683
684     pad->cd(2);
685     hSSTc->SetMinimum(6930);
686     hSSTc->SetMaximum(7030);
687     hSSTc->Draw("lego2z");
688
689     pad->cd(3);
690
691     pad->cd(4);
692     hSSTcp->SetMinimum(9600);
693     hSSTcp->SetMaximum(9850);
694     hSSTcp->Draw("surf2z");
695   }
696
697   pad->Update();
698   if( saveps ){
699     ps->Close();
700   }
701
702   TFile *hFile = new TFile(Form("spCH%d_PL2G_ISF.root",chId+1),"RECREATE");
703   gDeltaDiffCenXDetElem0->Write("gDeltaDiffCenXDetElem0");
704   gDeltaDiffCenYDetElem0->Write("gDeltaDiffCenYDetElem0");
705   gDeltaDiffCenZDetElem0->Write("gDeltaDiffCenZDetElem0");
706   gDeltaDiffPsiDetElem0->Write("gDeltaDiffPsiDetElem0");
707   gDeltaDiffThtDetElem0->Write("gDeltaDiffThtDetElem0");
708   gDeltaDiffPhiDetElem0->Write("gDeltaDiffPhiDetElem0");
709   myDetElemDeltaDiffCenX->Write();
710   myDetElemDeltaDiffCenY->Write();
711   myDetElemDeltaDiffCenZ->Write();
712   myDetElemDeltaDiffRotX->Write();
713   myDetElemDeltaDiffRotY->Write();
714   myDetElemDeltaDiffRotZ->Write();
715   hFile->Close();
716   hFile->Delete();
717 }