]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/MUONSurveyCh4.C
Possibility to vary the cos(PA) cut
[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 MUONSurveyCh4.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.Print();
209     myChamber->SetLocalTransformation(new TGeoCombiTrans(trfThChamber),kTRUE);
210
211     // Set Chamber plane function
212     cout << "Setting plane for Chamber" << iCh+1 << " ..." << endl;
213     myChamber->SetPlane(Form("fChamber%d",iCh+1));
214     myChamber->SetPlaneParameters(0.,0.,0.);
215
216     // Fit a plane to sticker targets
217     Double_t lCChi2 = myChamber->FitPlane();
218     cout << "... chi2 = " << lCChi2 << " ..." << endl; 
219     
220     // Get best transformation from fitting method 
221     // (local to global transformation)
222     cout << "Trying fitting method for chamber " << iCh << endl;
223     myChamber->SurveyToAlign();    
224     //    myChamber->SurveyToAlign(myChamber->GetPlane()->GetParameter(0),myChamber->GetPlane()->GetParameter(1),myChamber->GetPlane()->GetParError(0),myChamber->GetPlane()->GetParError(1));    
225     
226     for (int iDetElemI=0; iDetElemI<nDetElemsI; iDetElemI++){
227       myDetElem =  myChamber->GetDetElem(iDetElemI);
228       Int_t iDetElem = iDetElemsOfChamber[iCh][iDetElemI];
229
230       trfThDetElem.Clear();
231       trfThDetElem = TGeoCombiTrans(*transform->GetDetElement(iDetElemToDetElemId[iDetElem])->GetLocalTransformation());
232       trfThDetElem.Print();
233       myDetElem->SetLocalTransformation(new TGeoCombiTrans(trfThDetElem),kTRUE);
234
235       for (int iCor=0; iCor<3; iCor++){
236         lCenDetElem[iDetElem][iCor] = 0;
237         lRotDetElem[iDetElem][iCor] = 0;
238       }
239
240       if (bMonitor){
241         // MONITOR: Draw graph with sticker targets for plane fitting
242         myDetElem->DrawSTargets();
243         gPad->Update();
244       }
245
246       // Get DetElem transformation. 
247       // Psi and Theta are obtained by fitting a plane to the sticker targets.
248       // Then Xc, Yc, Zc and Phi are obtained by solving the equations to the ref.
249       // syst. transformation of the button targets
250
251       // Set DetElem plane function
252       cout << "Setting plane for DetElem" << iDetElem+1 << " ..." << endl;
253       myDetElem->SetPlane(Form("fDetElem%d",iDetElem+1));
254       myDetElem->SetPlaneParameters(0.,0.,3.)
255 ;
256       // Fit a plane to sticker targets
257       Double_t lChi2 = myDetElem->FitPlane();
258       cout << "... chi2 = " << lChi2 << " ..." << endl; 
259       
260       lRotDetElem[iDetElem][0]=TMath::ATan(myDetElem->GetPlane()->GetParameter(0));
261       lRotDetElem[iDetElem][1]=TMath::ATan(myDetElem->GetPlane()->GetParameter(1));
262       
263       // Calculate Mean transformation using previous plane fit 
264       // and pairs of button targets
265       myDetElem->CalculateMeanTransf(lCenDetElem[iDetElem],lRotDetElem[iDetElem]);
266     
267       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;          
268       
269      
270       // Get best transformation from fitting method 
271       // (local to global transformation)
272       cout << "Trying fitting method for DetElemId " << iDetElemToDetElemId[iDetElem] << endl;
273       //      myDetElem->SurveyToAlign();     
274       myDetElem->SurveyToAlign(lRotDetElem[iDetElem][0],lRotDetElem[iDetElem][1],myDetElem->GetPlane()->GetParError(0),myDetElem->GetPlane()->GetParError(1));    
275     }
276   }
277
278   // Print found transformation of Detection Element (plane fit + loop)        
279   for (int iCh =0; iCh<=1; iCh++) {
280     for (int iDetElemI=0; iDetElemI<nDetElemsI; iDetElemI++){
281       Int_t iDetElem = iDetElemsOfChamber[iCh][iDetElemI];
282       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;          
283     }
284   }
285
286   // Print Theoretical transformation of Detection Element
287   for (int iCh =0; iCh<=1; iCh++) {
288     myChamber = (AliMUONSurveyChamber*)myChamberArray->At(iCh);
289     for (int iDetElemI=0; iDetElemI<nDetElemsI; iDetElemI++){
290       myChamber->GetDetElem(iDetElemI)->PrintLocalTrf();
291     }
292   }
293
294   // Print found delta transformation of Detection Element
295   for (int iCh =0; iCh<=1; iCh++) {
296     myChamber = (AliMUONSurveyChamber*)myChamberArray->At(iCh);      
297     for (int iDetElemI=0; iDetElemI<nDetElemsI; iDetElemI++){
298       myChamber->GetDetElem(iDetElemI)->PrintAlignTrf();
299     }
300   }
301
302   //
303   // Compare transformations to expected ones 
304   //
305   Int_t iDetElemToPos[18] = {0, 1, 2, 3};
306  
307   TGraphErrors *gDeltaDiffCenXDetElem0 = new TGraphErrors(nDetElems);
308   TGraphErrors *gDeltaDiffCenYDetElem0 = new TGraphErrors(nDetElems);
309   TGraphErrors *gDeltaDiffCenZDetElem0 = new TGraphErrors(nDetElems);
310   TGraphErrors *gDeltaDiffPsiDetElem0 = new TGraphErrors(nDetElems);
311   TGraphErrors *gDeltaDiffThtDetElem0 = new TGraphErrors(nDetElems);
312   TGraphErrors *gDeltaDiffPhiDetElem0 = new TGraphErrors(nDetElems);
313
314   for (int iCh =0; iCh<=1; iCh++) {
315     myChamber = (AliMUONSurveyChamber*)myChamberArray->At(iCh);      
316     for (int iDetElemI=0; iDetElemI<nDetElemsI; iDetElemI++){
317       myChamber->GetDetElem(iDetElemI)->GetAlignTrf()->Print();
318     }
319   }
320
321   for (int iCh =0; iCh<=1; iCh++) {
322     myChamber = (AliMUONSurveyChamber*)myChamberArray->At(iCh);
323     // Store delta transformations to use for CDB entry creation
324     dtrfDetElem[nDetElems+iCh].Clear();
325     dtrfDetElem[nDetElems+iCh] = *(myChamber->GetAlignTrf());
326
327     // Those are for checks and visualizations
328     localTrfDetElem[nDetElems+iCh].Clear();
329     localTrfDetElem[nDetElems+iCh] = (*(myChamber->GetLocalTrf())*(*(myChamber->GetAlignTrf())));
330     localTrfDetElem[nDetElems+iCh].Print();
331     lDiffCenDetElem0[nDetElems+iCh] = (Double_t*)localTrfDetElem[nDetElems+iCh].GetTranslation();
332     AliMUONSurveyUtil::MatrixToAngles(localTrfDetElem[nDetElems+iCh].GetRotationMatrix(),lDiffRotDetElem0[nDetElems+iCh]);
333     
334     localTrfThDetElem[nDetElems+iCh].Clear();
335     localTrfThDetElem[nDetElems+iCh] = (*(myChamber->GetLocalTrf()));
336     localTrfThDetElem[nDetElems+iCh].Print();
337     lDiffThCenDetElem0[nDetElems+iCh] = (Double_t*)localTrfThDetElem[nDetElems+iCh].GetTranslation();
338     AliMUONSurveyUtil::MatrixToAngles(localTrfThDetElem[nDetElems+iCh].GetRotationMatrix(),lDiffThRotDetElem0[nDetElems+iCh]);
339
340     for (int iCor=0; iCor<3; iCor++){
341       lDeltaDiffCenDetElem0[nDetElems+iCh][iCor] = lDiffCenDetElem0[nDetElems+iCh][iCor]-lDiffThCenDetElem0[nDetElems+iCh][iCor];
342       lDeltaDiffRotDetElem0[nDetElems+iCh][iCor] = lDiffRotDetElem0[nDetElems+iCh][iCor]-lDiffThRotDetElem0[nDetElems+iCh][iCor];
343       if (lDeltaDiffRotDetElem0[nDetElems+iCh][iCor]>TMath::Pi()) lDeltaDiffRotDetElem0[nDetElems+iCh][iCor]-=TMath::TwoPi();
344       if (lDeltaDiffRotDetElem0[nDetElems+iCh][iCor]<-TMath::Pi()) lDeltaDiffRotDetElem0[nDetElems+iCh][iCor]+=TMath::TwoPi();
345     }      
346
347     for (int iDetElemI=0; iDetElemI<nDetElemsI; iDetElemI++){
348       Int_t iDetElem = iDetElemsOfChamber[iCh][iDetElemI];
349       myDetElem =  myChamber->GetDetElem(iDetElemI);
350       // Store delta transformations to use for CDB entry creation
351       dtrfDetElem[iDetElem].Clear();
352       dtrfDetElem[iDetElem] = *(myDetElem->GetAlignTrf());
353
354       // Those are for checks and visualizations
355       localTrfDetElem[iDetElem].Clear();
356       localTrfDetElem[iDetElem] = (*(myDetElem->GetLocalTrf())*(*(myDetElem->GetAlignTrf())));
357       //      localTrfDetElem[iDetElem] = (*(myDetElem->GetBaseTrf())*(*(myDetElem->GetAlignTrf())));
358       localTrfDetElem[iDetElem].Print();
359       lDiffCenDetElem0[iDetElem] = (Double_t*)localTrfDetElem[iDetElem].GetTranslation();
360       AliMUONSurveyUtil::MatrixToAngles(localTrfDetElem[iDetElem].GetRotationMatrix(),lDiffRotDetElem0[iDetElem]);
361 //       lDiffCenDetElem0[iDetElem] = lCenDetElem[iDetElem];
362 //       lDiffRotDetElem0[iDetElem] = lRotDetElem[iDetElem];
363
364       localTrfThDetElem[iDetElem].Clear();
365       localTrfThDetElem[iDetElem] = (*(myDetElem->GetLocalTrf()));
366       //      localTrfThDetElem[iDetElem] = (*(myDetElem->GetBaseTrf()));
367       localTrfThDetElem[iDetElem].Print();
368       lDiffThCenDetElem0[iDetElem] = (Double_t*)localTrfThDetElem[iDetElem].GetTranslation();
369       AliMUONSurveyUtil::MatrixToAngles(localTrfThDetElem[iDetElem].GetRotationMatrix(),lDiffThRotDetElem0[iDetElem]);
370       
371       for (int iCor=0; iCor<3; iCor++){
372         lDeltaDiffCenDetElem0[iDetElem][iCor] = lDiffCenDetElem0[iDetElem][iCor]-lDiffThCenDetElem0[iDetElem][iCor];
373         lDeltaDiffRotDetElem0[iDetElem][iCor] = lDiffRotDetElem0[iDetElem][iCor]-lDiffThRotDetElem0[iDetElem][iCor];
374         if (lDeltaDiffRotDetElem0[iDetElem][iCor]>TMath::Pi()) lDeltaDiffRotDetElem0[iDetElem][iCor]-=TMath::TwoPi();
375         if (lDeltaDiffRotDetElem0[iDetElem][iCor]<-TMath::Pi()) lDeltaDiffRotDetElem0[iDetElem][iCor]+=TMath::TwoPi();
376       }
377
378       gDeltaDiffCenXDetElem0->SetPoint(iDetElem,1e1*lDeltaDiffCenDetElem0[iDetElem][0],iDetElemToPos[iDetElem]+1);
379       gDeltaDiffCenYDetElem0->SetPoint(iDetElem,1e1*lDeltaDiffCenDetElem0[iDetElem][1],iDetElemToPos[iDetElem]+1);
380       gDeltaDiffCenZDetElem0->SetPoint(iDetElem,1e1*lDeltaDiffCenDetElem0[iDetElem][2],iDetElemToPos[iDetElem]+1);
381       gDeltaDiffPsiDetElem0->SetPoint(iDetElem,1e3*lDeltaDiffRotDetElem0[iDetElem][0],iDetElemToPos[iDetElem]+1);
382       gDeltaDiffThtDetElem0->SetPoint(iDetElem,1e3*lDeltaDiffRotDetElem0[iDetElem][1],iDetElemToPos[iDetElem]+1);
383       gDeltaDiffPhiDetElem0->SetPoint(iDetElem,1e3*lDeltaDiffRotDetElem0[iDetElem][2],iDetElemToPos[iDetElem]+1);
384       gDeltaDiffCenXDetElem0->SetPointError(iDetElem,1e1*myDetElem->GetFitter()->GetParError(0),0.);
385       gDeltaDiffCenYDetElem0->SetPointError(iDetElem,1e1*myDetElem->GetFitter()->GetParError(1),0.);
386       gDeltaDiffCenZDetElem0->SetPointError(iDetElem,1e1*myDetElem->GetFitter()->GetParError(2),0.);
387       gDeltaDiffPsiDetElem0->SetPointError(iDetElem,1e3*myDetElem->GetFitter()->GetParError(3),0.);
388       gDeltaDiffThtDetElem0->SetPointError(iDetElem,1e3*myDetElem->GetFitter()->GetParError(4),0.);
389       gDeltaDiffPhiDetElem0->SetPointError(iDetElem,1e3*myDetElem->GetFitter()->GetParError(5),0.);  
390     }
391   }
392
393   // Apply the found misalignments to the geometry
394   AliMUONGeometryTransformer *newTransform = AliMUONSurveyUtil::ReAlign(transform,chId,nDetElems,iDetElemPseudoIdToDetElem,dtrfDetElem,true); 
395   newTransform->WriteTransformations(Form("transform2ReAlignSurveyCh%d.dat",chId+1));
396
397   if(bOCDB){
398     // Generate realigned data in local cdb
399     const TClonesArray* array = newTransform->GetMisAlignmentData();
400     
401     // Set the alignment resolution in the align objects for this chamber
402     Double_t chResX = myChamberI->GetAlignResX();
403     Double_t chResY = myChamberI->GetAlignResY();
404     Double_t deResX = (myChamberI->GetMeanDetElemAlignResX()+myChamberO->GetMeanDetElemAlignResX())/2.;
405     Double_t deResY = (myChamberI->GetMeanDetElemAlignResY()+myChamberO->GetMeanDetElemAlignResY())/2.;
406     printf("Chamber alignment resolutions: resX=%f , resY=%f\n",chResX,chResY); 
407     printf("Detection Elements alignment resolutions: resX=%f , resY=%f\n",deResX,deResY); 
408     chResX = TMath::Sqrt(0.1*0.1+chResX*chResX);
409     chResY = TMath::Sqrt(0.1*0.1+chResY*chResY);
410     AliMUONSurveyUtil::SetAlignmentResolution(array,chId,chResX,chResY,deResX,deResY);
411     
412     // CDB manager
413     AliCDBManager* cdbManager = AliCDBManager::Instance();
414     cdbManager->SetDefaultStorage(Form("local://ReAlignSurveyCh%dCDB",chId+1));
415     
416     AliCDBMetaData* cdbData = new AliCDBMetaData();
417     cdbData->SetResponsible("Dimuon Offline project");
418     cdbData->SetComment("MUON alignment objects with survey realignment");
419     AliCDBId id("MUON/Align/Data", 0, AliCDBRunRange::Infinity()); 
420     cdbManager->Put(const_cast<TClonesArray*>(array), id, cdbData);
421   }
422
423   for(Int_t iCor=0; iCor<3; iCor++){
424     for(Int_t iDetElem=0; iDetElem<nDetElems; iDetElem++){
425       cout << lDeltaDiffCenDetElem0[iDetElem][iCor] << " " << lDiffCenDetElem0[iDetElem][iCor] << " " << lDiffThCenDetElem0[iDetElem][iCor] << endl;
426     }
427     cout << endl;
428   }
429   cout << endl;
430   for(Int_t iCor=0; iCor<3; iCor++){
431     for(Int_t iDetElem=0; iDetElem<nDetElems; iDetElem++){
432       cout << lDeltaDiffRotDetElem0[iDetElem][iCor] << " " << lDiffRotDetElem0[iDetElem][iCor] << " " << lDiffThRotDetElem0[iDetElem][iCor] << endl;
433     }
434     cout << endl;
435   }
436
437   TH1F *myDetElemDeltaDiffCenX = new TH1F("myDetElemDeltaDiffCenX","myDetElemDeltaDiffCenX",100,-10,10);
438   myDetElemDeltaDiffCenX->SetMaximum(nDetElems+1);
439   myDetElemDeltaDiffCenX->SetMinimum(0);
440   TH1F *myDetElemDeltaDiffCenY = new TH1F("myDetElemDeltaDiffCenY","myDetElemDeltaDiffCenY",100,-10,10);
441   myDetElemDeltaDiffCenY->SetMaximum(nDetElems+1);
442   myDetElemDeltaDiffCenY->SetMinimum(0);
443   TH1F *myDetElemDeltaDiffCenZ = new TH1F("myDetElemDeltaDiffCenZ","myDetElemDeltaDiffCenZ",100,-20,20);
444   myDetElemDeltaDiffCenZ->SetMaximum(nDetElems+1);
445   myDetElemDeltaDiffCenZ->SetMinimum(0);
446
447   TH1F *myDetElemDeltaDiffRotX = new TH1F("myDetElemDeltaDiffRotX","myDetElemDeltaDiffRotX",100,-15,15);
448   myDetElemDeltaDiffRotX->SetMaximum(nDetElems+1);
449   myDetElemDeltaDiffRotX->SetMinimum(0);
450   TH1F *myDetElemDeltaDiffRotY = new TH1F("myDetElemDeltaDiffRotY","myDetElemDeltaDiffRotY",100,-15,15);
451   myDetElemDeltaDiffRotY->SetMaximum(nDetElems+1);
452   myDetElemDeltaDiffRotY->SetMinimum(0);
453   TH1F *myDetElemDeltaDiffRotZ = new TH1F("myDetElemDeltaDiffRotZ","myDetElemDeltaDiffRotZ",100,-5,5);
454   myDetElemDeltaDiffRotZ->SetMaximum(nDetElems+1);
455   myDetElemDeltaDiffRotZ->SetMinimum(0);
456
457   //
458   // ******** Starting plots 
459   //
460   TCanvas *canvas;
461   TPad *pad;
462   TPaveLabel *theTitle;
463   gStyle->SetPalette(1);
464
465   TPostScript *ps = 0;
466   
467   if( saveps ){
468     ps = new TPostScript(Form("surveyChamber%d",chId+1),filetype); 
469     ps->NewPage();
470   }
471
472   // Observed misalignments
473   str = Form("Chamber %d",chId+1);
474   TCanvas *cvn2 = new TCanvas("cvn2",str,cWidth,cHeight);
475   canvas = cvn2;
476   canvas->Range(0,0,21,29);
477   
478   title = Form(" MisAlignments Chamber %d - PL2G - In Single Frame ",chId+1);
479   TPaveLabel *theTitle2 = new TPaveLabel(3,27.0,18,28.5,title,"br");
480   theTitle = theTitle2;
481   theTitle->SetFillColor(18);
482   theTitle->SetTextFont(32);
483   theTitle->SetTextSize(0.4);
484   theTitle->SetTextColor(1);
485   theTitle->Draw();
486  
487   TPad *pad2 = new TPad("pad2","pad2",0.01,0.01,0.98,0.91,0);
488   pad = pad2;
489   pad->Draw();
490   TLine *ch0Line = new TLine(0,1,0,2);
491   TLine *ch1Line = new TLine(0,1,0,2);
492   ch1Line->SetLineStyle(2);
493   pad->Divide(3,2);
494   
495   pad->cd(1);
496   myDetElemDeltaDiffCenX->Draw();
497   myDetElemDeltaDiffCenX->SetXTitle("#Delta[xc_{i}^{m}-xc_{i}^{th}] (mm)");
498   myDetElemDeltaDiffCenX->SetYTitle("DetElem arbitrary ordering");
499   gDeltaDiffCenXDetElem0->SetMarkerStyle(20);
500   gDeltaDiffCenXDetElem0->Draw("P");
501   ch0Line->DrawLine(1e1*lDeltaDiffCenDetElem0[nDetElems+0][0],0.5,1e1*lDeltaDiffCenDetElem0[nDetElems+0][0],1.5);
502   ch0Line->DrawLine(1e1*lDeltaDiffCenDetElem0[nDetElems+0][0],3.5,1e1*lDeltaDiffCenDetElem0[nDetElems+0][0],4.5);
503   ch1Line->DrawLine(1e1*lDeltaDiffCenDetElem0[nDetElems+1][0],1.5,1e1*lDeltaDiffCenDetElem0[nDetElems+1][0],3.5);
504
505   pad->cd(2);
506   myDetElemDeltaDiffCenY->Draw();
507   myDetElemDeltaDiffCenY->SetXTitle("#Delta[yc_{i}^{m}-yc_{i}^{th}] (mm)");
508   myDetElemDeltaDiffCenY->SetYTitle("DetElem arbitrary ordering");
509   gDeltaDiffCenYDetElem0->SetMarkerStyle(20);
510   gDeltaDiffCenYDetElem0->Draw("P");
511   ch0Line->DrawLine(1e1*lDeltaDiffCenDetElem0[nDetElems+0][1],0.5,1e1*lDeltaDiffCenDetElem0[nDetElems+0][1],1.5);
512   ch0Line->DrawLine(1e1*lDeltaDiffCenDetElem0[nDetElems+0][1],3.5,1e1*lDeltaDiffCenDetElem0[nDetElems+0][1],4.5);
513   ch1Line->DrawLine(1e1*lDeltaDiffCenDetElem0[nDetElems+1][1],1.5,1e1*lDeltaDiffCenDetElem0[nDetElems+1][1],3.5);
514
515   pad->cd(3);
516   myDetElemDeltaDiffCenZ->Draw();
517   myDetElemDeltaDiffCenZ->SetXTitle("#Delta[zc_{i}^{m}-zc_{i}^{th}] (mm)");
518   myDetElemDeltaDiffCenZ->SetYTitle("DetElem arbitrary ordering");
519   gDeltaDiffCenZDetElem0->SetMarkerStyle(20);
520   gDeltaDiffCenZDetElem0->Draw("P");
521   ch0Line->DrawLine(1e1*lDeltaDiffCenDetElem0[nDetElems+0][2],0.5,1e1*lDeltaDiffCenDetElem0[nDetElems+0][2],1.5);
522   ch0Line->DrawLine(1e1*lDeltaDiffCenDetElem0[nDetElems+0][2],3.5,1e1*lDeltaDiffCenDetElem0[nDetElems+0][2],4.5);
523   ch1Line->DrawLine(1e1*lDeltaDiffCenDetElem0[nDetElems+1][2],1.5,1e1*lDeltaDiffCenDetElem0[nDetElems+1][2],3.5);
524
525   pad->cd(4);
526   myDetElemDeltaDiffRotX->Draw();
527   myDetElemDeltaDiffRotX->SetXTitle("#Delta[#psi_{i}^{m}-#psi_{i}^{th}] (mrad)");
528   myDetElemDeltaDiffRotX->SetYTitle("DetElem arbitrary ordering");
529   gDeltaDiffPsiDetElem0->SetMarkerStyle(20);
530   gDeltaDiffPsiDetElem0->Draw("P");
531   ch0Line->DrawLine(1e3*lDeltaDiffRotDetElem0[nDetElems+0][0],0.5,1e3*lDeltaDiffRotDetElem0[nDetElems+0][0],1.5);
532   ch0Line->DrawLine(1e3*lDeltaDiffRotDetElem0[nDetElems+0][0],3.5,1e3*lDeltaDiffRotDetElem0[nDetElems+0][0],4.5);
533   ch1Line->DrawLine(1e3*lDeltaDiffRotDetElem0[nDetElems+1][0],1.5,1e3*lDeltaDiffRotDetElem0[nDetElems+1][0],3.5);
534
535   pad->cd(5);
536   myDetElemDeltaDiffRotY->Draw();
537   myDetElemDeltaDiffRotY->SetXTitle("#Delta[#theta_{i}^{m}-#theta_{i}^{th}] (mrad)");
538   myDetElemDeltaDiffRotY->SetYTitle("DetElem arbitrary ordering");
539   gDeltaDiffThtDetElem0->SetMarkerStyle(20);
540   gDeltaDiffThtDetElem0->Draw("P");
541   ch0Line->DrawLine(1e3*lDeltaDiffRotDetElem0[nDetElems+0][1],0.5,1e3*lDeltaDiffRotDetElem0[nDetElems+0][1],1.5);
542   ch0Line->DrawLine(1e3*lDeltaDiffRotDetElem0[nDetElems+0][1],3.5,1e3*lDeltaDiffRotDetElem0[nDetElems+0][1],4.5);
543   ch1Line->DrawLine(1e3*lDeltaDiffRotDetElem0[nDetElems+1][1],1.5,1e3*lDeltaDiffRotDetElem0[nDetElems+1][1],3.5);
544
545   pad->cd(6);
546   myDetElemDeltaDiffRotZ->Draw();
547   myDetElemDeltaDiffRotZ->SetXTitle("#Delta[#phi_{i}^{m}-#phi_{i}^{th}] (mrad)");
548   myDetElemDeltaDiffRotZ->SetYTitle("DetElem arbitrary ordering");
549   gDeltaDiffPhiDetElem0->SetMarkerStyle(20);
550   gDeltaDiffPhiDetElem0->Draw("P");
551   ch0Line->DrawLine(1e3*lDeltaDiffRotDetElem0[nDetElems+0][2],0.5,1e3*lDeltaDiffRotDetElem0[nDetElems+0][2],1.5);
552   ch0Line->DrawLine(1e3*lDeltaDiffRotDetElem0[nDetElems+0][2],3.5,1e3*lDeltaDiffRotDetElem0[nDetElems+0][2],4.5);
553   ch1Line->DrawLine(1e3*lDeltaDiffRotDetElem0[nDetElems+1][2],1.5,1e3*lDeltaDiffRotDetElem0[nDetElems+1][2],3.5);
554
555   pad->Update();
556  
557   if(bMonitor){    
558     // MONITOR: Histograms for monitoring
559     TH2F *hCPSTa = new TH2F("hCPSTa","hCPSTa",100,-200,200,100,-200,200);
560     TH2F *hCPSTc = new TH2F("hCPSTc","hCPSTc",100,-200,200,100,-200,200);
561     TH2F *hSSTa = new TH2F("hSSTa","hSSTa",100,-200,200,100,-200,200);
562     TH2F *hSSTc = new TH2F("hSSTc","hSSTc",100,-200,200,100,-200,200);
563     TH2F *hSSTap = new TH2F("hSSTap","hSSTap",800,-200,200,800,-200,200);
564     TH2F *hSSTcp = new TH2F("hSSTcp","hSSTcp",800,-200,200,800,-200,200);
565     
566     // MONITOR: Fill histograms with chambers and slats sticker target positions
567     for (int iCh =0; iCh<=1; iCh++) {
568       myChamber = (AliMUONSurveyChamber*)myChamberArray->At(iCh);            
569       cout << "Filling histograms of sticker target for chamber" << iCh+1 << " ..." << endl;
570       //    myChamber->FillCPSTHistograms(TString("C"),hCPSTc);
571       myChamber->FillDESTHistograms(TString("T1_46"),hSSTc);
572
573       for (int iDetElemI=0; iDetElemI<nDetElemsI; iDetElemI++){
574         //      Int_t iDetElem = iDetElemsOfChamber[iCh][iDetElemI];
575         myDetElem =  myChamber->GetDetElem(iDetElemI);
576         // MONITOR: Fill slat plane for fit monitor.
577         Double_t pl[3] = {0};
578         Double_t pg[3] = {0};
579         AliSurveyPoint *pointSBT0 = myDetElem->GetLButtonTarget(0);
580         AliSurveyPoint *pointSBT1 = myDetElem->GetLButtonTarget(2);
581         if(pointSBT0&&pointSBT1) {
582           if (pointSBT0->GetX()>pointSBT1->GetX()){
583             pointSBT0=myDetElem->GetLButtonTarget(2);
584             pointSBT1=myDetElem->GetLButtonTarget(0);
585           }
586           Double_t lX = pointSBT0->GetX();
587           while(lX<pointSBT1->GetX()){
588             Double_t lY = pointSBT0->GetY()-50;
589             while(lY<pointSBT0->GetY()+50){
590               pl[0] = lX;  pl[1] = lY;  pl[2] = 0.;
591               (TGeoCombiTrans((*(myDetElem->GetBaseTrf()))*(*(myDetElem->GetAlignTrf())))).LocalToMaster(pl,pg);
592               if(myDetElem->GetGButtonTarget(0)->GetPointName().Contains("A")){
593                 if (hSSTap->GetBinContent(hSSTap->FindBin(pg[0],pg[1]))==0)
594                   hSSTap->Fill(pg[0],pg[1],-pg[2]);
595               }
596               else {
597                 if (hSSTcp->GetBinContent(hSSTcp->FindBin(pg[0],pg[1]))==0)
598                   hSSTcp->Fill(pg[0],pg[1],-pg[2]);
599               }
600               lY+=hSSTap->GetYaxis()->GetBinWidth(1);
601             }
602             lX+=hSSTap->GetXaxis()->GetBinWidth(1);
603           }
604         }
605       }
606     }
607
608     if( saveps ){
609       ps->NewPage();
610     }
611
612     // View from side A
613     str = Form("Chamber %d - side A",chId+1);
614     TCanvas *cvn0 = new TCanvas("cvn0",str,cWidth,cHeight);
615     canvas = cvn0;
616     canvas->Range(0,0,21,29);
617   
618     title = Form(" Deformations of chamber %d - View from side A ",chId+1);
619     TPaveLabel *theTitle0 = new TPaveLabel(3,27.0,18,28.5,title,"br");
620     theTitle = theTitle0;
621     theTitle->SetFillColor(18);
622     theTitle->SetTextFont(32);
623     theTitle->SetTextSize(0.4);
624     theTitle->SetTextColor(1);
625     theTitle->Draw();
626  
627     TPad *pad0 = new TPad("pad0","pad0",0.01,0.01,0.98,0.91,0);
628     pad = pad0;
629     pad->Draw();
630     pad->Divide(2,2);
631
632     pad->cd(1);
633     hCPSTa->SetMinimum(950);
634     hCPSTa->SetMaximum(975);
635     hCPSTa->Draw("lego2z");
636
637     pad->cd(2);
638     hSSTa->SetMinimum(950);
639     hSSTa->SetMaximum(980);
640     hSSTa->Draw("lego2z");
641
642     pad->cd(3);
643
644     pad->cd(4);
645     hSSTap->SetMinimum(950);
646     hSSTap->SetMaximum(980);
647     hSSTap->Draw("surf2z");
648
649     pad->Update();
650     if(saveps){
651       ps->NewPage();
652     }
653
654     // Inv Mass, Multiplicity
655     str = Form("Chamber %d - side C",chId+1);
656     TCanvas *cvn1 = new TCanvas("cvn1",str,cWidth,cHeight);
657     canvas = cvn1;
658     canvas->Range(0,0,21,29);
659
660     title = Form(" Deformations of chamber %d - View from side C ",chId+1);
661     TPaveLabel *theTitle1 = new TPaveLabel(3,27.0,18,28.5,title,"br");
662     theTitle = theTitle1;
663     theTitle->SetFillColor(18);
664     theTitle->SetTextFont(32);
665     theTitle->SetTextSize(0.4);
666     theTitle->SetTextColor(1);
667     theTitle->Draw();
668  
669     TPad *pad1 = new TPad("pad1","pad1",0.01,0.01,0.98,0.91,0);
670     pad = pad1;
671     pad->Draw();
672     pad->Divide(2,2);
673
674     Double_t lMin, lMax;
675
676     pad->cd(1);
677 //     hCPSTc->SetMinimum(955);
678 //     hCPSTc->SetMaximum(985);
679 //     hCPSTc->Draw("lego2z");
680
681     pad->cd(2);
682     lMin = hSSTc->GetMinimum(0);
683     hSSTc->SetMinimum(TMath::Floor(lMin));
684     lMax = hSSTc->GetMaximum();
685     hSSTc->SetMaximum(TMath::Ceil(lMax));
686     hSSTc->Draw("lego2z");
687
688     pad->cd(3);
689
690     pad->cd(4);
691     lMin = hSSTcp->GetMinimum(0);
692     hSSTcp->SetMinimum(TMath::Floor(lMin));
693     lMax = hSSTcp->GetMaximum();
694     hSSTcp->SetMaximum(TMath::Ceil(lMax));
695     hSSTcp->Draw("surf2z");
696   }
697
698   pad->Update();
699   if( saveps ){
700     ps->Close();
701   }
702
703   TFile *hFile = new TFile(Form("spCH%d_PL2G_ISF.root",chId+1),"RECREATE");
704   gDeltaDiffCenXDetElem0->Write("gDeltaDiffCenXDetElem0");
705   gDeltaDiffCenYDetElem0->Write("gDeltaDiffCenYDetElem0");
706   gDeltaDiffCenZDetElem0->Write("gDeltaDiffCenZDetElem0");
707   gDeltaDiffPsiDetElem0->Write("gDeltaDiffPsiDetElem0");
708   gDeltaDiffThtDetElem0->Write("gDeltaDiffThtDetElem0");
709   gDeltaDiffPhiDetElem0->Write("gDeltaDiffPhiDetElem0");
710   myDetElemDeltaDiffCenX->Write();
711   myDetElemDeltaDiffCenY->Write();
712   myDetElemDeltaDiffCenZ->Write();
713   myDetElemDeltaDiffRotX->Write();
714   myDetElemDeltaDiffRotY->Write();
715   myDetElemDeltaDiffRotZ->Write();
716   hFile->Close();
717   hFile->Delete();
718 }