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