]> git.uio.no Git - u/mrichter/AliRoot.git/blame_incremental - TPC/CalibMacros/RegisterCorrection.C
Update master to aliroot
[u/mrichter/AliRoot.git] / TPC / CalibMacros / RegisterCorrection.C
... / ...
CommitLineData
1/// \file RegisterCorrection.C
2///
3/// \author marian.ivanov@cern.ch
4///
5/// Register primitive corrections: base functions for minimization:
6/// Id numbers are associated to given primitive corrections.
7/// See comments in the function headers.
8/// Used only for residuals minimization not in the reconstruction.
9/// File with all primitives expected to be in the current directory:
10/// filenames:
11/// TPCCorrectionPrimitives.root - Alignment, Quadrants, 2D symentrid and rod misalingment
12/// TPCCorrectionPrimitivesROC.root - ExB distortion due to the common ROC misalignment
13/// TPCCorrectionPrimitivesSector.root - ExB distortion due to the one sector misalingment (sector 0)
14/// TPCCorrectionPrimitivesFieldCage.root - ExB distortion due to the Field cage missalignment
15///
16///
17/// RegisterCorrection(); - Reserved id's 0 -999
18///
19/// RegisterAliTPCFCVoltError3D(); - Reserved id's 0 -99
20/// RegisterAliTPCBoundaryVoltError(); - Reserved id's 100-199
21/// RegisterAliTPCCalibGlobalMisalignment(); - Reserved id's 200-499
22/// RegisterAliTPCExBBShape(); - Reserved id's 500-600
23/// RegisterAliTPCExBTwist(); - Reserved id's 600-700
24/// RegisterAliTPCROCVoltError3D() - Reserved id's 700-800
25/// RegisterAliTPCROCVoltError3DSector() - Reserved id's 800-900
26///
27/// ~~~
28/// .x ~/rootlogon.C
29/// .L $ALICE_ROOT/TPC/CalibMacros/RegisterCorrection.C+
30/// RegisterCorrection();
31/// ~~~
32///
33/// Example use:
34///
35/// ~~~
36/// .x ~/rootlogon.C
37/// .L $ALICE_ROOT/PWGPP/CalibMacros/CPass0/ConfigCalibTrain.C
38/// ConfigCalibTrain(119037,"local:///cvmfs/alice.gsi.de/alice/data/2010/OCDB/")
39///
40/// .L $ALICE_ROOT/TPC/CalibMacros/RegisterCorrection.C+
41/// RegisterCorrection(0);
42/// ~~~
43///
44/// See example usage of correction primitive/derivatives in file
45
46#if !defined(__CINT__) || defined(__MAKECINT__)
47#include "TFile.h"
48#include "TObjArray.h"
49#include "TTreeStream.h"
50#include "TMath.h"
51#include "TGraph.h"
52#include "TRandom.h"
53#include "TTree.h"
54#include "TF1.h"
55#include "TH2F.h"
56#include "TH3F.h"
57#include "TAxis.h"
58#include "TPad.h"
59#include "TCanvas.h"
60#include "TStyle.h"
61#include "AliTPCParamSR.h"
62#include "TDatabasePDG.h"
63#include "AliTPCCorrection.h"
64#include "AliTPCFCVoltError3D.h"
65#include "AliTPCROCVoltError3D.h"
66#include "AliTPCComposedCorrection.h"
67#include "AliTPCBoundaryVoltError.h"
68#include "AliTPCCalibGlobalMisalignment.h"
69#include "AliTPCExBBShape.h"
70#include "AliTPCExBTwist.h"
71#include "AliTPCCorrectionDrift.h"
72#include "AliMagF.h"
73#include "AliCDBEntry.h"
74#include "AliTPCROC.h"
75#include <TStatToolkit.h>
76#include "TCut.h"
77#include "TGraphErrors.h"
78#include "AliTrackerBase.h"
79#include "TGeoGlobalMagField.h"
80#include "TROOT.h"
81#include "TLinearFitter.h"
82#include "TStopwatch.h"
83#include "AliTPCCorrectionFit.h"
84#endif
85
86
87
88
89TFile *fCorrections=0; //file with corrections
90//
91//models E field distortion AliTPCFCVoltError3D
92//
93AliTPCFCVoltError3D *rotOFC =0; // fit models
94AliTPCFCVoltError3D *rodOFC1 =0;
95AliTPCFCVoltError3D *rodOFC2 =0;
96AliTPCFCVoltError3D *rotIFC =0;
97AliTPCFCVoltError3D *rodIFC1 =0;
98AliTPCFCVoltError3D *rodIFC2 =0;
99//
100AliTPCFCVoltError3D *rodIFCShift=0; //common IFC shift
101AliTPCFCVoltError3D *rodIFCSin=0; //common IFC shift -sinus
102AliTPCFCVoltError3D *rodOFCSin=0; //common OFC shift -sinus
103AliTPCFCVoltError3D *rodOFCShift=0; //common OFC shift
104AliTPCFCVoltError3D *rodIFCCos=0; //common IFC shift -cos
105AliTPCFCVoltError3D *rodOFCCos=0; //common OFC shift -cos
106//
107TObjArray *rodFCSideRadiusType=0; // field cage shift, A/C side, OFC/IFC, Fourier(const,sin,cos)
108//
109
110
111
112AliTPCROCVoltError3D *rocRotgXA=0; // roc rotation A side - inclination in X
113AliTPCROCVoltError3D *rocRotgYA=0; // roc rotation A side - inclination in Y
114AliTPCROCVoltError3D *rocRotgXC=0; // roc rotation C side - inclination in X
115AliTPCROCVoltError3D *rocRotgYC=0; // roc rotation C side - inclination in Y
116//
117AliTPCROCVoltError3D *rocDzIROCA=0; // roc shift A side - in Z
118AliTPCROCVoltError3D *rocDzIROCC=0; // roc shift C side - in Z
119AliTPCROCVoltError3D *rocRotIROCA=0; // roc rotation A side - in Z
120AliTPCROCVoltError3D *rocRotIROCC=0; // roc rotation C side - in Z
121AliTPCROCVoltError3D *rocDzUDA=0; // roc shift Up-Down A side
122AliTPCROCVoltError3D *rocDzUDC=0; // roc shift Up-Down C side
123AliTPCROCVoltError3D *rocRotUDA=0; // roc rotation updown A side
124AliTPCROCVoltError3D *rocRotUDC=0; // roc rotation updown C side
125
126//
127AliTPCROCVoltError3D *rocShiftIROCA0=0; // IROC shift A0 side
128AliTPCROCVoltError3D *rocRotIROCA0=0; // IROC rot A0 side
129AliTPCROCVoltError3D *rocShiftOROCA0=0; // OROC shift A0 side
130AliTPCROCVoltError3D *rocRotOROCA0=0; // OROC rot A0 side
131AliTPCROCVoltError3D *rocShiftIROCC0=0; // IROC shift C0 side
132AliTPCROCVoltError3D *rocRotIROCC0=0; // IROC rot C0 side
133AliTPCROCVoltError3D *rocShiftOROCC0=0; // OROC shift C0 side
134AliTPCROCVoltError3D *rocRotOROCC0=0; // OROC rot C0 side
135//
136//
137//
138AliTPCBoundaryVoltError *boundaryVoltErrorA[8]; // boundary voltage error A side
139AliTPCBoundaryVoltError *boundaryVoltErrorC[8]; // boundary voltage error C side
140AliTPCExBBShape *exbShape = 0; // nominal correctin
141AliTPCExBBShape *exbShapeT1X = 0; // nominal +deltaT1=0.1
142AliTPCExBBShape *exbShapeT2X = 0; // nominal +deltaT2=0.1
143AliTPCExBTwist *twistX = 0;
144AliTPCExBTwist *twistY = 0;
145//
146AliTPCCalibGlobalMisalignment *alignRot0=0;
147AliTPCCalibGlobalMisalignment *alignRot1=0;
148AliTPCCalibGlobalMisalignment *alignRot2=0;
149AliTPCCalibGlobalMisalignment *alignTrans0=0;
150AliTPCCalibGlobalMisalignment *alignTrans1=0;
151AliTPCCalibGlobalMisalignment *alignTrans2=0;
152//
153AliTPCCalibGlobalMisalignment *alignTrans0D[4]={0}; // delta alignemnt 4 quadrants
154AliTPCCalibGlobalMisalignment *alignTrans1D[4]={0}; //
155AliTPCCalibGlobalMisalignment *alignTrans2D[4]={0}; //
156AliTPCCalibGlobalMisalignment *alignRot0D[4]={0};
157AliTPCCalibGlobalMisalignment *alignRot1D[4]={0};
158AliTPCCalibGlobalMisalignment *alignRot2D[4]={0};
159
160
161
162
163AliTPCCorrectionDrift *calibDrift[7]={0};
164
165
166//
167void RegisterAliTPCFCVoltError3D();
168void RegisterAliTPCCalibGlobalMisalignment();
169void RegisterAliTPCBoundaryVoltError();
170void RegisterAliTPCExBShape();
171void RegisterAliTPCExBTwist();
172void RegisterAliTPCROCVoltError3D();
173void RegisterAliTPCROCVoltError3DSector();
174void RegisterAliTPCCorrectionDrift();
175void RegisterAliTPCFCVoltError3DRodFCSideRadiusType();
176//
177
178void RegisterCorrection(Int_t type=0){
179 /// check the presence of corrections in file
180 ///
181 /// gROOT->Macro("ConfigCalibTrain.C(119037)");
182
183 if (type==1) return RegisterAliTPCROCVoltError3D(); // 3D distortion due misalignemnt of FC ....
184 if (type==2) return RegisterAliTPCROCVoltError3DSector(); // 3D distortion due misalingment of Sectors
185 if (type==4) return RegisterAliTPCFCVoltError3DRodFCSideRadiusType();
186 //
187 if (type==3) { // 2D distortions
188 fCorrections = TFile::Open("TPCCorrectionPrimitives.root","recreate");
189 RegisterAliTPCCalibGlobalMisalignment();
190 RegisterAliTPCBoundaryVoltError();
191 RegisterAliTPCFCVoltError3D();
192 RegisterAliTPCExBShape();
193 RegisterAliTPCExBTwist();
194 RegisterAliTPCCorrectionDrift();
195 fCorrections->Close();
196 delete fCorrections;
197 return;
198 }
199 fCorrections = TFile::Open("TPCCorrectionPrimitives.root");
200 AliTPCComposedCorrection *corrField3D = (AliTPCComposedCorrection*) fCorrections->Get("TPCFCVoltError3D");
201 // if not make new file
202 if (!corrField3D) fCorrections = TFile::Open("TPCCorrectionPrimitives.root","update");
203 if (type==0) {
204 RegisterAliTPCROCVoltError3D();
205 RegisterAliTPCROCVoltError3DSector();
206 RegisterAliTPCFCVoltError3DRodFCSideRadiusType();
207 }
208 RegisterAliTPCCalibGlobalMisalignment();
209 RegisterAliTPCBoundaryVoltError();
210 RegisterAliTPCFCVoltError3D();
211 RegisterAliTPCExBShape();
212 RegisterAliTPCExBTwist();
213 RegisterAliTPCCorrectionDrift();
214 if (fCorrections) fCorrections->Close();
215 //
216}
217
218
219
220void RegisterAliTPCFCVoltError3D(){
221 /// Load the models from the file
222 /// Or create it
223 /// Register functions with following IDs:
224 /// IMPORTANT: The nominal shift is in mm
225 ///
226 /// rotOFC - 0
227 /// rodOFC1 - 1
228 /// rodOFC2 - 2
229 /// rotIFC - 3
230 /// rodIFC1 - 4
231 /// rodIFC2 - 5
232 /// rodIFCShift - 6
233 /// rodIFCSin - 7
234 /// rodIFCCos - 8
235 /// rodOFCShift - 9
236 /// rodOFCSin - 10
237 /// rodOFCCos - 11
238
239 printf("RegisterAliTPCFCVoltError3D()");
240 Int_t volt = 40; // 40 V ~ 1mm
241 AliTPCComposedCorrection *corrField3D = (AliTPCComposedCorrection*) fCorrections->Get("TPCFCVoltError3D");
242 if (corrField3D) { // load form file
243 corrField3D->Print();
244 TCollection *iter = corrField3D->GetCorrections();
245 rotOFC = (AliTPCFCVoltError3D*)iter->FindObject("rotOFC");
246 rodOFC1 = (AliTPCFCVoltError3D*)iter->FindObject("rodOFC1");
247 rodOFC2 = (AliTPCFCVoltError3D*)iter->FindObject("rodOFC2");
248 rotIFC = (AliTPCFCVoltError3D*)iter->FindObject("rotIFC");
249 rodIFC1 = (AliTPCFCVoltError3D*)iter->FindObject("rodIFC1");
250 rodIFC2 = (AliTPCFCVoltError3D*)iter->FindObject("rodIFC2");
251 //
252 rodIFCShift = (AliTPCFCVoltError3D*)iter->FindObject("rodIFCShift");
253 rodOFCShift = (AliTPCFCVoltError3D*)iter->FindObject("rodOFCShift");
254 rodIFCSin = (AliTPCFCVoltError3D*)iter->FindObject("rodIFCSin");
255 rodOFCSin = (AliTPCFCVoltError3D*)iter->FindObject("rodOFCSin");
256 rodIFCCos = (AliTPCFCVoltError3D*)iter->FindObject("rodIFCCos");
257 rodOFCCos = (AliTPCFCVoltError3D*)iter->FindObject("rodOFCCos");
258 } else {
259 // OFC
260 rotOFC = new AliTPCFCVoltError3D();
261 rotOFC->SetOmegaTauT1T2(0,1,1);
262 rotOFC->SetRotatedClipVoltA(1,volt);
263 rotOFC->SetRotatedClipVoltC(1,volt);
264 //
265 rodOFC1 = new AliTPCFCVoltError3D();
266 rodOFC1->SetOmegaTauT1T2(0,1,1);
267 rodOFC1->SetRodVoltShiftA(18,volt);
268 rodOFC1->SetRodVoltShiftC(18,volt);
269 //
270 rodOFC2 = new AliTPCFCVoltError3D();
271 rodOFC2->SetOmegaTauT1T2(0,1,1);
272 rodOFC2->SetCopperRodShiftA(18,volt);
273 rodOFC2->SetCopperRodShiftC(18,volt);
274 // IFC
275 rotIFC = new AliTPCFCVoltError3D();
276 rotIFC->SetOmegaTauT1T2(0,1,1);
277 rotIFC->SetRotatedClipVoltA(0,volt);
278 rotIFC->SetRotatedClipVoltC(0,volt);
279 //
280 rodIFC1 = new AliTPCFCVoltError3D();
281 rodIFC1->SetOmegaTauT1T2(0,1,1);
282 rodIFC1->SetRodVoltShiftA(0,volt);
283 rodIFC1->SetRodVoltShiftC(0,volt);
284 //
285 rodIFC2 = new AliTPCFCVoltError3D();
286 rodIFC2->SetOmegaTauT1T2(0,1,1);
287 rodIFC2->SetCopperRodShiftA(0,volt);
288 rodIFC2->SetCopperRodShiftC(0,volt);
289 //
290 rodIFCShift = new AliTPCFCVoltError3D();
291 rodIFCSin = new AliTPCFCVoltError3D();
292 rodIFCCos = new AliTPCFCVoltError3D();
293 rodOFCShift = new AliTPCFCVoltError3D();
294 rodOFCSin = new AliTPCFCVoltError3D();
295 rodOFCCos = new AliTPCFCVoltError3D();
296 for (Int_t isec=0; isec<18; isec++){
297 Double_t phi=TMath::Pi()*isec/9.;
298 rodIFCShift->SetOmegaTauT1T2(0,1,1);
299 rodIFCShift->SetRodVoltShiftA(isec,volt);
300 rodIFCShift->SetRodVoltShiftC(isec,volt);
301 rodIFCSin->SetOmegaTauT1T2(0,1,1);
302 rodIFCSin->SetRodVoltShiftA(isec,volt*TMath::Sin(phi));
303 rodIFCSin->SetRodVoltShiftC(isec,volt*TMath::Sin(phi));
304 rodIFCCos->SetOmegaTauT1T2(0,1,1);
305 rodIFCCos->SetRodVoltShiftA(isec,volt*TMath::Cos(phi));
306 rodIFCCos->SetRodVoltShiftC(isec,volt*TMath::Cos(phi));
307 //
308 rodOFCShift->SetOmegaTauT1T2(0,1,1);
309 rodOFCShift->SetRodVoltShiftA(18+isec,volt);
310 rodOFCShift->SetRodVoltShiftC(18+isec,volt);
311 rodOFCSin->SetOmegaTauT1T2(0,1,1);
312 rodOFCSin->SetRodVoltShiftA(18+isec,volt*TMath::Sin(phi));
313 rodOFCSin->SetRodVoltShiftC(18+isec,volt*TMath::Sin(phi));
314 rodOFCCos->SetOmegaTauT1T2(0,1,1);
315 rodOFCCos->SetRodVoltShiftA(18+isec,volt*TMath::Cos(phi));
316 rodOFCCos->SetRodVoltShiftC(18+isec,volt*TMath::Cos(phi));
317 }
318 //
319 //
320 // Initialization of the lookup tables
321 //
322 printf(" ------- OFC rotated clip:\n"); rotOFC->InitFCVoltError3D();
323 printf(" ------- OFC rod & strip:\n"); rodOFC1->InitFCVoltError3D();
324 printf(" ------- OFC copper rod:\n"); rodOFC2->InitFCVoltError3D();
325 printf(" ------- IFC rotated clip:\n"); rotIFC->InitFCVoltError3D();
326 printf(" ------- IFC rod & strip:\n"); rodIFC1->InitFCVoltError3D();
327 printf(" ------- IFC copper rod:\n"); rodIFC2->InitFCVoltError3D();
328
329 printf(" ------- IFC rod & strip shift:\n"); rodIFCShift->InitFCVoltError3D();
330 printf(" ------- IFC rod & strip sin:\n"); rodIFCSin->InitFCVoltError3D();
331 printf(" ------- IFC rod & strip cos:\n"); rodIFCCos->InitFCVoltError3D();
332 printf(" ------- OFC rod & strip shift:\n"); rodOFCShift->InitFCVoltError3D();
333 printf(" ------- OFC rod & strip sin:\n"); rodOFCSin->InitFCVoltError3D();
334 printf(" ------- OFC rod & strip cos:\n"); rodOFCCos->InitFCVoltError3D();
335
336 // give names
337 rotOFC->SetName("rotOFC");rotOFC->SetTitle("rotOFC");
338 rodOFC1->SetName("rodOFC1");rodOFC1->SetTitle("rodOFC1");
339 rodOFC2->SetName("rodOFC2");rodOFC2->SetTitle("rodOFC2");
340 rotIFC->SetName("rotIFC");rotIFC->SetTitle("rotIFC");
341 rodIFC1->SetName("rodIFC1");rodIFC1->SetTitle("rodIFC1");
342 rodIFC2->SetName("rodIFC2");rodIFC2->SetTitle("rodIFC2");
343 //
344 rodIFCShift->SetName("rodIFCShift");rodIFCShift->SetTitle("rodIFCShift");
345 rodIFCSin->SetName("rodIFCSin");rodIFCSin->SetTitle("rodIFCSin");
346 rodIFCCos->SetName("rodIFCCos");rodIFCCos->SetTitle("rodIFCCos");
347 //
348 rodOFCShift->SetName("rodOFCShift");rodOFCShift->SetTitle("rodOFCShift");
349 rodOFCSin->SetName("rodOFCSin");rodOFCSin->SetTitle("rodOFCSin");
350 rodOFCCos->SetName("rodOFCCos");rodOFCCos->SetTitle("rodOFCCos");
351 //
352 // save in file
353 corrField3D = new AliTPCComposedCorrection();
354 TObjArray *cs = new TObjArray();
355 cs->Add(rotIFC); cs->Add(rotOFC);
356 cs->Add(rodIFC1); cs->Add(rodOFC1);
357 cs->Add(rodIFC2); cs->Add(rodOFC2);
358 cs->Add(rodIFCShift); cs->Add(rodIFCSin); cs->Add(rodIFCCos);
359 cs->Add(rodOFCShift); cs->Add(rodOFCSin); cs->Add(rodOFCCos);
360 //
361 corrField3D->SetCorrections(cs);
362 corrField3D->SetOmegaTauT1T2(0,1.,1.);
363 corrField3D->Print();
364 fCorrections->cd();
365 corrField3D->Write("TPCFCVoltError3D");
366 }
367 //
368 AliTPCCorrection::AddVisualCorrection(rotOFC,0);
369 AliTPCCorrection::AddVisualCorrection(rodOFC1,1);
370 AliTPCCorrection::AddVisualCorrection(rodOFC2,2);
371 AliTPCCorrection::AddVisualCorrection(rotIFC,3);
372 AliTPCCorrection::AddVisualCorrection(rodIFC1,4);
373 AliTPCCorrection::AddVisualCorrection(rodIFC2,5);
374 // common corrections
375 //
376 AliTPCCorrection::AddVisualCorrection(rodIFCShift,6);
377 AliTPCCorrection::AddVisualCorrection(rodIFCSin,7);
378 AliTPCCorrection::AddVisualCorrection(rodIFCCos,8);
379 //
380 AliTPCCorrection::AddVisualCorrection(rodOFCShift,9);
381 AliTPCCorrection::AddVisualCorrection(rodOFCSin,10);
382 AliTPCCorrection::AddVisualCorrection(rodOFCCos,11);
383}
384
385
386void RegisterAliTPCFCVoltError3DRodFCSideRadiusType(){
387 /// Load the models from the file
388 /// Or create it
389 /// Register functions with following IDs:
390 /// IMPORTANT: The nominal shift is in mm
391 ///
392 /// naming convention:
393 /// rodFCSide%dRadius%dType%d
394
395 ::Info("RegisterAliTPCFCVoltError3DRodFCSideRadiusType()","Start");
396 Int_t volt = 40; // 40 V ~ 1mm
397 TFile * fCorrectionsRodFCSideRadiusType = TFile::Open("TPCCorrectionPrimitivesFieldCage.root","update");
398 if (!fCorrectionsRodFCSideRadiusType) fCorrectionsRodFCSideRadiusType = TFile::Open("TPCCorrectionPrimitivesFieldCage.root","recreate");
399 rodFCSideRadiusType= new TObjArray(100);
400 //
401 AliTPCComposedCorrection *corrFieldCage = (AliTPCComposedCorrection*) fCorrectionsRodFCSideRadiusType->Get("TPCRodFCSideRadiusType");
402 if (corrFieldCage) { // load from file
403 TCollection *corrections = corrFieldCage->GetCorrections();
404 for (Int_t itype=0; itype<3; itype++){
405 for (Int_t iside=0; iside<2; iside++){
406 for (Int_t ifc=0; ifc<2; ifc++){
407 Int_t id=4*itype+2*iside+ifc;
408 AliTPCFCVoltError3D *corr = (AliTPCFCVoltError3D *)corrections->FindObject(TString::Format("rodFCSide%dRadius%dType%d",iside,ifc,itype).Data());
409 rodFCSideRadiusType->AddLast(corr);
410 ::Info("Read AliTPCFCVoltError3DRodFCSideRadiusType", TString::Format("Init: rodFCSide%dRadius%dType%d\tId=%d",iside,ifc,itype,900+id).Data());
411 AliTPCCorrection::AddVisualCorrection(corr,900+id);
412 }
413 }
414 }
415 corrFieldCage->Print();
416 } else {
417 for (Int_t itype=0; itype<3; itype++){
418 for (Int_t iside=0; iside<2; iside++){
419 for (Int_t ifc=0; ifc<2; ifc++){
420 Int_t id=4*itype+2*iside+ifc;
421 AliTPCFCVoltError3D *corr = new AliTPCFCVoltError3D();
422 corr->SetName(TString::Format("rodFCSide%dRadius%dType%d",iside,ifc,itype).Data());
423 corr->SetTitle(TString::Format("rodFCSide%dRadius%dType%d",iside,ifc,itype).Data());
424 //
425 for (Int_t isec=0; isec<18; isec++){
426 Double_t phi=TMath::Pi()*isec/9.;
427 Int_t sectorOffset=(ifc==0) ? 0:18;
428 corr->SetOmegaTauT1T2(0,1,1);
429 Double_t sectorVoltage=0;
430 if (itype==0) sectorVoltage=volt;
431 if (itype==1) sectorVoltage=volt*TMath::Sin(phi);
432 if (itype==2) sectorVoltage=volt*TMath::Cos(phi);
433 if (iside==0){
434 corr->SetRodVoltShiftA(isec+sectorOffset,sectorVoltage);
435 }
436 if (iside==1){
437 corr->SetRodVoltShiftC(isec+sectorOffset,sectorVoltage);
438 }
439 }
440 rodFCSideRadiusType->AddLast(corr);
441 ::Info("Generate AliTPCFCVoltError3DRodFCSideRadiusType", TString::Format("Init: rodFCSide%dRadius%dType%d",iside,ifc,itype).Data());
442 corr->InitFCVoltError3D();
443 AliTPCCorrection::AddVisualCorrection(corr,900+id);
444 }
445 }
446 }
447 //
448 corrFieldCage = new AliTPCComposedCorrection();
449 corrFieldCage->SetCorrections(rodFCSideRadiusType);
450 corrFieldCage->SetOmegaTauT1T2(0,1.,1.);
451
452 corrFieldCage->Print();
453 fCorrectionsRodFCSideRadiusType->cd();
454 corrFieldCage->Write("TPCRodFCSideRadiusType");
455 }
456 fCorrectionsRodFCSideRadiusType->Close();
457 ::Info("RegisterAliTPCFCVoltError3DRodFCSideRadiusType()","End");
458}
459
460
461void RegisterAliTPCCalibGlobalMisalignment(){
462 /// Register primitive alignment components.
463 /// Linear conbination of primitev forulas used for fit
464 /// The nominal delta 1 mm in shift and 1 mrad in rotation
465 /// Primitive formulas registeren in AliTPCCoreection::AddvisualCorrection
466 /// 20 - deltaX
467 /// 21 - deltaY
468 /// 22 - deltaZ
469 /// 23 - rot0 (phi)
470 /// 24 - rot1 (theta)
471 /// 25 - rot2
472
473 printf("RegisterAliTPCCalibGlobalMisalignment()\n");
474 TGeoHMatrix matrixX;
475 TGeoHMatrix matrixY;
476 TGeoHMatrix matrixZ;
477 TGeoRotation rot0;
478 TGeoRotation rot1;
479 TGeoRotation rot2; //transformation matrices
480 TGeoRotation rot90; //transformation matrices
481 matrixX.SetDx(0.1); matrixY.SetDy(0.1); matrixZ.SetDz(0.1); //1 mm translation
482 rot0.SetAngles(0.001*TMath::RadToDeg(),0,0);
483 rot1.SetAngles(0,0.001*TMath::RadToDeg(),0);
484 rot2.SetAngles(0,0,0.001*TMath::RadToDeg());
485 //how to get rot02 ?
486 rot90.SetAngles(0,90,0);
487 rot2.MultiplyBy(&rot90,kTRUE);
488 rot90.SetAngles(0,-90,0);
489 rot2.MultiplyBy(&rot90,kFALSE);
490 //
491 alignRot0 =new AliTPCCalibGlobalMisalignment;
492 alignRot0->SetAlignGlobal(&rot0);
493 alignRot0->SetName("alignRot0");
494 alignRot1=new AliTPCCalibGlobalMisalignment;
495 alignRot1->SetAlignGlobal(&rot1);
496 alignRot1->SetName("alignRot1");
497 alignRot2=new AliTPCCalibGlobalMisalignment;
498 alignRot2->SetAlignGlobal(&rot2);
499 alignRot2->SetName("alignRot2");
500 //
501 alignTrans0 =new AliTPCCalibGlobalMisalignment;
502 alignTrans0->SetAlignGlobal(&matrixX);
503 alignTrans0->SetName("alignTrans0");
504 alignTrans1=new AliTPCCalibGlobalMisalignment;
505 alignTrans1->SetAlignGlobal(&matrixY);
506 alignTrans1->SetName("alignTrans1");
507 alignTrans2=new AliTPCCalibGlobalMisalignment;
508 alignTrans2->SetAlignGlobal(&matrixZ);
509 alignTrans2->SetName("alignTrans2");
510 //
511
512
513 alignRot0D[0] =new AliTPCCalibGlobalMisalignment;
514 alignRot0D[0]->SetAlignGlobalDelta(&rot0);
515 alignRot0D[0]->SetName("alignRot0D");
516 alignRot1D[0]=new AliTPCCalibGlobalMisalignment;
517 alignRot1D[0]->SetAlignGlobalDelta(&rot1);
518 alignRot1D[0]->SetName("alignRot1D");
519 alignRot2D[0]=new AliTPCCalibGlobalMisalignment;
520 alignRot2D[0]->SetAlignGlobalDelta(&rot2);
521 alignRot2D[0]->SetName("alignRot2D");
522 //
523 alignTrans0D[0] =new AliTPCCalibGlobalMisalignment;
524 alignTrans0D[0]->SetAlignGlobalDelta(&matrixX);
525 alignTrans0D[0]->SetName("alignTrans0D");
526 alignTrans1D[0]=new AliTPCCalibGlobalMisalignment;
527 alignTrans1D[0]->SetAlignGlobalDelta(&matrixY);
528 alignTrans1D[0]->SetName("alignTrans1D");
529 alignTrans2D[0]=new AliTPCCalibGlobalMisalignment;
530 alignTrans2D[0]->SetAlignGlobalDelta(&matrixZ);
531 alignTrans2D[0]->SetName("alignTrans2D");
532
533 TObjArray * arrayDX = new TObjArray(72);
534 TObjArray * arrayDY = new TObjArray(72);
535 TObjArray * arrayDPhi = new TObjArray(72);
536 //
537 // Up down A side
538 //
539 for (Int_t isec=0; isec<72; isec++){ //A side
540 TGeoHMatrix *matrixDX = new TGeoHMatrix;
541 TGeoHMatrix *matrixDY = new TGeoHMatrix;
542 TGeoRotation *matrixDPhi= new TGeoRotation;
543 arrayDX->AddAt(matrixDX,isec);
544 arrayDY->AddAt(matrixDY,isec);
545 arrayDPhi->AddAt(matrixDPhi,isec);
546 if (isec%36<18) matrixDX->SetDx(isec%18<9?0.1:-0.1);
547 if (isec%36<18) matrixDY->SetDy(isec%18<9?0.1:-0.1);
548 if (isec%36<18) matrixDPhi->SetAngles((isec%18<9?0.001:-0.001)*TMath::RadToDeg(),0,0);
549 }
550 alignTrans0D[1] =new AliTPCCalibGlobalMisalignment;
551 alignTrans0D[1]->SetName("alignTrans0UDA");
552 alignTrans0D[1]->SetAlignSectors((TObjArray *)(arrayDX->Clone()));
553 alignTrans1D[1] =new AliTPCCalibGlobalMisalignment;
554 alignTrans1D[1]->SetName("alignTrans1UDA");
555 alignTrans1D[1]->SetAlignSectors((TObjArray *)(arrayDY->Clone()));
556 alignRot0D[1] =new AliTPCCalibGlobalMisalignment;
557 alignRot0D[1]->SetName("alignRot0UDA");
558 alignRot0D[1]->SetAlignSectors((TObjArray *)(arrayDPhi->Clone()));
559 //
560 // Uu-down C side
561 //
562 for (Int_t isec=0; isec<72; isec++){ //A side
563 TGeoHMatrix *matrixDX = new TGeoHMatrix;
564 TGeoHMatrix *matrixDY = new TGeoHMatrix;
565 TGeoRotation *matrixDPhi= new TGeoRotation;
566 arrayDX->AddAt(matrixDX,isec);
567 arrayDY->AddAt(matrixDY,isec);
568 arrayDPhi->AddAt(matrixDPhi,isec);
569 if (isec%36>=18) matrixDX->SetDx(isec%18<9?0.1:-0.1);
570 if (isec%36>=18) matrixDY->SetDy(isec%18<9?0.1:-0.1);
571 if (isec%36>=18) matrixDPhi->SetAngles((isec%18<9?0.001:-0.001)*TMath::RadToDeg(),0,0);
572 }
573 alignTrans0D[2] =new AliTPCCalibGlobalMisalignment;
574 alignTrans0D[2]->SetName("alignTrans0UDC");
575 alignTrans0D[2]->SetAlignSectors((TObjArray *)(arrayDX->Clone()));
576 alignTrans1D[2] =new AliTPCCalibGlobalMisalignment;
577 alignTrans1D[2]->SetName("alignTrans1UDC");
578 alignTrans1D[2]->SetAlignSectors((TObjArray *)(arrayDY->Clone()));
579 alignRot0D[2] =new AliTPCCalibGlobalMisalignment;
580 alignRot0D[2]->SetName("alignRot0UDC");
581 alignRot0D[2]->SetAlignSectors((TObjArray *)(arrayDPhi->Clone()));
582 //
583
584
585 //
586 AliTPCCorrection::AddVisualCorrection(alignTrans0 ,200);
587 AliTPCCorrection::AddVisualCorrection(alignTrans1 ,201);
588 AliTPCCorrection::AddVisualCorrection(alignTrans2 ,202);
589 AliTPCCorrection::AddVisualCorrection(alignRot0 ,203);
590 AliTPCCorrection::AddVisualCorrection(alignRot1 ,204);
591 AliTPCCorrection::AddVisualCorrection(alignRot2 ,205);
592
593}
594
595
596void RegisterAliTPCBoundaryVoltError(){
597 /// Register phi symetric E filed distortions
598 /// 100-108 - A side 0 Field
599 /// 110-118 - C side 0 Field
600 /// 120-128 - A side +0.5 Field
601 /// 130-138 - C side +0.5 Field
602 /// 140-148 - A side -0.5 Field
603 /// 150-158 - C side -0.5 Field
604
605 Double_t vdrift = 2.64; // [cm/us] // to be updated: per second (ideally)
606 Double_t ezField = 400; // [V/cm] // to be updated: never (hopefully)
607 Double_t T1 = 1.0;
608 Double_t T2 = 1.0;
609 Double_t wtP = -10.0 * (0.5*10) * vdrift / ezField ;
610 Double_t wtM = -10.0 * (0.5*10) * vdrift / -ezField ;
611
612 printf("RegisterAliTPCBoundaryVoltError()\n");
613 AliTPCComposedCorrection *corrField2D = (AliTPCComposedCorrection*) fCorrections->Get("TPCFCVoltError2D");
614 //
615 if (!corrField2D){
616 TObjArray *array=new TObjArray(16);
617 Double_t val = 40.; // 1mm
618 Float_t bound0[8] = { 0, 0,0,0,0,0,0,0};
619 Float_t boundAi[8] = { 0, 0,0,0,0,0,0,0};
620 Float_t boundCi[8] = { 0, 0,0,0,0,0,0,0};
621 for (Int_t ipar=0; ipar<8; ipar++){
622 //
623 boundaryVoltErrorA[ipar] = new AliTPCBoundaryVoltError;
624 boundaryVoltErrorC[ipar] = new AliTPCBoundaryVoltError;
625 boundaryVoltErrorA[ipar]->SetName(Form("BoundaryVoltErrorAsidePar%d",ipar));
626 boundaryVoltErrorA[ipar]->SetTitle(Form("BoundaryVoltErrorAsidePar%d",ipar));
627 boundaryVoltErrorC[ipar]->SetName(Form("BoundaryVoltErrorCsidePar%d",ipar));
628 boundaryVoltErrorC[ipar]->SetTitle(Form("BoundaryVoltErrorCsidePar%d",ipar));
629 for (Int_t jpar=0; jpar<8; jpar++) if (ipar!=jpar){
630 boundAi[jpar]=0;
631 boundCi[jpar]=0;
632 }
633 boundAi[ipar]=val;
634 boundCi[ipar]=val;
635 //
636 boundaryVoltErrorA[ipar]->SetBoundariesA(boundAi);
637 boundaryVoltErrorA[ipar]->SetBoundariesC(bound0);
638 boundaryVoltErrorA[ipar]->InitBoundaryVoltErrorDistortion();
639 boundaryVoltErrorA[ipar]->SetOmegaTauT1T2(0.,1,1);
640 //
641 Float_t tempboundAi[8] = { 0, 0,0,0,0,0,-boundCi[6],-boundCi[7]};
642 boundaryVoltErrorC[ipar]->SetBoundariesA(tempboundAi);
643 boundaryVoltErrorC[ipar]->SetBoundariesC(boundCi);
644
645 boundaryVoltErrorC[ipar]->InitBoundaryVoltErrorDistortion();
646 boundaryVoltErrorC[ipar]->SetOmegaTauT1T2(0.,1,1);
647 array->AddAt(boundaryVoltErrorA[ipar],ipar);
648 array->AddAt(boundaryVoltErrorC[ipar],ipar+8);
649 boundaryVoltErrorA[ipar]->Print();
650 boundaryVoltErrorC[ipar]->Print();
651 AliTPCCorrection::AddVisualCorrection(boundaryVoltErrorA[ipar], 100+ipar);
652 AliTPCCorrection::AddVisualCorrection(boundaryVoltErrorC[ipar], 150+ipar);
653 }
654 corrField2D = new AliTPCComposedCorrection;
655 corrField2D->SetCorrections(array);
656 corrField2D->SetOmegaTauT1T2(0,1.,1.);
657 corrField2D->Print();
658 fCorrections->cd();
659 corrField2D->SetName("TPCFCVoltError2D");
660 corrField2D->SetTitle("TPCFCVoltError2D");
661 corrField2D->Write("TPCFCVoltError2D");
662 }else{
663 TObjArray *array = (TObjArray*)corrField2D->GetCorrections();
664 for (Int_t ipar=0; ipar<8; ipar++){
665 boundaryVoltErrorA[ipar] = (AliTPCBoundaryVoltError*) array->At(ipar);
666 boundaryVoltErrorC[ipar] = (AliTPCBoundaryVoltError*) array->At(ipar+8);
667 }
668 }
669 //
670 // Register correction
671 for (Int_t ipar=0; ipar<8; ipar++){
672 AliTPCCorrection::AddVisualCorrection(boundaryVoltErrorA[ipar], 100+ipar);
673 AliTPCCorrection::AddVisualCorrection(boundaryVoltErrorC[ipar], 110+ipar);
674 //
675 // correction for +-0.5 T setting
676 AliTPCCorrection *corrField =0;
677 corrField=(AliTPCCorrection *)boundaryVoltErrorA[ipar]->Clone();
678 corrField->SetOmegaTauT1T2(wtP,T1,T2);
679 AliTPCCorrection::AddVisualCorrection(corrField,120+ipar);
680
681 corrField=(AliTPCCorrection *)boundaryVoltErrorC[ipar]->Clone();
682 corrField->SetOmegaTauT1T2(wtP,T1,T2);
683 AliTPCCorrection::AddVisualCorrection(corrField,130+ipar);
684 // correction for +-0.5 T setting
685 corrField=(AliTPCCorrection *)boundaryVoltErrorA[ipar]->Clone();
686 corrField->SetOmegaTauT1T2(wtM,T1,T2);
687 AliTPCCorrection::AddVisualCorrection(corrField,140+ipar);
688
689 corrField=(AliTPCCorrection *)boundaryVoltErrorC[ipar]->Clone();
690 corrField->SetOmegaTauT1T2(wtM,T1,T2);
691 AliTPCCorrection::AddVisualCorrection(corrField,150+ipar);
692 }
693
694}
695
696
697
698void RegisterAliTPCExBShape(){
699 ///
700
701 AliMagF *magF = new AliMagF("mag","mag");
702
703 exbShape = new AliTPCExBBShape;
704 exbShape->SetBField(magF);
705 exbShape->SetName("TPCExBShape");
706 exbShape->SetTitle("TPCExBShape");
707 exbShape->SetOmegaTauT1T2(0,1.,1.);
708 exbShape->Print();
709 AliTPCCorrection::AddVisualCorrection(exbShape,500);
710 exbShapeT1X = new AliTPCExBBShape;
711 exbShapeT1X->SetBField(magF);
712 exbShapeT1X->SetName("TPCExbShapeT1X");
713 exbShapeT1X->SetTitle("TPCExbShapeT1X");
714 exbShapeT1X->SetOmegaTauT1T2(0,1.2,1.);
715 exbShapeT1X->Print();
716 AliTPCCorrection::AddVisualCorrection(exbShapeT1X,501);
717 exbShapeT2X = new AliTPCExBBShape;
718 exbShapeT2X->SetBField(magF);
719 exbShapeT2X->SetName("TPCExbShapeT2X");
720 exbShapeT2X->SetTitle("TPCExbShapeT2X");
721 exbShapeT2X->SetOmegaTauT1T2(0,1.0,1.2);
722 exbShapeT2X->Print();
723 AliTPCCorrection::AddVisualCorrection(exbShapeT2X,502);
724}
725
726
727void RegisterAliTPCExBTwist(){
728 ///
729
730 twistX = new AliTPCExBTwist;
731 twistY = new AliTPCExBTwist;
732 twistX->SetXTwist(0.001); // 1 mrad twist in x
733 twistX->SetName("ExBTwistX");
734 twistX->SetTitle("ExBTwistX");
735 twistY->SetYTwist(0.001); // 1 mrad twist in y
736 twistY->SetName("ExBTwistY");
737 twistY->SetTitle("ExBTwistY");
738 twistX->SetOmegaTauT1T2(0,1.,1.);
739 twistY->SetOmegaTauT1T2(0,1.,1.);
740 AliTPCCorrection::AddVisualCorrection(twistX,600);
741 AliTPCCorrection::AddVisualCorrection(twistY,601);
742}
743
744void RegisterAliTPCCorrectionDrift(){
745 /// Drift distortion/correction
746
747 for (Int_t idrift=0; idrift<7; idrift++) {
748 calibDrift[idrift]=new AliTPCCorrectionDrift;
749 }
750 calibDrift[0]->SetName("driftT0");
751 calibDrift[0]->SetTitle("driftT0");
752 calibDrift[0]->fZ0Aside=0.1;
753 calibDrift[0]->fZ0Cside=0.1;
754 calibDrift[1]->SetName("driftScale0");
755 calibDrift[1]->SetTitle("driftScale0");
756 calibDrift[1]->fVScale0=0.001;
757 calibDrift[2]->SetName("driftScaleR");
758 calibDrift[2]->SetTitle("driftScaleR");
759 calibDrift[2]->fVScaleR=0.001;
760 calibDrift[3]->SetName("driftScaleX");
761 calibDrift[3]->SetTitle("driftScaleX");
762 calibDrift[3]->fVScaleX=0.001;
763 calibDrift[4]->SetName("driftScaleY");
764 calibDrift[4]->SetTitle("driftScaleY");
765 calibDrift[4]->fVScaleY=0.001;
766 //
767 calibDrift[5]->SetName("driftIROCDZ");
768 calibDrift[5]->SetTitle("driftIROCDZ");
769 calibDrift[5]->fIROCZ0=0.1; //delta Z for IROCORC
770 //
771 calibDrift[6]->SetName("driftOROCDT");
772 calibDrift[6]->SetTitle("driftOROCDT");
773 calibDrift[6]->fOROCDZ=0.001; //delta Z for IROCORC
774
775
776}
777
778
779void RegisterAliTPCROCVoltError3D(){
780 /// ROC rotation transformation
781 /// 700 -709 - 0 field
782 /// 710 -719 - +0.5 field
783 /// 720 -729 - -0.5 field
784
785 Double_t vdrift = 2.64; // [cm/us] // to be updated: per second (ideally)
786 Double_t ezField = 400; // [V/cm] // to be updated: never (hopefully)
787 Double_t T1 = 1.0;
788 Double_t T2 = 1.0;
789 Double_t wtP = -10.0 * (0.5*10) * vdrift / ezField ;
790 Double_t wtM = -10.0 * (0.5*10) * vdrift / -ezField ;
791
792 //
793 rocRotgXA=0; // roc rotation A side - inclination in X
794 rocRotgYA=0; // roc rotation A side - inclination in Y
795 rocRotgXC=0; // roc rotation C side - inclination in X
796 rocRotgYC=0; // roc rotation C side - inclination in Y
797 rocDzIROCA=0; // roc shift A side - in Z
798 rocDzIROCC=0; // roc shift C side - in Z
799 //
800 rocDzUDA=0; // roc shift up-down - A side
801 rocDzUDC=0; // roc shift up-down - C side
802 //
803 rocRotIROCA=0; // roc rot IROC A side - in Z
804 rocRotIROCC=0; // roc rot OROC C side - in Z
805 rocRotUDA=0; // roc rot updown A side
806 rocRotUDC=0; // roc rot updown C side
807 //
808 printf("RegisterAliTPCROCVoltError3D()");
809 Double_t kAngle=0.001;
810 // reference in lx
811 AliTPCROC * rocInfo = AliTPCROC::Instance();
812 Double_t lxRef = (rocInfo->GetPadRowRadii(0,62)+rocInfo->GetPadRowRadii(36,0))/2;
813 //
814 TMatrixD matrix(72,3);
815
816 AliTPCComposedCorrection *corrField3D = 0;
817 TFile *fCorrectionsROC=0;
818 fCorrectionsROC = TFile::Open("TPCCorrectionPrimitivesROC.root");
819 corrField3D = ( AliTPCComposedCorrection *)fCorrectionsROC->Get("TPCROCVoltError3DRotationgXgY");
820 //
821 if (!corrField3D){
822 fCorrectionsROC = TFile::Open("TPCCorrectionPrimitivesROC.root","recreate");
823 }
824 if (corrField3D) { // load from file
825 corrField3D->Print();
826 TCollection *iter = corrField3D->GetCorrections();
827
828 rocRotgXA=(AliTPCROCVoltError3D*)iter->FindObject("rocRotgXA");
829
830 rocRotgYA=(AliTPCROCVoltError3D*)iter->FindObject("rocRotgYA");
831
832 rocRotgXC=(AliTPCROCVoltError3D*)iter->FindObject("rocRotgXC");
833
834 rocRotgYC=(AliTPCROCVoltError3D*)iter->FindObject("rocRotgYC");
835
836 rocDzIROCA=(AliTPCROCVoltError3D*)iter->FindObject("rocDzIROCA");
837
838 rocDzIROCC=(AliTPCROCVoltError3D*)iter->FindObject("rocDzIROCC");
839 //
840 rocDzUDA=(AliTPCROCVoltError3D*)iter->FindObject("rocDzUDA");
841
842 rocDzUDC=(AliTPCROCVoltError3D*)iter->FindObject("rocDzUDC");
843
844 rocRotIROCA=(AliTPCROCVoltError3D*)iter->FindObject("rocRotIROCA");
845
846 rocRotIROCC=(AliTPCROCVoltError3D*)iter->FindObject("rocRotIROCC");
847
848 rocRotUDA=(AliTPCROCVoltError3D*)iter->FindObject("rocRotUDA");
849
850 rocRotUDC=(AliTPCROCVoltError3D*)iter->FindObject("rocRotUDC");
851
852 } else {
853 corrField3D = new AliTPCComposedCorrection;
854 rocRotgXA=new AliTPCROCVoltError3D;
855 rocRotgYA=new AliTPCROCVoltError3D;
856 rocRotgXC=new AliTPCROCVoltError3D;
857 rocRotgYC=new AliTPCROCVoltError3D;
858 rocDzIROCA=new AliTPCROCVoltError3D;
859 rocDzIROCC=new AliTPCROCVoltError3D;
860 rocDzUDA=new AliTPCROCVoltError3D;
861 rocDzUDC=new AliTPCROCVoltError3D;
862 rocRotIROCA=new AliTPCROCVoltError3D;
863 rocRotIROCC=new AliTPCROCVoltError3D;
864 rocRotUDA=new AliTPCROCVoltError3D;
865 rocRotUDC=new AliTPCROCVoltError3D;
866 //
867 for (Int_t isec=0; isec<72; isec++){
868 Double_t secAlpha = TMath::DegToRad()*(10.+20.*(((Int_t)isec)%18));
869 matrix(isec,0)=0; matrix(isec,1)=0; matrix(isec,2)=0;
870 if (isec%36<18){
871 matrix(isec,0) = kAngle*TMath::Cos(secAlpha)*lxRef;
872 matrix(isec,1) = kAngle*TMath::Cos(secAlpha);
873 matrix(isec,2) = -kAngle*TMath::Sin(secAlpha);
874 }
875 }
876 rocRotgXA->SetROCData(&matrix);
877 //
878 for (Int_t isec=0; isec<72; isec++){
879 Double_t secAlpha = TMath::DegToRad()*(10.+20.*(((Int_t)isec)%18));
880 matrix(isec,0)=0; matrix(isec,1)=0; matrix(isec,2)=0;
881 if (isec%36<18){
882 matrix(isec,0) = kAngle*TMath::Sin(secAlpha)*lxRef;
883 matrix(isec,1) = kAngle*TMath::Sin(secAlpha);
884 matrix(isec,2) = kAngle*TMath::Cos(secAlpha);
885 }
886 }
887 rocRotgYA->SetROCData(&matrix);
888 //
889 for (Int_t isec=0; isec<72; isec++){
890 Double_t secAlpha = TMath::DegToRad()*(10.+20.*(((Int_t)isec)%18));
891 matrix(isec,0)=0; matrix(isec,1)=0; matrix(isec,2)=0;
892 if (isec%36>=18){
893 matrix(isec,0) = kAngle*TMath::Cos(secAlpha)*lxRef;
894 matrix(isec,1) = kAngle*TMath::Cos(secAlpha);
895 matrix(isec,2) = -kAngle*TMath::Sin(secAlpha);
896 }
897 }
898 rocRotgXC->SetROCData(&matrix);
899 //
900 for (Int_t isec=0; isec<72; isec++){
901 Double_t secAlpha = TMath::DegToRad()*(10.+20.*(((Int_t)isec)%18));
902 matrix(isec,0)=0; matrix(isec,1)=0; matrix(isec,2)=0;
903 if (isec%36>=18){
904 matrix(isec,0) = kAngle*TMath::Sin(secAlpha)*lxRef;
905 matrix(isec,1) = kAngle*TMath::Sin(secAlpha);
906 matrix(isec,2) = kAngle*TMath::Cos(secAlpha);
907 }
908 }
909 rocRotgYC->SetROCData(&matrix);
910
911 //
912 //
913 for (Int_t isec=0; isec<72; isec++){
914 matrix(isec,0)=0; matrix(isec,1)=0; matrix(isec,2)=0;
915 if (isec<18){
916 matrix(isec,0) = 0.1; // 1 mm
917 matrix(isec,1) = 0;
918 matrix(isec,2) = 0;
919 }
920 }
921 rocDzIROCA->SetROCData(&matrix);
922 //
923 for (Int_t isec=0; isec<72; isec++){
924 matrix(isec,0)=0; matrix(isec,1)=0; matrix(isec,2)=0;
925 if (isec>=18 && isec<36){
926 matrix(isec,0) = 0.1; // 1 mm
927 matrix(isec,1) = 0;
928 matrix(isec,2) = 0;
929 }
930 }
931 rocDzIROCC->SetROCData(&matrix);
932
933 //
934 //
935 for (Int_t isec=0; isec<72; isec++){
936 matrix(isec,0)=0; matrix(isec,1)=0; matrix(isec,2)=0;
937 if (isec%36<18){
938 matrix(isec,0) = (isec%18<9)? 0.05:-0.05; // 1 mm
939 matrix(isec,1) = 0;
940 matrix(isec,2) = 0;
941 }
942 }
943 rocDzUDA->SetROCData(&matrix);
944 //
945 for (Int_t isec=0; isec<72; isec++){
946 matrix(isec,0)=0; matrix(isec,1)=0; matrix(isec,2)=0;
947 if (isec%36>=18){
948 matrix(isec,0) = (isec%18<9)?0.05:-0.05; // 1 mm
949 matrix(isec,1) = 0;
950 matrix(isec,2) = 0;
951 }
952 }
953 rocDzUDC->SetROCData(&matrix);
954 //
955 //
956 for (Int_t isec=0; isec<72; isec++){
957 matrix(isec,0)=0; matrix(isec,1)=0; matrix(isec,2)=0;
958 if (isec<18){
959 matrix(isec,0) = 0;
960 matrix(isec,1) = kAngle;
961 matrix(isec,2) = 0;
962 }
963 }
964 //
965 rocRotIROCA->SetROCData(&matrix);
966 //
967 for (Int_t isec=0; isec<72; isec++){
968 matrix(isec,0)=0; matrix(isec,1)=0; matrix(isec,2)=0;
969 if (isec>=18 && isec<36){
970 matrix(isec,0) = 0;
971 matrix(isec,1) = kAngle;
972 matrix(isec,2) = 0;
973 }
974 }
975 rocRotIROCC->SetROCData(&matrix);
976 //
977 //
978 for (Int_t isec=0; isec<72; isec++){
979 matrix(isec,0)=0; matrix(isec,1)=0; matrix(isec,2)=0;
980 if (isec%36<18){
981 matrix(isec,0) = 0;
982 matrix(isec,1) = (isec%18<9)?kAngle:-kAngle;
983 matrix(isec,2) = 0;
984 }
985 }
986 //
987 rocRotUDA->SetROCData(&matrix);
988 //
989 for (Int_t isec=0; isec<72; isec++){
990 matrix(isec,0)=0; matrix(isec,1)=0; matrix(isec,2)=0;
991 if (isec%36>=18){
992 matrix(isec,0) = 0;
993 matrix(isec,1) = (isec%18<9)?kAngle:-kAngle;
994 matrix(isec,2) = 0;
995 }
996 }
997 rocRotUDC->SetROCData(&matrix);
998
999
1000 //
1001 rocRotgXA->SetElectronArrivalCorrection(kFALSE);
1002 rocRotgYA->SetElectronArrivalCorrection(kFALSE);
1003 rocRotgXC->SetElectronArrivalCorrection(kFALSE);
1004 rocRotgYC->SetElectronArrivalCorrection(kFALSE);
1005 rocDzIROCA->SetElectronArrivalCorrection(kFALSE);
1006 rocDzIROCC->SetElectronArrivalCorrection(kFALSE);
1007 rocDzUDA->SetElectronArrivalCorrection(kFALSE);
1008 rocDzUDC->SetElectronArrivalCorrection(kFALSE);
1009 rocRotIROCA->SetElectronArrivalCorrection(kFALSE);
1010 rocRotIROCC->SetElectronArrivalCorrection(kFALSE);
1011 rocRotUDA->SetElectronArrivalCorrection(kFALSE);
1012 rocRotUDC->SetElectronArrivalCorrection(kFALSE);
1013
1014 /* // verification plots
1015 rocRotgXA.CreateHistoOfZAlignment(0,500,500)->Draw("surf2");
1016 rocRotgYA.CreateHistoOfZAlignment(0,500,500)->Draw("surf2");
1017 rocRotgXC.CreateHistoOfZAlignment(1,500,500)->Draw("surf2");
1018 rocRotgYC.CreateHistoOfZAlignment(1,500,500)->Draw("surf2");
1019 */
1020
1021 //
1022 rocRotgXA->SetName("rocRotgXA");rocRotgXA->SetTitle("rocRotgXA");
1023 rocRotgYA->SetName("rocRotgYA");rocRotgYA->SetTitle("rocRotgYA");
1024 rocRotgXC->SetName("rocRotgXC");rocRotgXC->SetTitle("rocRotgXC");
1025 rocRotgYC->SetName("rocRotgYC");rocRotgYC->SetTitle("rocRotgYC");
1026 rocDzIROCA->SetName("rocDzIROCA");rocDzIROCA->SetTitle("rocDzIROCA");
1027 rocDzIROCC->SetName("rocDzIROCC");rocDzIROCC->SetTitle("rocDzIROCC");
1028 rocDzUDA->SetName("rocDzUDA");rocDzUDA->SetTitle("rocDzUDA");
1029 rocDzUDC->SetName("rocDzUDC");rocDzUDC->SetTitle("rocDzUDC");
1030 rocRotIROCA->SetName("rocRotIROCA");rocRotIROCA->SetTitle("rocRotIROCA");
1031 rocRotIROCC->SetName("rocRotIROCC");rocRotIROCC->SetTitle("rocRotIROCC");
1032 rocRotUDA->SetName("rocRotUDA");rocRotUDA->SetTitle("rocRotUDA");
1033 rocRotUDC->SetName("rocRotUDC");rocRotUDC->SetTitle("rocRotUDC");
1034 //
1035 //
1036 TObjArray *cs = new TObjArray();
1037 cs->Add(rocRotgXA);
1038 cs->Add(rocRotgYA);
1039 cs->Add(rocRotgXC);
1040 cs->Add(rocRotgYC);
1041 cs->Add(rocDzIROCA);
1042 cs->Add(rocDzIROCC);
1043 cs->Add(rocDzUDA);
1044 cs->Add(rocDzUDC);
1045 cs->Add(rocRotIROCA);
1046 cs->Add(rocRotIROCC);
1047 cs->Add(rocRotUDA);
1048 cs->Add(rocRotUDC);
1049 corrField3D->SetCorrections(cs);
1050 corrField3D->SetOmegaTauT1T2(0,1.,1.);
1051 corrField3D->Print();
1052 fCorrectionsROC->cd();
1053 corrField3D->Init();
1054 corrField3D->Print("da");
1055 fCorrectionsROC->cd();
1056 corrField3D->Write("TPCROCVoltError3DRotationgXgY");
1057 }
1058 AliTPCCorrection::AddVisualCorrection(rocRotgXA,701);
1059 AliTPCCorrection::AddVisualCorrection(rocRotgYA,702);
1060 AliTPCCorrection::AddVisualCorrection(rocRotgXC,703);
1061 AliTPCCorrection::AddVisualCorrection(rocRotgYC,704);
1062 AliTPCCorrection::AddVisualCorrection(rocDzIROCA,705);
1063 AliTPCCorrection::AddVisualCorrection(rocDzIROCC,706);
1064 AliTPCCorrection::AddVisualCorrection(rocDzUDA,709);
1065 AliTPCCorrection::AddVisualCorrection(rocDzUDC,710);
1066 AliTPCCorrection::AddVisualCorrection(rocRotIROCA,707);
1067 AliTPCCorrection::AddVisualCorrection(rocRotIROCC,708);
1068 AliTPCCorrection::AddVisualCorrection(rocRotUDA,711);
1069 AliTPCCorrection::AddVisualCorrection(rocRotUDC,712);
1070
1071 AliTPCCorrection *corrPlus =0;
1072 //
1073 corrPlus=(AliTPCCorrection *)rocRotgXA->Clone();
1074 corrPlus->SetOmegaTauT1T2(wtP,T1,T2);
1075 AliTPCCorrection::AddVisualCorrection(corrPlus,711);
1076 //
1077 corrPlus=(AliTPCCorrection *)rocRotgYA->Clone();
1078 corrPlus->SetOmegaTauT1T2(wtP,T1,T2);
1079 AliTPCCorrection::AddVisualCorrection(corrPlus,712);
1080 //
1081 corrPlus=(AliTPCCorrection *)rocRotgXC->Clone();
1082 corrPlus->SetOmegaTauT1T2(wtP,T1,T2);
1083 AliTPCCorrection::AddVisualCorrection(corrPlus,713);
1084 //
1085 corrPlus=(AliTPCCorrection *)rocRotgYC->Clone();
1086 corrPlus->SetOmegaTauT1T2(wtP,T1,T2);
1087 AliTPCCorrection::AddVisualCorrection(corrPlus,714);
1088 //
1089 corrPlus=(AliTPCCorrection *)rocDzIROCA->Clone();
1090 corrPlus->SetOmegaTauT1T2(wtP,T1,T2);
1091 AliTPCCorrection::AddVisualCorrection(corrPlus,715);
1092 //
1093 corrPlus=(AliTPCCorrection *)rocDzIROCC->Clone();
1094 corrPlus->SetOmegaTauT1T2(wtP,T1,T2);
1095 AliTPCCorrection::AddVisualCorrection(corrPlus,716);
1096 //
1097 corrPlus=(AliTPCCorrection *)rocRotIROCA->Clone();
1098 corrPlus->SetOmegaTauT1T2(wtP,T1,T2);
1099 AliTPCCorrection::AddVisualCorrection(corrPlus,717);
1100 //
1101 corrPlus=(AliTPCCorrection *)rocDzIROCC->Clone();
1102 corrPlus->SetOmegaTauT1T2(wtP,T1,T2);
1103 AliTPCCorrection::AddVisualCorrection(corrPlus,718);
1104 //
1105 //
1106 AliTPCCorrection *corrMinus =0;
1107 //
1108 corrMinus=(AliTPCCorrection *)rocRotgXA->Clone();
1109 corrMinus->SetOmegaTauT1T2(wtM,T1,T2);
1110 AliTPCCorrection::AddVisualCorrection(corrMinus,721);
1111 //
1112 corrMinus=(AliTPCCorrection *)rocRotgYA->Clone();
1113 corrMinus->SetOmegaTauT1T2(wtM,T1,T2);
1114 AliTPCCorrection::AddVisualCorrection(corrMinus,722);
1115 //
1116 corrMinus=(AliTPCCorrection *)rocRotgXC->Clone();
1117 corrMinus->SetOmegaTauT1T2(wtM,T1,T2);
1118 AliTPCCorrection::AddVisualCorrection(corrMinus,723);
1119 //
1120 corrMinus=(AliTPCCorrection *)rocRotgYC->Clone();
1121 corrMinus->SetOmegaTauT1T2(wtM,T1,T2);
1122 AliTPCCorrection::AddVisualCorrection(corrMinus,724);
1123 //
1124 corrMinus=(AliTPCCorrection *)rocDzIROCA->Clone();
1125 corrMinus->SetOmegaTauT1T2(wtM,T1,T2);
1126 AliTPCCorrection::AddVisualCorrection(corrMinus,725);
1127 //
1128 corrMinus=(AliTPCCorrection *)rocDzIROCC->Clone();
1129 corrMinus->SetOmegaTauT1T2(wtM,T1,T2);
1130 AliTPCCorrection::AddVisualCorrection(corrMinus,726);
1131 //
1132 corrMinus=(AliTPCCorrection *)rocRotIROCA->Clone();
1133 corrMinus->SetOmegaTauT1T2(wtM,T1,T2);
1134 AliTPCCorrection::AddVisualCorrection(corrMinus,727);
1135 //
1136 corrMinus=(AliTPCCorrection *)rocDzIROCC->Clone();
1137 corrMinus->SetOmegaTauT1T2(wtM,T1,T2);
1138 AliTPCCorrection::AddVisualCorrection(corrMinus,728);
1139 //
1140
1141 fCorrectionsROC->Close();
1142 delete fCorrectionsROC;
1143}
1144
1145
1146void RegisterAliTPCROCVoltError3DSector(){
1147 /// ROC rotation and shift transformation
1148 /// 800-819 - 0.0 Field
1149 /// 820-839 - +0.5 Field
1150 /// 840-859 - +0.5 Field
1151
1152 rocShiftIROCA0=0; // IROC shift A0 side
1153 rocRotIROCA0=0; // IROC rot A0 side
1154 rocShiftOROCA0=0; // OROC shift A0 side
1155 rocRotOROCA0=0; // OROC rot A0 side
1156 rocShiftIROCC0=0; // IROC shift C0 side
1157 rocRotIROCC0=0; // IROC rot C0 side
1158 rocShiftOROCC0=0; // OROC shift C0 side
1159 rocRotOROCC0=0; // OROC rot C0 side
1160 Double_t vdrift = 2.64; // [cm/us] // to be updated: per second (ideally)
1161 Double_t ezField = 400; // [V/cm] // to be updated: never (hopefully)
1162 Double_t T1 = 1.0;
1163 Double_t T2 = 1.0;
1164
1165 // Double_t wt = -10.0 * (bzField*10) * vdrift / ezField ;
1166
1167 //
1168 //
1169 printf("RegisterAliTPCROCVoltError3DSector()");
1170 Double_t kAngle=0.001;
1171 Double_t kDz=0.1;
1172 // reference in lx
1173 //
1174 TMatrixD matrix(72,3);
1175
1176 AliTPCComposedCorrection *corrField3DSector = 0;
1177 TFile *fCorrectionsROC=0;
1178 fCorrectionsROC = TFile::Open("TPCCorrectionPrimitivesSector.root");
1179 corrField3DSector = ( AliTPCComposedCorrection *)fCorrectionsROC->Get("TPCROCVoltError3DSector");
1180 //
1181 if (!corrField3DSector){
1182 fCorrectionsROC = TFile::Open("TPCCorrectionPrimitivesSector.root","recreate");
1183 }
1184 if (corrField3DSector) { // load from file
1185 corrField3DSector->Print();
1186 TCollection *iter = corrField3DSector->GetCorrections();
1187 //
1188 rocShiftIROCA0=(AliTPCROCVoltError3D*)iter->FindObject("rocShiftIROCA0"); // IROC shift A0 side
1189 rocRotIROCA0=(AliTPCROCVoltError3D*)iter->FindObject("rocRotIROCA0"); // IROC rot A0 side
1190 rocShiftOROCA0=(AliTPCROCVoltError3D*)iter->FindObject("rocShiftOROCA0"); // OROC shift A0 side
1191 rocRotOROCA0=(AliTPCROCVoltError3D*)iter->FindObject("rocRotOROCA0"); // OROC rot A0 side
1192 rocShiftIROCC0=(AliTPCROCVoltError3D*)iter->FindObject("rocShiftIROCC0"); // IROC shift C0 side
1193 rocRotIROCC0=(AliTPCROCVoltError3D*)iter->FindObject("rocRotIROCC0"); // IROC rot C0 side
1194 rocShiftOROCC0=(AliTPCROCVoltError3D*)iter->FindObject("rocShiftOROCC0"); // OROC shift C0 side
1195 rocRotOROCC0=(AliTPCROCVoltError3D*)iter->FindObject("rocRotOROCC0"); // OROC rot C0 side
1196
1197 } else {
1198 corrField3DSector = new AliTPCComposedCorrection;
1199 rocShiftIROCA0=new AliTPCROCVoltError3D; // IROC shift A0 side
1200 rocRotIROCA0=new AliTPCROCVoltError3D; // IROC rot A0 side
1201 rocShiftOROCA0=new AliTPCROCVoltError3D; // OROC shift A0 side
1202 rocRotOROCA0=new AliTPCROCVoltError3D; // OROC rot A0 side
1203 rocShiftIROCC0=new AliTPCROCVoltError3D; // IROC shift C0 side
1204 rocRotIROCC0=new AliTPCROCVoltError3D; // IROC rot C0 side
1205 rocShiftOROCC0=new AliTPCROCVoltError3D; // OROC shift C0 side
1206 rocRotOROCC0=new AliTPCROCVoltError3D; // OROC rot C0 side
1207 //
1208 matrix.Zero(); matrix(0,0)=kDz;
1209 rocShiftIROCA0->SetROCData(&matrix);
1210 matrix.Zero(); matrix(0,1)=kAngle;
1211 rocRotIROCA0->SetROCData(&matrix);
1212 matrix.Zero(); matrix(36,0)=kDz;
1213 rocShiftOROCA0->SetROCData(&matrix);
1214 matrix.Zero(); matrix(36,1)=kAngle;
1215 rocRotOROCA0->SetROCData(&matrix);
1216
1217 matrix.Zero(); matrix(18,0)=kDz;
1218 rocShiftIROCC0->SetROCData(&matrix);
1219 matrix.Zero(); matrix(18,1)=kAngle;
1220 rocRotIROCC0->SetROCData(&matrix);
1221 matrix.Zero(); matrix(36+18,0)=kDz;
1222 rocShiftOROCC0->SetROCData(&matrix);
1223 matrix.Zero(); matrix(36+18,1)=kAngle;
1224 rocRotOROCC0->SetROCData(&matrix);
1225 //
1226 rocShiftIROCA0->SetElectronArrivalCorrection(kFALSE); // IROC shift A0 side
1227 rocRotIROCA0->SetElectronArrivalCorrection(kFALSE); // IROC rot A0 side
1228 rocShiftOROCA0->SetElectronArrivalCorrection(kFALSE); // OROC shift A0 side
1229 rocRotOROCA0->SetElectronArrivalCorrection(kFALSE); // OROC rot A0 side
1230 rocShiftIROCC0->SetElectronArrivalCorrection(kFALSE); // IROC shift C0 side
1231 rocRotIROCC0->SetElectronArrivalCorrection(kFALSE); // IROC rot C0 side
1232 rocShiftOROCC0->SetElectronArrivalCorrection(kFALSE); // OROC shift C0 side
1233 rocRotOROCC0->SetElectronArrivalCorrection(kFALSE); // OROC rot C0 side
1234
1235 /* // verification plots
1236 */
1237 //
1238 rocShiftIROCA0->SetName("rocShiftIROCA0");rocShiftIROCA0->SetTitle("rocShiftIROCA0");
1239 rocRotIROCA0->SetName("rocRotIROCA0");rocRotIROCA0->SetTitle("rocRotIROCA0");
1240 rocShiftOROCA0->SetName("rocShiftOROCA0"); rocShiftOROCA0->SetTitle("rocShiftOROCA0");
1241 rocRotOROCA0->SetName("rocRotOROCA0");rocRotOROCA0->SetTitle("rocRotOROCA0");
1242 //
1243 rocShiftIROCC0->SetName("rocShiftIROCC0");rocShiftIROCC0->SetTitle("rocShiftIROCC0");
1244 rocRotIROCC0->SetName("rocRotIROCC0");rocRotIROCC0->SetTitle("rocRotIROCC0");
1245 rocShiftOROCC0->SetName("rocShiftOROCC0");rocShiftOROCC0->SetTitle("rocShiftOROCC0");
1246 rocRotOROCC0->SetName("rocRotOROCC0");rocRotOROCC0->SetTitle("rocRotOROCC0");
1247 //
1248 //
1249 TObjArray *cs = new TObjArray();
1250 cs->Add(rocShiftIROCA0); // IROC shift A0 side
1251 cs->Add(rocRotIROCA0); // IROC rot A0 side
1252 cs->Add(rocShiftOROCA0); // OROC shift A0 side
1253 cs->Add(rocRotOROCA0); // OROC rot A0 side
1254 cs->Add(rocShiftIROCC0); // IROC shift C0 side
1255 cs->Add(rocRotIROCC0); // IROC rot C0 side
1256 cs->Add(rocShiftOROCC0); // OROC shift C0 side
1257 cs->Add(rocRotOROCC0); // OROC rot C0 side
1258 //
1259 corrField3DSector->SetCorrections(cs);
1260 corrField3DSector->SetOmegaTauT1T2(0,T1,T2);
1261 corrField3DSector->Print();
1262 //
1263
1264 fCorrectionsROC->cd();
1265 corrField3DSector->Init();
1266 corrField3DSector->Print("da");
1267 fCorrectionsROC->cd();
1268 corrField3DSector->Write("TPCROCVoltError3DSector");
1269 }
1270 AliTPCCorrection::AddVisualCorrection(rocShiftIROCA0,800); // IROC shift A0 side
1271 AliTPCCorrection::AddVisualCorrection(rocRotIROCA0,801); // IROC rot A0 side
1272 AliTPCCorrection::AddVisualCorrection(rocShiftOROCA0,802); // OROC shift A0 side
1273 AliTPCCorrection::AddVisualCorrection(rocRotOROCA0,803); // OROC rot A0 side
1274 AliTPCCorrection::AddVisualCorrection(rocShiftIROCC0,804); // IROC shift C0 side
1275 AliTPCCorrection::AddVisualCorrection(rocRotIROCC0,805); // IROC rot C0 side
1276 AliTPCCorrection::AddVisualCorrection(rocShiftOROCC0,806); // OROC shift C0 side
1277 AliTPCCorrection::AddVisualCorrection(rocRotOROCC0,807); // OROC rot C0 side
1278 //
1279 // Register correction for plus setting
1280 AliTPCCorrection *corr =0;
1281 Double_t wtp = -10.0 * (0.5*10) * vdrift / ezField ;
1282 Double_t wtm = -10.0 * (0.5*10) * vdrift / -ezField ;
1283
1284 corr=(AliTPCCorrection *)rocShiftIROCA0->Clone();
1285 corr->SetOmegaTauT1T2(wtp,T1,T2);
1286 AliTPCCorrection::AddVisualCorrection(corr,820); // IROC shift A0 side + Plus field
1287 //
1288 corr=(AliTPCCorrection *)rocRotIROCA0->Clone();
1289 corr->SetOmegaTauT1T2(wtp,T1,T2);
1290 AliTPCCorrection::AddVisualCorrection(corr,821); // IROC rot A0 side
1291 //
1292 corr=(AliTPCCorrection *)rocShiftOROCA0->Clone();
1293 corr->SetOmegaTauT1T2(wtp,T1,T2);
1294 AliTPCCorrection::AddVisualCorrection(corr,822); // OROC shift A0 side
1295 //
1296 corr=(AliTPCCorrection *)rocRotOROCA0->Clone();
1297 corr->SetOmegaTauT1T2(wtp,T1,T2);
1298 AliTPCCorrection::AddVisualCorrection(corr,823); // OROC rot A0 side
1299 corr=(AliTPCCorrection *)rocShiftIROCC0->Clone();
1300 corr->SetOmegaTauT1T2(wtp,T1,T2);
1301 AliTPCCorrection::AddVisualCorrection(corr,824); // IROC shift C0 side + Plus field
1302 //
1303 corr=(AliTPCCorrection *)rocRotIROCC0->Clone();
1304 corr->SetOmegaTauT1T2(wtp,T1,T2);
1305 AliTPCCorrection::AddVisualCorrection(corr,825); // IROC rot C0 side
1306 //
1307 corr=(AliTPCCorrection *)rocShiftOROCC0->Clone();
1308 corr->SetOmegaTauT1T2(wtp,T1,T2);
1309 AliTPCCorrection::AddVisualCorrection(corr,826); // OROC shift C0 side
1310 //
1311 corr=(AliTPCCorrection *)rocRotOROCC0->Clone();
1312 corr->SetOmegaTauT1T2(wtp,T1,T2);
1313 AliTPCCorrection::AddVisualCorrection(corr,827); // OROC rot C0 side
1314 //
1315 corr=(AliTPCCorrection *)rocShiftIROCA0->Clone();
1316 corr->SetOmegaTauT1T2(wtm,T1,T2);
1317 AliTPCCorrection::AddVisualCorrection(corr,840); // IROC shift A0 side + Plus field
1318 //
1319 corr=(AliTPCCorrection *)rocRotIROCA0->Clone();
1320 corr->SetOmegaTauT1T2(wtm,T1,T2);
1321 AliTPCCorrection::AddVisualCorrection(corr,841); // IROC rot A0 side
1322 //
1323 corr=(AliTPCCorrection *)rocShiftOROCA0->Clone();
1324 corr->SetOmegaTauT1T2(wtm,T1,T2);
1325 AliTPCCorrection::AddVisualCorrection(corr,842); // OROC shift A0 side
1326 //
1327 corr=(AliTPCCorrection *)rocRotOROCA0->Clone();
1328 corr->SetOmegaTauT1T2(wtm,T1,T2);
1329 AliTPCCorrection::AddVisualCorrection(corr,843); // OROC rot A0 side
1330 corr=(AliTPCCorrection *)rocShiftIROCC0->Clone();
1331 corr->SetOmegaTauT1T2(wtm,T1,T2);
1332 AliTPCCorrection::AddVisualCorrection(corr,844); // IROC shift C0 side + Plus field
1333 //
1334 corr=(AliTPCCorrection *)rocRotIROCC0->Clone();
1335 corr->SetOmegaTauT1T2(wtm,T1,T2);
1336 AliTPCCorrection::AddVisualCorrection(corr,845); // IROC rot C0 side
1337 //
1338 corr=(AliTPCCorrection *)rocShiftOROCC0->Clone();
1339 corr->SetOmegaTauT1T2(wtm,T1,T2);
1340 AliTPCCorrection::AddVisualCorrection(corr,846); // OROC shift C0 side
1341 //
1342 corr=(AliTPCCorrection *)rocRotOROCC0->Clone();
1343 corr->SetOmegaTauT1T2(wtm,T1,T2);
1344 AliTPCCorrection::AddVisualCorrection(corr,847); // OROC rot C0 side
1345 //
1346 fCorrectionsROC->Close();
1347 delete fCorrectionsROC;
1348}
1349
1350
1351
1352
1353
1354AliTPCComposedCorrection * MakeComposedCorrectionExB(){
1355 /// make composed corection for ExB scanning
1356
1357 RegisterCorrection();
1358 //
1359 //
1360 TObjArray * corr = new TObjArray; // primitive corrections - for fitting
1361 TObjArray * testCorr = new TObjArray; // test corrections - to be used as benchmark for fitting
1362 //
1363 corr->AddLast(twistX);
1364 corr->AddLast(twistY);
1365 corr->AddLast(exbShape);
1366 corr->AddLast(exbShapeT1X);
1367 corr->AddLast(exbShapeT2X);
1368 //
1369 for (Int_t i=0; i<8; i++){
1370 corr->AddLast(boundaryVoltErrorA[i]);
1371 corr->AddLast(boundaryVoltErrorC[i]);
1372 }
1373 //
1374 // ROD alignment
1375 //
1376 corr->AddLast(rodIFCCos);
1377 corr->AddLast(rodIFCSin);
1378 corr->AddLast(rodOFCSin);
1379 corr->AddLast(rodOFCCos);
1380 //alignment
1381 corr->AddLast(alignTrans0);
1382 corr->AddLast(alignTrans1);
1383 corr->AddLast(alignTrans2);
1384 corr->AddLast(alignRot0);
1385 corr->AddLast(alignRot1);
1386 corr->AddLast(alignRot2);
1387 //
1388 corr->AddLast(alignTrans0D[0]);
1389 corr->AddLast(alignTrans1D[0]);
1390 corr->AddLast(alignTrans2D[0]);
1391 corr->AddLast(alignRot0D[0]);
1392 corr->AddLast(alignRot1D[0]);
1393 corr->AddLast(alignRot2D[0]);
1394 corr->AddLast(alignTrans0D[1]);
1395 corr->AddLast(alignTrans1D[1]);
1396 corr->AddLast(alignRot0D[1]);
1397 corr->AddLast(alignTrans0D[2]);
1398 corr->AddLast(alignTrans1D[2]);
1399 corr->AddLast(alignRot0D[2]);
1400 //
1401 // z alignment + E field distortion due z misalignment
1402 //
1403 corr->AddLast(rocRotgXA); // A side C side
1404 corr->AddLast(rocRotgYA);
1405 corr->AddLast(rocRotgXC);
1406 corr->AddLast(rocRotgYC);
1407 corr->AddLast(rocDzIROCA);
1408 corr->AddLast(rocDzIROCC);
1409 corr->AddLast(rocRotIROCA);
1410 corr->AddLast(rocRotIROCC);
1411 corr->AddLast(rocDzUDA);
1412 corr->AddLast(rocDzUDC);
1413 corr->AddLast(rocRotUDA);
1414 corr->AddLast(rocRotUDC);
1415 //
1416 //
1417 corr->AddLast(calibDrift[0]);
1418 corr->AddLast(calibDrift[1]);
1419 corr->AddLast(calibDrift[2]);
1420 corr->AddLast(calibDrift[3]);
1421 corr->AddLast(calibDrift[4]);
1422 corr->AddLast(calibDrift[5]);
1423 corr->AddLast(calibDrift[6]);
1424 //
1425 // setup test correction
1426 //
1427 testCorr->AddLast(rodIFCCos);
1428 testCorr->AddLast(rodIFCSin);
1429 testCorr->AddLast(rodOFCSin);
1430 testCorr->AddLast(rodOFCCos);
1431 //
1432 testCorr->AddLast(twistX);
1433 testCorr->AddLast(twistY);
1434 testCorr->AddLast(alignTrans0);
1435 testCorr->AddLast(alignTrans1);
1436 testCorr->AddLast(alignTrans2);
1437 testCorr->AddLast(alignRot0);
1438 testCorr->AddLast(alignRot1);
1439 testCorr->AddLast(alignRot2);
1440 testCorr->AddLast(rocRotgXA); // A side C side
1441 testCorr->AddLast(rocRotgYA);
1442 testCorr->AddLast(rocRotgXC);
1443 testCorr->AddLast(rocRotgYC);
1444 testCorr->AddLast(rocDzIROCA);
1445 testCorr->AddLast(rocDzIROCC);
1446 testCorr->AddLast(rocRotIROCA);
1447 testCorr->AddLast(rocRotIROCC);
1448 Int_t entries=testCorr->GetEntries();
1449 TVectorD weights(entries);
1450 for (Int_t i=0; i<entries; i++) weights[i]=1;
1451 //
1452 AliTPCComposedCorrection *composedTest= new AliTPCComposedCorrection ;
1453 composedTest->SetName("FitSample");
1454 composedTest->SetTitle("FitSample");
1455 composedTest->SetCorrections(testCorr);
1456 composedTest->SetWeights(&weights);
1457 //
1458 corr->AddLast(composedTest);
1459 AliTPCComposedCorrection *cc= new AliTPCComposedCorrection ;
1460 cc->SetCorrections((TObjArray*)(corr));
1461 //cc->Init();
1462 cc->Print("DA"); // Print used correction classes
1463 cc->SetName("ComposedExB");
1464 TFile fexb("RegisterCorrectionExB.root","recreate");
1465 cc->Write("ComposedExB");
1466 fexb.Close();
1467 return cc;
1468}
1469
1470
1471AliTPCComposedCorrection * GetCorrectionFromFile(){
1472 /// Getthe appropariate correction form the closest file
1473
1474 TFile * fexb= TFile::Open("RegisterCorrectionExB.root");
1475 if (!fexb) fexb= TFile::Open("../RegisterCorrectionExB.root");
1476 if (!fexb) fexb= TFile::Open("../../RegisterCorrectionExB.root");
1477 if (!fexb) fexb= TFile::Open("../../../RegisterCorrectionExB.root");
1478 //
1479 if (!fexb) return 0;
1480 TFile * fitter= TFile::Open("fitCorrection.root");
1481
1482 if (!fitter) fitter= TFile::Open("../fitCorrection.root");
1483 if (!fitter) fitter= TFile::Open("../../fitCorrection.root");
1484 if (!fitter) fitter= TFile::Open("../../../fitCorrection.root");
1485 //
1486 AliTPCComposedCorrection *cc= (AliTPCComposedCorrection*) fexb->Get("ComposedExB");
1487 {if (!cc){
1488 printf("File or correction RegisterCorrectionExB.root doees not exist or corrupted\n\n\n");
1489 return 0;
1490 }}
1491 TObjArray * corr = (TObjArray*)(cc->GetCorrections());
1492 // TObjArray * corrLocal =new TObjArray;
1493 // TObjArray * corrGlobal =new TObjArray;
1494 //
1495 if (fitter){
1496 if (fitter->GetKey("FitBoundary")) corr->AddLast(fitter->Get("FitBoundary"));
1497 if (fitter->GetKey("FitExBTwist")) corr->AddLast(fitter->Get("FitExBTwist"));
1498 if (fitter->GetKey("FitAlignGlobal")) corr->AddLast(fitter->Get("FitAlignGlobal"));
1499 if (fitter->GetKey("FitRodAlignGloba")) corr->AddLast(fitter->Get("FitRodAlignGlobal"));
1500 if (fitter->GetKey("FitRocAlignGlobal")) corr->AddLast(fitter->Get("FitRocAlignGlobal"));
1501 if (fitter->GetKey("FitRocAlignZ")) corr->AddLast(fitter->Get("FitRocAlignZ"));
1502 if (fitter->GetKey("FitAlignLocal")) corr->AddLast(fitter->Get("FitAlignLocal"));
1503 if (fitter->GetKey("FitAlignTPC")) corr->AddLast(fitter->Get("FitAlignTPC"));
1504 if (fitter->GetKey("FitAlignTOF")) corr->AddLast(fitter->Get("FitAlignTOF"));
1505 if (fitter->GetKey("FitAlignTRD")) corr->AddLast(fitter->Get("FitAlignTRD"));
1506 }
1507 return cc;
1508}
1509
1510void TestParExample(){
1511 /// dz shift example: AliTPCCorrection::AddVisualCorrection(rocDzIROCA,705);
1512 /// => parabolic fit and helix fit agrees - once significant ammount of points used
1513 /// 160 point - agreement ~2%;
1514 /// 80 points - agreement ~5%
1515
1516 AliTPCCorrection* corr = AliTPCCorrection::GetVisualCorrection(705);
1517 corr->SetOmegaTauT1T2(0.33,1,1);
1518 TF1 f705Par("f705","AliTPCCorrectionFit::EvalAtPar(0,0,85, x,0.1,705,0,80)",0,360);
1519 f705Par.SetLineColor(2);f705Par.SetNpx(500);
1520 TF1 f705Helix("f705Helix","AliTPCCorrectionFit::EvalAtHelix(0,0,85,x,0.1,705,0,80)",0,360);
1521 f705Helix.SetLineColor(4);f705Helix.SetNpx(500);
1522 f705Helix.Draw();
1523 f705Par.Draw("same");
1524
1525}
1526
1527
1528
1529void TestFitSpeed(Int_t nEvals){
1530 /// test speed of helix fit/ resp. parabolic fir
1531
1532 TStopwatch timerh;
1533 ::Info("TestFitSpeed","Helix fit");
1534 for (Int_t i=0;i<nEvals; i++) AliTPCCorrectionFit::EvalAtPar(0,0,85,gRandom->Rndm(),0.1,705,0,80);
1535 timerh.Print();
1536 TStopwatch timerP;
1537 ::Info("TestFitSpeed","Parabolicfit");
1538 for (Int_t i=0;i<nEvals; i++) AliTPCCorrectionFit::EvalAtPar(0,0,85,gRandom->Rndm(),0.1,705,0,80);
1539 timerP.Print();
1540 /*
1541 For the test system CPU comsumption is the same:
1542 I-TestFitSpeed: Helix fit
1543 Real time 0:00:03, CP time 3.310
1544 I-TestFitSpeed: Parabolicfit
1545 Real time 0:00:03, CP time 3.280
1546 */
1547}
1548