1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
16 // _________________________________________________________________
19 // <h2> AliTPCCorrection class </h2>
21 // The AliTPCCorrection class provides a general framework to deal with space point distortions.
22 // An correction class which inherits from here is for example AliTPCExBBShape or AliTPCExBTwist. <br>
23 // General virtual functions are (for example) CorrectPoint(x,roc) where x is the vector of initial
24 // positions in cartesian coordinates and roc represents the read-out chamber number according to
25 // the offline numbering convention. The vector x is overwritten with the corrected coordinates. <br>
26 // An alternative usage would be CorrectPoint(x,roc,dx), which leaves the vector x untouched, but
27 // returns the distortions via the vector dx. <br>
28 // This class is normally used via the general class AliTPCComposedCorrection.
30 // Furthermore, the class contains basic geometrical descriptions like field cage radii
31 // (fgkIFCRadius, fgkOFCRadius) and length (fgkTPCZ0) plus the voltages. Also, the definitions
32 // of size and widths of the fulcrums building the grid of the final look-up table, which is
33 // then interpolated, is defined in kNX and fgkXList).
35 // All physics-model classes below are derived from this class in order to not duplicate code
36 // and to allow a uniform treatment of all physics models.
38 // <h3> Poisson solver </h3>
39 // A numerical solver of the Poisson equation (relaxation technique) is implemented for 2-dimensional
40 // geometries (r,z) as well as for 3-dimensional problems (r,$\phi$,z). The corresponding function
41 // names are PoissonRelaxation?D. The relevant function arguments are the arrays of the boundary and
42 // initial conditions (ArrayofArrayV, ArrayofChargeDensities) as well as the grid granularity which
43 // is used during the calculation. These inputs can be chosen according to the needs of the physical
44 // effect which is supposed to be simulated. In the 3D version, different symmetry conditions can be set
45 // in order to reduce the calculation time (used in AliTPCFCVoltError3D).
47 // <h3> Unified plotting functionality </h3>
48 // Generic plot functions were implemented. They return a histogram pointer in the chosen plane of
49 // the TPC drift volume with a selectable grid granularity and the magnitude of the correction vector.
50 // For example, the function CreateHistoDZinXY(z,nx,ny) returns a 2-dimensional histogram which contains
51 // the longitudinal corrections $dz$ in the (x,y)-plane at the given z position with the granularity of
52 // nx and ny. The magnitude of the corrections is defined by the class from which this function is called.
53 // In the same manner, standard plots for the (r,$\phi$)-plane and for the other corrections like $dr$ and $rd\phi$ are available
55 // Note: This class is normally used via the class AliTPCComposedCorrection
58 // Begin_Macro(source)
60 // gROOT->SetStyle("Plain"); gStyle->SetPalette(1);
61 // TCanvas *c2 = new TCanvas("cAliTPCCorrection","cAliTPCCorrection",700,1050); c2->Divide(2,3);
62 // AliTPCROCVoltError3D roc; // EXAMPLE PLOTS - SEE BELOW
63 // roc.SetROCDataFileName("$ALICE_ROOT/TPC/Calib/maps/TPCROCdzSurvey.root");
64 // roc.SetOmegaTauT1T2(0,1,1); // B=0
65 // Float_t z0 = 1; // at +1 cm -> A side
66 // c2->cd(1); roc.CreateHistoDRinXY(1.,300,300)->Draw("cont4z");
67 // c2->cd(3);roc.CreateHistoDRPhiinXY(1.,300,300)->Draw("cont4z");
68 // c2->cd(5);roc.CreateHistoDZinXY(1.,300,300)->Draw("cont4z");
70 // c2->cd(2);roc.CreateHistoDRinZR(phi0)->Draw("surf2");
71 // c2->cd(4);roc.CreateHistoDRPhiinZR(phi0)->Draw("surf2");
72 // c2->cd(6);roc.CreateHistoDZinZR(phi0)->Draw("surf2");
79 // Date: 27/04/2010 <br>
80 // Authors: Magnus Mager, Stefan Rossegger, Jim Thomas
82 // _________________________________________________________________
85 #include "Riostream.h"
90 #include <TTreeStream.h>
93 #include <TTimeStamp.h>
94 #include <AliCDBStorage.h>
96 #include <AliCDBMetaData.h>
98 #include "AliTPCParamSR.h"
100 #include "AliTPCCorrection.h"
103 #include "AliExternalTrackParam.h"
104 #include "AliTrackPointArray.h"
105 #include "TDatabasePDG.h"
106 #include "AliTrackerBase.h"
107 #include "AliTPCROC.h"
108 #include "THnSparse.h"
110 #include "AliTPCLaserTrack.h"
111 #include "AliESDVertex.h"
112 #include "AliVertexerTracks.h"
113 #include "TDatabasePDG.h"
117 #include "TDatabasePDG.h"
119 #include "AliTPCTransform.h"
120 #include "AliTPCcalibDB.h"
121 #include "AliTPCExB.h"
123 #include "AliTPCRecoParam.h"
124 #include "TLinearFitter.h"
125 #include <AliSysInfo.h>
127 ClassImp(AliTPCCorrection)
130 TObjArray *AliTPCCorrection::fgVisualCorrection=0;
131 // instance of correction for visualization
134 // FIXME: the following values should come from the database
135 const Double_t AliTPCCorrection::fgkTPCZ0 = 249.7; // nominal gating grid position
136 const Double_t AliTPCCorrection::fgkIFCRadius= 83.5; // radius which renders the "18 rod manifold" best -> compare calc. of Jim Thomas
137 // compare gkIFCRadius= 83.05: Mean Radius of the Inner Field Cage ( 82.43 min, 83.70 max) (cm)
138 const Double_t AliTPCCorrection::fgkOFCRadius= 254.5; // Mean Radius of the Outer Field Cage (252.55 min, 256.45 max) (cm)
139 const Double_t AliTPCCorrection::fgkZOffSet = 0.2; // Offset from CE: calculate all distortions closer to CE as if at this point
140 const Double_t AliTPCCorrection::fgkCathodeV = -100000.0; // Cathode Voltage (volts)
141 const Double_t AliTPCCorrection::fgkGG = -70.0; // Gating Grid voltage (volts)
143 const Double_t AliTPCCorrection::fgkdvdE = 0.0024; // [cm/V] drift velocity dependency on the E field (from Magboltz for NeCO2N2 at standard environment)
145 const Double_t AliTPCCorrection::fgkEM = -1.602176487e-19/9.10938215e-31; // charge/mass in [C/kg]
146 const Double_t AliTPCCorrection::fgke0 = 8.854187817e-12; // vacuum permittivity [A·s/(V·m)]
149 AliTPCCorrection::AliTPCCorrection()
150 : TNamed("correction_unity","unity"),fILow(0),fJLow(0),fKLow(0), fT1(1), fT2(1), fIsLocal(kFALSE)
153 // default constructor
155 if (!fgVisualCorrection) fgVisualCorrection= new TObjArray;
157 InitLookUpfulcrums();
161 AliTPCCorrection::AliTPCCorrection(const char *name,const char *title)
162 : TNamed(name,title),fILow(0),fJLow(0),fKLow(0), fT1(1), fT2(1), fIsLocal(kFALSE)
165 // default constructor, that set the name and title
167 if (!fgVisualCorrection) fgVisualCorrection= new TObjArray;
169 InitLookUpfulcrums();
173 AliTPCCorrection::~AliTPCCorrection() {
175 // virtual destructor
179 void AliTPCCorrection::CorrectPoint(Float_t x[], Short_t roc) {
181 // Corrects the initial coordinates x (cartesian coordinates)
182 // according to the given effect (inherited classes)
183 // roc represents the TPC read out chamber (offline numbering convention)
186 GetCorrection(x,roc,dx);
187 for (Int_t j=0;j<3;++j) x[j]+=dx[j];
190 void AliTPCCorrection::CorrectPoint(const Float_t x[], Short_t roc,Float_t xp[]) {
192 // Corrects the initial coordinates x (cartesian coordinates) and stores the new
193 // (distorted) coordinates in xp. The distortion is set according to the given effect (inherited classes)
194 // roc represents the TPC read out chamber (offline numbering convention)
197 GetCorrection(x,roc,dx);
198 for (Int_t j=0;j<3;++j) xp[j]=x[j]+dx[j];
201 void AliTPCCorrection::DistortPoint(Float_t x[], Short_t roc) {
203 // Distorts the initial coordinates x (cartesian coordinates)
204 // according to the given effect (inherited classes)
205 // roc represents the TPC read out chamber (offline numbering convention)
208 GetDistortion(x,roc,dx);
209 for (Int_t j=0;j<3;++j) x[j]+=dx[j];
212 void AliTPCCorrection::DistortPointLocal(Float_t x[], Short_t roc) {
214 // Distorts the initial coordinates x (cartesian coordinates)
215 // according to the given effect (inherited classes)
216 // roc represents the TPC read out chamber (offline numbering convention)
218 Float_t gxyz[3]={0,0,0};
219 Double_t alpha = TMath::TwoPi()*(roc%18+0.5)/18;
220 Double_t ca=TMath::Cos(alpha), sa= TMath::Sin(alpha);
221 gxyz[0]= ca*x[0]+sa*x[1];
222 gxyz[1]= -sa*x[0]+ca*x[1];
224 DistortPoint(gxyz,roc);
225 x[0]= ca*gxyz[0]-sa*gxyz[1];
226 x[1]= +sa*gxyz[0]+ca*gxyz[1];
229 void AliTPCCorrection::CorrectPointLocal(Float_t x[], Short_t roc) {
231 // Distorts the initial coordinates x (cartesian coordinates)
232 // according to the given effect (inherited classes)
233 // roc represents the TPC read out chamber (offline numbering convention)
235 Float_t gxyz[3]={0,0,0};
236 Double_t alpha = TMath::TwoPi()*(roc%18+0.5)/18;
237 Double_t ca=TMath::Cos(alpha), sa= TMath::Sin(alpha);
238 gxyz[0]= ca*x[0]+sa*x[1];
239 gxyz[1]= -sa*x[0]+ca*x[1];
241 CorrectPoint(gxyz,roc);
242 x[0]= ca*gxyz[0]-sa*gxyz[1];
243 x[1]= sa*gxyz[0]+ca*gxyz[1];
247 void AliTPCCorrection::DistortPoint(const Float_t x[], Short_t roc,Float_t xp[]) {
249 // Distorts the initial coordinates x (cartesian coordinates) and stores the new
250 // (distorted) coordinates in xp. The distortion is set according to the given effect (inherited classes)
251 // roc represents the TPC read out chamber (offline numbering convention)
254 GetDistortion(x,roc,dx);
255 for (Int_t j=0;j<3;++j) xp[j]=x[j]+dx[j];
258 void AliTPCCorrection::GetCorrection(const Float_t /*x*/[], Short_t /*roc*/,Float_t dx[]) {
260 // This function delivers the correction values dx in respect to the inital coordinates x
261 // roc represents the TPC read out chamber (offline numbering convention)
262 // Note: The dx is overwritten by the inherited effectice class ...
264 for (Int_t j=0;j<3;++j) { dx[j]=0.; }
267 void AliTPCCorrection::GetDistortion(const Float_t x[], Short_t roc,Float_t dx[]) {
269 // This function delivers the distortion values dx in respect to the inital coordinates x
270 // roc represents the TPC read out chamber (offline numbering convention)
272 GetCorrection(x,roc,dx);
273 for (Int_t j=0;j<3;++j) dx[j]=-dx[j];
276 void AliTPCCorrection::GetCorrectionDz(const Float_t x[], Short_t roc,Float_t dx[], Float_t delta) {
277 // author: marian.ivanov@cern.ch
279 // In this (virtual)function calculates the dx'/dz, dy'/dz and dz'/dz at given point (x,y,z)
280 // Generic implementation. Better precision can be acchieved knowing the internal structure
281 // of underlying trasnformation. Derived classes can reimplement it.
282 // To calculate correction is fitted in small neighberhood:
283 // (x+-delta,y+-delta,z+-delta) where delta is an argument
286 // x[] - space point corrdinate
287 // roc - readout chamber identifier (important e.g to do not miss the side of detector)
288 // delta - define the size of neighberhood
290 // dx[] - array {dx'/dz, dy'/dz , dz'/dz }
292 // if (fIsLocal){ //standard implemenation provides the correction/distortion integrated over full drift length
295 // GetCorrection(xyz,roc,dxyz);
297 static TLinearFitter fitx(2,"pol1");
298 static TLinearFitter fity(2,"pol1");
299 static TLinearFitter fitz(2,"pol1");
305 //adjust limits around CE to stay on one side
308 if ((x[2]+zmin*delta)<0){
320 if ((x[2]+zmax*delta)>0){
330 for (Int_t xdelta=-1; xdelta<=1; xdelta++)
331 for (Int_t ydelta=-1; ydelta<=1; ydelta++){
332 // for (Int_t zdelta=-1; zdelta<=1; zdelta++){
333 // for (Int_t xdelta=-2; xdelta<=0; xdelta++)
334 // for (Int_t ydelta=-2; ydelta<=0; ydelta++){
335 for (Int_t zdelta=zmin; zdelta<=zmax; zdelta++){
336 //TODO: what happens if x[2] is on the A-Side, but x[2]+zdelta*delta
337 // will be on the C-Side?
338 Float_t xyz[3]={x[0]+xdelta*delta, x[1]+ydelta*delta, x[2]+zdelta*delta};
340 GetCorrection(xyz,roc,dxyz);
341 Double_t adelta=zdelta*delta;
342 fitx.AddPoint(&adelta, dxyz[0]);
343 fity.AddPoint(&adelta, dxyz[1]);
344 fitz.AddPoint(&adelta, dxyz[2]);
350 dx[0] = fitx.GetParameter(1);
351 dx[1] = fity.GetParameter(1);
352 dx[2] = fitz.GetParameter(1);
355 void AliTPCCorrection::GetDistortionDz(const Float_t x[], Short_t roc,Float_t dx[], Float_t delta) {
356 // author: marian.ivanov@cern.ch
358 // In this (virtual)function calculates the dx'/dz, dy'/dz and dz'/dz at given point (x,y,z)
359 // Generic implementation. Better precision can be acchieved knowing the internal structure
360 // of underlying trasnformation. Derived classes can reimplement it.
361 // To calculate distortion is fitted in small neighberhood:
362 // (x+-delta,y+-delta,z+-delta) where delta is an argument
365 // x[] - space point corrdinate
366 // roc - readout chamber identifier (important e.g to do not miss the side of detector)
367 // delta - define the size of neighberhood
369 // dx[] - array {dx'/dz, dy'/dz , dz'/dz }
371 static TLinearFitter fitx(2,"pol1");
372 static TLinearFitter fity(2,"pol1");
373 static TLinearFitter fitz(2,"pol1");
380 //adjust limits around CE to stay on one side
383 if ((x[2]+zmin*delta)<0){
389 if ((x[2]+zmax*delta)>0){
395 //TODO: in principle one shuld check that x[2]+zdelta*delta does not get 'out of' bounds,
396 // so close to the CE it doesn't change the sign, since then the corrections will be wrong ...
397 for (Int_t xdelta=-1; xdelta<=1; xdelta++)
398 for (Int_t ydelta=-1; ydelta<=1; ydelta++){
399 for (Int_t zdelta=zmin; zdelta<=zmax; zdelta++){
400 //TODO: what happens if x[2] is on the A-Side, but x[2]+zdelta*delta
401 // will be on the C-Side?
402 //TODO: For the C-Side, does this have the correct sign?
403 Float_t xyz[3]={x[0]+xdelta*delta, x[1]+ydelta*delta, x[2]+zdelta*delta};
405 GetDistortion(xyz,roc,dxyz);
406 Double_t adelta=zdelta*delta;
407 fitx.AddPoint(&adelta, dxyz[0]);
408 fity.AddPoint(&adelta, dxyz[1]);
409 fitz.AddPoint(&adelta, dxyz[2]);
415 dx[0] = fitx.GetParameter(1);
416 dx[1] = fity.GetParameter(1);
417 dx[2] = fitz.GetParameter(1);
420 void AliTPCCorrection::GetCorrectionIntegralDz(const Float_t x[], Short_t roc,Float_t dx[], Float_t delta){
422 // Integrate 3D distortion along drift lines starting from the roc plane
423 // to the expected z position of the point, this assumes that dz is small
424 // and the error propagating to z' instead of the correct z is negligible
425 // To define the drift lines virtual function AliTPCCorrection::GetCorrectionDz is used
428 // x[] - space point corrdinate
429 // roc - readout chamber identifier (important e.g to do not miss the side of detector)
430 // delta - define the size of neighberhood
432 // dx[] - array { integral(dx'/dz), integral(dy'/dz) , integral(dz'/dz) }
434 Float_t zroc= ((roc%36)<18) ? fgkTPCZ0:-fgkTPCZ0;
435 Double_t zdrift = TMath::Abs(x[2]-zroc);
436 Int_t nsteps = Int_t(zdrift/delta)+1;
439 Float_t xyz[3]={x[0],x[1],zroc};
440 Float_t dxyz[3]={x[0],x[1],x[2]};
441 Short_t side=(roc/18)%2;
442 Float_t sign=1-2*side;
444 for (Int_t i=0;i<nsteps; i++){
445 //propagate backwards, therefore opposite signs
446 Float_t deltaZ=delta*(-sign);
447 // if (xyz[2]+deltaZ>fgkTPCZ0) deltaZ=TMath::Abs(xyz[2]-fgkTPCZ0);
448 // if (xyz[2]-deltaZ<-fgkTPCZ0) deltaZ=TMath::Abs(xyz[2]-fgkTPCZ0);
449 // protect again integrating through the CE
451 if (xyz[2]+deltaZ<0) deltaZ=-xyz[2]+1e-20;
453 if (xyz[2]+deltaZ>0) deltaZ=xyz[2]-+1e-20;
455 // since at larger drift (smaller z) the corrections are larger (absolute, but negative)
456 // the slopes will be positive.
457 // but since we chose deltaZ opposite sign the singn of the corretion should be fine
459 Float_t xyz2[3]={xyz[0],xyz[1],static_cast<Float_t>(xyz[2]+deltaZ/2.)};
460 GetCorrectionDz(xyz2,roc,dxyz,delta/2.);
461 xyz[0]+=deltaZ*dxyz[0];
462 xyz[1]+=deltaZ*dxyz[1];
464 sumdz+=deltaZ*dxyz[2];
469 dx[2]= sumdz; //TODO: is sumdz correct?
472 void AliTPCCorrection::GetDistortionIntegralDz(const Float_t x[], Short_t roc,Float_t dx[], Float_t delta){
474 // Integrate 3D distortion along drift lines
475 // To define the drift lines virtual function AliTPCCorrection::GetCorrectionDz is used
478 // x[] - space point corrdinate
479 // roc - readout chamber identifier (important e.g to do not miss the side of detector)
480 // delta - define the size of neighberhood
482 // dx[] - array { integral(dx'/dz), integral(dy'/dz) , integral(dz'/dz) }
484 Float_t zroc= ((roc%36)<18) ? fgkTPCZ0:-fgkTPCZ0;
485 Double_t zdrift = TMath::Abs(x[2]-zroc);
486 Int_t nsteps = Int_t(zdrift/delta)+1;
489 Float_t xyz[3]={x[0],x[1],x[2]};
490 Float_t dxyz[3]={x[0],x[1],x[2]};
491 Float_t sign=((roc%36)<18) ? 1.:-1.;
493 for (Int_t i=0;i<nsteps; i++){
494 Float_t deltaZ=delta;
495 if (xyz[2]+deltaZ>fgkTPCZ0) deltaZ=TMath::Abs(xyz[2]-zroc);
496 if (xyz[2]-deltaZ<-fgkTPCZ0) deltaZ=TMath::Abs(xyz[2]-zroc);
497 // since at larger drift (smaller z) the distortions are larger
498 // the slopes will be negative.
499 // and since we are moving towards the read-out plane the deltaZ for
500 // weighting the dK/dz should have the opposite sign
502 Float_t xyz2[3]={xyz[0],xyz[1],static_cast<Float_t>(xyz[2]+deltaZ/2.)};
503 GetDistortionDz(xyz2,roc,dxyz,delta/2.);
504 xyz[0]+=-deltaZ*dxyz[0];
505 xyz[1]+=-deltaZ*dxyz[1];
506 xyz[2]+=deltaZ; //TODO: Should this also be corrected for the dxyz[2]
507 sumdz+=-deltaZ*dxyz[2];
512 dx[2]= sumdz; //TODO: is sumdz correct?
517 void AliTPCCorrection::Init() {
519 // Initialization funtion (not used at the moment)
523 void AliTPCCorrection::Update(const TTimeStamp &/*timeStamp*/) {
529 void AliTPCCorrection::Print(Option_t* /*option*/) const {
531 // Print function to check which correction classes are used
532 // option=="d" prints details regarding the setted magnitude
533 // option=="a" prints the C0 and C1 coefficents for calibration purposes
535 printf("TPC spacepoint correction: \"%s\"\n",GetTitle());
538 void AliTPCCorrection:: SetOmegaTauT1T2(Float_t /*omegaTau*/,Float_t t1,Float_t t2) {
540 // Virtual funtion to pass the wt values (might become event dependent) to the inherited classes
541 // t1 and t2 represent the "effective omegaTau" corrections and were measured in a dedicated
546 //SetOmegaTauT1T2(omegaTau, t1, t2);
549 TH2F* AliTPCCorrection::CreateHistoDRinXY(Float_t z,Int_t nx,Int_t ny) {
551 // Simple plot functionality.
552 // Returns a 2d hisogram which represents the corrections in radial direction (dr)
553 // in respect to position z within the XY plane.
554 // The histogramm has nx times ny entries.
556 AliTPCParam* tpcparam = new AliTPCParamSR;
558 TH2F *h=CreateTH2F("dr_xy",GetTitle(),"x [cm]","y [cm]","dr [cm]",
559 nx,-250.,250.,ny,-250.,250.);
562 Int_t roc=z>0.?0:18; // FIXME
563 for (Int_t iy=1;iy<=ny;++iy) {
564 x[1]=h->GetYaxis()->GetBinCenter(iy);
565 for (Int_t ix=1;ix<=nx;++ix) {
566 x[0]=h->GetXaxis()->GetBinCenter(ix);
567 GetCorrection(x,roc,dx);
568 Float_t r0=TMath::Sqrt((x[0] )*(x[0] )+(x[1] )*(x[1] ));
569 if (tpcparam->GetPadRowRadii(0,0)<=r0 && r0<=tpcparam->GetPadRowRadii(36,95)) {
570 Float_t r1=TMath::Sqrt((x[0]+dx[0])*(x[0]+dx[0])+(x[1]+dx[1])*(x[1]+dx[1]));
571 h->SetBinContent(ix,iy,r1-r0);
574 h->SetBinContent(ix,iy,0.);
581 TH2F* AliTPCCorrection::CreateHistoDRPhiinXY(Float_t z,Int_t nx,Int_t ny) {
583 // Simple plot functionality.
584 // Returns a 2d hisogram which represents the corrections in rphi direction (drphi)
585 // in respect to position z within the XY plane.
586 // The histogramm has nx times ny entries.
589 AliTPCParam* tpcparam = new AliTPCParamSR;
591 TH2F *h=CreateTH2F("drphi_xy",GetTitle(),"x [cm]","y [cm]","drphi [cm]",
592 nx,-250.,250.,ny,-250.,250.);
595 Int_t roc=z>0.?0:18; // FIXME
596 for (Int_t iy=1;iy<=ny;++iy) {
597 x[1]=h->GetYaxis()->GetBinCenter(iy);
598 for (Int_t ix=1;ix<=nx;++ix) {
599 x[0]=h->GetXaxis()->GetBinCenter(ix);
600 GetCorrection(x,roc,dx);
601 Float_t r0=TMath::Sqrt((x[0] )*(x[0] )+(x[1] )*(x[1] ));
602 if (tpcparam->GetPadRowRadii(0,0)<=r0 && r0<=tpcparam->GetPadRowRadii(36,95)) {
603 Float_t phi0=TMath::ATan2(x[1] ,x[0] );
604 Float_t phi1=TMath::ATan2(x[1]+dx[1],x[0]+dx[0]);
606 Float_t dphi=phi1-phi0;
607 if (dphi<TMath::Pi()) dphi+=TMath::TwoPi();
608 if (dphi>TMath::Pi()) dphi-=TMath::TwoPi();
610 h->SetBinContent(ix,iy,r0*dphi);
613 h->SetBinContent(ix,iy,0.);
620 TH2F* AliTPCCorrection::CreateHistoDZinXY(Float_t z,Int_t nx,Int_t ny) {
622 // Simple plot functionality.
623 // Returns a 2d hisogram which represents the corrections in longitudinal direction (dz)
624 // in respect to position z within the XY plane.
625 // The histogramm has nx times ny entries.
628 AliTPCParam* tpcparam = new AliTPCParamSR;
630 TH2F *h=CreateTH2F("dz_xy",GetTitle(),"x [cm]","y [cm]","dz [cm]",
631 nx,-250.,250.,ny,-250.,250.);
634 Int_t roc=z>0.?0:18; // FIXME
635 for (Int_t iy=1;iy<=ny;++iy) {
636 x[1]=h->GetYaxis()->GetBinCenter(iy);
637 for (Int_t ix=1;ix<=nx;++ix) {
638 x[0]=h->GetXaxis()->GetBinCenter(ix);
639 GetCorrection(x,roc,dx);
640 Float_t r0=TMath::Sqrt((x[0] )*(x[0] )+(x[1] )*(x[1] ));
641 if (tpcparam->GetPadRowRadii(0,0)<=r0 && r0<=tpcparam->GetPadRowRadii(36,95)) {
642 h->SetBinContent(ix,iy,dx[2]);
645 h->SetBinContent(ix,iy,0.);
652 TH2F* AliTPCCorrection::CreateHistoDRinZR(Float_t phi,Int_t nz,Int_t nr) {
654 // Simple plot functionality.
655 // Returns a 2d hisogram which represents the corrections in r direction (dr)
656 // in respect to angle phi within the ZR plane.
657 // The histogramm has nx times ny entries.
659 TH2F *h=CreateTH2F("dr_zr",GetTitle(),"z [cm]","r [cm]","dr [cm]",
660 nz,-250.,250.,nr,85.,250.);
662 for (Int_t ir=1;ir<=nr;++ir) {
663 Float_t radius=h->GetYaxis()->GetBinCenter(ir);
664 x[0]=radius*TMath::Cos(phi);
665 x[1]=radius*TMath::Sin(phi);
666 for (Int_t iz=1;iz<=nz;++iz) {
667 x[2]=h->GetXaxis()->GetBinCenter(iz);
668 Int_t roc=x[2]>0.?0:18; // FIXME
669 GetCorrection(x,roc,dx);
670 Float_t r0=TMath::Sqrt((x[0] )*(x[0] )+(x[1] )*(x[1] ));
671 Float_t r1=TMath::Sqrt((x[0]+dx[0])*(x[0]+dx[0])+(x[1]+dx[1])*(x[1]+dx[1]));
672 h->SetBinContent(iz,ir,r1-r0);
679 TH2F* AliTPCCorrection::CreateHistoDRPhiinZR(Float_t phi,Int_t nz,Int_t nr) {
681 // Simple plot functionality.
682 // Returns a 2d hisogram which represents the corrections in rphi direction (drphi)
683 // in respect to angle phi within the ZR plane.
684 // The histogramm has nx times ny entries.
686 TH2F *h=CreateTH2F("drphi_zr",GetTitle(),"z [cm]","r [cm]","drphi [cm]",
687 nz,-250.,250.,nr,85.,250.);
689 for (Int_t iz=1;iz<=nz;++iz) {
690 x[2]=h->GetXaxis()->GetBinCenter(iz);
691 Int_t roc=x[2]>0.?0:18; // FIXME
692 for (Int_t ir=1;ir<=nr;++ir) {
693 Float_t radius=h->GetYaxis()->GetBinCenter(ir);
694 x[0]=radius*TMath::Cos(phi);
695 x[1]=radius*TMath::Sin(phi);
696 GetCorrection(x,roc,dx);
697 Float_t r0=TMath::Sqrt((x[0] )*(x[0] )+(x[1] )*(x[1] ));
698 Float_t phi0=TMath::ATan2(x[1] ,x[0] );
699 Float_t phi1=TMath::ATan2(x[1]+dx[1],x[0]+dx[0]);
701 Float_t dphi=phi1-phi0;
702 if (dphi<TMath::Pi()) dphi+=TMath::TwoPi();
703 if (dphi>TMath::Pi()) dphi-=TMath::TwoPi();
705 h->SetBinContent(iz,ir,r0*dphi);
711 TH2F* AliTPCCorrection::CreateHistoDZinZR(Float_t phi,Int_t nz,Int_t nr) {
713 // Simple plot functionality.
714 // Returns a 2d hisogram which represents the corrections in longitudinal direction (dz)
715 // in respect to angle phi within the ZR plane.
716 // The histogramm has nx times ny entries.
718 TH2F *h=CreateTH2F("dz_zr",GetTitle(),"z [cm]","r [cm]","dz [cm]",
719 nz,-250.,250.,nr,85.,250.);
721 for (Int_t ir=1;ir<=nr;++ir) {
722 Float_t radius=h->GetYaxis()->GetBinCenter(ir);
723 x[0]=radius*TMath::Cos(phi);
724 x[1]=radius*TMath::Sin(phi);
725 for (Int_t iz=1;iz<=nz;++iz) {
726 x[2]=h->GetXaxis()->GetBinCenter(iz);
727 Int_t roc=x[2]>0.?0:18; // FIXME
728 GetCorrection(x,roc,dx);
729 h->SetBinContent(iz,ir,dx[2]);
737 TH2F* AliTPCCorrection::CreateTH2F(const char *name,const char *title,
738 const char *xlabel,const char *ylabel,const char *zlabel,
739 Int_t nbinsx,Double_t xlow,Double_t xup,
740 Int_t nbinsy,Double_t ylow,Double_t yup) {
742 // Helper function to create a 2d histogramm of given size
748 while (gDirectory->FindObject(hname.Data())) {
755 TH2F *h=new TH2F(hname.Data(),title,
758 h->GetXaxis()->SetTitle(xlabel);
759 h->GetYaxis()->SetTitle(ylabel);
760 h->GetZaxis()->SetTitle(zlabel);
765 // Simple Interpolation functions: e.g. with bi(tri)cubic interpolations (not yet in TH2 and TH3)
767 void AliTPCCorrection::Interpolate2DEdistortion( Int_t order, Double_t r, Double_t z,
768 const Double_t er[kNZ][kNR], Double_t &erValue ) {
770 // Interpolate table - 2D interpolation
772 Double_t saveEr[5] = {0,0,0,0,0};
774 Search( kNZ, fgkZList, z, fJLow ) ;
775 Search( kNR, fgkRList, r, fKLow ) ;
776 if ( fJLow < 0 ) fJLow = 0 ; // check if out of range
777 if ( fKLow < 0 ) fKLow = 0 ;
778 if ( fJLow + order >= kNZ - 1 ) fJLow = kNZ - 1 - order ;
779 if ( fKLow + order >= kNR - 1 ) fKLow = kNR - 1 - order ;
781 for ( Int_t j = fJLow ; j < fJLow + order + 1 ; j++ ) {
782 saveEr[j-fJLow] = Interpolate( &fgkRList[fKLow], &er[j][fKLow], order, r ) ;
784 erValue = Interpolate( &fgkZList[fJLow], saveEr, order, z ) ;
788 void AliTPCCorrection::Interpolate3DEdistortion( Int_t order, Double_t r, Float_t phi, Double_t z,
789 const Double_t er[kNZ][kNPhi][kNR], const Double_t ephi[kNZ][kNPhi][kNR], const Double_t ez[kNZ][kNPhi][kNR],
790 Double_t &erValue, Double_t &ephiValue, Double_t &ezValue) {
792 // Interpolate table - 3D interpolation
795 Double_t saveEr[5]= {0,0,0,0,0};
796 Double_t savedEr[5]= {0,0,0,0,0} ;
798 Double_t saveEphi[5]= {0,0,0,0,0};
799 Double_t savedEphi[5]= {0,0,0,0,0} ;
801 Double_t saveEz[5]= {0,0,0,0,0};
802 Double_t savedEz[5]= {0,0,0,0,0} ;
804 Search( kNZ, fgkZList, z, fILow ) ;
805 Search( kNPhi, fgkPhiList, z, fJLow ) ;
806 Search( kNR, fgkRList, r, fKLow ) ;
808 if ( fILow < 0 ) fILow = 0 ; // check if out of range
809 if ( fJLow < 0 ) fJLow = 0 ;
810 if ( fKLow < 0 ) fKLow = 0 ;
812 if ( fILow + order >= kNZ - 1 ) fILow = kNZ - 1 - order ;
813 if ( fJLow + order >= kNPhi - 1 ) fJLow = kNPhi - 1 - order ;
814 if ( fKLow + order >= kNR - 1 ) fKLow = kNR - 1 - order ;
816 for ( Int_t i = fILow ; i < fILow + order + 1 ; i++ ) {
817 for ( Int_t j = fJLow ; j < fJLow + order + 1 ; j++ ) {
818 saveEr[j-fJLow] = Interpolate( &fgkRList[fKLow], &er[i][j][fKLow], order, r ) ;
819 saveEphi[j-fJLow] = Interpolate( &fgkRList[fKLow], &ephi[i][j][fKLow], order, r ) ;
820 saveEz[j-fJLow] = Interpolate( &fgkRList[fKLow], &ez[i][j][fKLow], order, r ) ;
822 savedEr[i-fILow] = Interpolate( &fgkPhiList[fJLow], saveEr, order, phi ) ;
823 savedEphi[i-fILow] = Interpolate( &fgkPhiList[fJLow], saveEphi, order, phi ) ;
824 savedEz[i-fILow] = Interpolate( &fgkPhiList[fJLow], saveEz, order, phi ) ;
826 erValue = Interpolate( &fgkZList[fILow], savedEr, order, z ) ;
827 ephiValue = Interpolate( &fgkZList[fILow], savedEphi, order, z ) ;
828 ezValue = Interpolate( &fgkZList[fILow], savedEz, order, z ) ;
832 Double_t AliTPCCorrection::Interpolate2DTable( Int_t order, Double_t x, Double_t y,
833 Int_t nx, Int_t ny, const Double_t xv[], const Double_t yv[],
834 const TMatrixD &array ) {
836 // Interpolate table (TMatrix format) - 2D interpolation
839 static Int_t jlow = 0, klow = 0 ;
840 Double_t saveArray[5] = {0,0,0,0,0} ;
842 Search( nx, xv, x, jlow ) ;
843 Search( ny, yv, y, klow ) ;
844 if ( jlow < 0 ) jlow = 0 ; // check if out of range
845 if ( klow < 0 ) klow = 0 ;
846 if ( jlow + order >= nx - 1 ) jlow = nx - 1 - order ;
847 if ( klow + order >= ny - 1 ) klow = ny - 1 - order ;
849 for ( Int_t j = jlow ; j < jlow + order + 1 ; j++ )
851 Double_t *ajkl = &((TMatrixD&)array)(j,klow);
852 saveArray[j-jlow] = Interpolate( &yv[klow], ajkl , order, y ) ;
855 return( Interpolate( &xv[jlow], saveArray, order, x ) ) ;
859 Double_t AliTPCCorrection::Interpolate3DTable( Int_t order, Double_t x, Double_t y, Double_t z,
860 Int_t nx, Int_t ny, Int_t nz,
861 const Double_t xv[], const Double_t yv[], const Double_t zv[],
862 TMatrixD **arrayofArrays ) {
864 // Interpolate table (TMatrix format) - 3D interpolation
867 static Int_t ilow = 0, jlow = 0, klow = 0 ;
868 Double_t saveArray[5]= {0,0,0,0,0};
869 Double_t savedArray[5]= {0,0,0,0,0} ;
871 Search( nx, xv, x, ilow ) ;
872 Search( ny, yv, y, jlow ) ;
873 Search( nz, zv, z, klow ) ;
875 if ( ilow < 0 ) ilow = 0 ; // check if out of range
876 if ( jlow < 0 ) jlow = 0 ;
877 if ( klow < 0 ) klow = 0 ;
879 if ( ilow + order >= nx - 1 ) ilow = nx - 1 - order ;
880 if ( jlow + order >= ny - 1 ) jlow = ny - 1 - order ;
881 if ( klow + order >= nz - 1 ) klow = nz - 1 - order ;
883 for ( Int_t k = klow ; k < klow + order + 1 ; k++ )
885 TMatrixD &table = *arrayofArrays[k] ;
886 for ( Int_t i = ilow ; i < ilow + order + 1 ; i++ )
888 saveArray[i-ilow] = Interpolate( &yv[jlow], &table(i,jlow), order, y ) ;
890 savedArray[k-klow] = Interpolate( &xv[ilow], saveArray, order, x ) ;
892 return( Interpolate( &zv[klow], savedArray, order, z ) ) ;
896 Double_t AliTPCCorrection::Interpolate( const Double_t xArray[], const Double_t yArray[],
897 Int_t order, Double_t x ) {
899 // Interpolate function Y(x) using linear (order=1) or quadratic (order=2) interpolation.
903 if ( order == 2 ) { // Quadratic Interpolation = 2
904 y = (x-xArray[1]) * (x-xArray[2]) * yArray[0] / ( (xArray[0]-xArray[1]) * (xArray[0]-xArray[2]) ) ;
905 y += (x-xArray[2]) * (x-xArray[0]) * yArray[1] / ( (xArray[1]-xArray[2]) * (xArray[1]-xArray[0]) ) ;
906 y += (x-xArray[0]) * (x-xArray[1]) * yArray[2] / ( (xArray[2]-xArray[0]) * (xArray[2]-xArray[1]) ) ;
907 } else { // Linear Interpolation = 1
908 y = yArray[0] + ( yArray[1]-yArray[0] ) * ( x-xArray[0] ) / ( xArray[1] - xArray[0] ) ;
915 Float_t AliTPCCorrection::Interpolate2DTable( Int_t order, Double_t x, Double_t y,
916 Int_t nx, Int_t ny, const Double_t xv[], const Double_t yv[],
917 const TMatrixF &array ) {
919 // Interpolate table (TMatrix format) - 2D interpolation
920 // Float version (in order to decrease the OCDB size)
923 static Int_t jlow = 0, klow = 0 ;
924 Float_t saveArray[5] = {0.,0.,0.,0.,0.} ;
926 Search( nx, xv, x, jlow ) ;
927 Search( ny, yv, y, klow ) ;
928 if ( jlow < 0 ) jlow = 0 ; // check if out of range
929 if ( klow < 0 ) klow = 0 ;
930 if ( jlow + order >= nx - 1 ) jlow = nx - 1 - order ;
931 if ( klow + order >= ny - 1 ) klow = ny - 1 - order ;
933 for ( Int_t j = jlow ; j < jlow + order + 1 ; j++ )
935 Float_t *ajkl = &((TMatrixF&)array)(j,klow);
936 saveArray[j-jlow] = Interpolate( &yv[klow], ajkl , order, y ) ;
939 return( Interpolate( &xv[jlow], saveArray, order, x ) ) ;
943 Float_t AliTPCCorrection::Interpolate3DTable( Int_t order, Double_t x, Double_t y, Double_t z,
944 Int_t nx, Int_t ny, Int_t nz,
945 const Double_t xv[], const Double_t yv[], const Double_t zv[],
946 TMatrixF **arrayofArrays ) {
948 // Interpolate table (TMatrix format) - 3D interpolation
949 // Float version (in order to decrease the OCDB size)
952 static Int_t ilow = 0, jlow = 0, klow = 0 ;
953 Float_t saveArray[5]= {0.,0.,0.,0.,0.};
954 Float_t savedArray[5]= {0.,0.,0.,0.,0.} ;
956 Search( nx, xv, x, ilow ) ;
957 Search( ny, yv, y, jlow ) ;
958 Search( nz, zv, z, klow ) ;
960 if ( ilow < 0 ) ilow = 0 ; // check if out of range
961 if ( jlow < 0 ) jlow = 0 ;
962 if ( klow < 0 ) klow = 0 ;
964 if ( ilow + order >= nx - 1 ) ilow = nx - 1 - order ;
965 if ( jlow + order >= ny - 1 ) jlow = ny - 1 - order ;
966 if ( klow + order >= nz - 1 ) klow = nz - 1 - order ;
968 for ( Int_t k = klow ; k < klow + order + 1 ; k++ )
970 TMatrixF &table = *arrayofArrays[k] ;
971 for ( Int_t i = ilow ; i < ilow + order + 1 ; i++ )
973 saveArray[i-ilow] = Interpolate( &yv[jlow], &table(i,jlow), order, y ) ;
975 savedArray[k-klow] = Interpolate( &xv[ilow], saveArray, order, x ) ;
977 return( Interpolate( &zv[klow], savedArray, order, z ) ) ;
980 Float_t AliTPCCorrection::Interpolate( const Double_t xArray[], const Float_t yArray[],
981 Int_t order, Double_t x ) {
983 // Interpolate function Y(x) using linear (order=1) or quadratic (order=2) interpolation.
984 // Float version (in order to decrease the OCDB size)
988 if ( order == 2 ) { // Quadratic Interpolation = 2
989 y = (x-xArray[1]) * (x-xArray[2]) * yArray[0] / ( (xArray[0]-xArray[1]) * (xArray[0]-xArray[2]) ) ;
990 y += (x-xArray[2]) * (x-xArray[0]) * yArray[1] / ( (xArray[1]-xArray[2]) * (xArray[1]-xArray[0]) ) ;
991 y += (x-xArray[0]) * (x-xArray[1]) * yArray[2] / ( (xArray[2]-xArray[0]) * (xArray[2]-xArray[1]) ) ;
992 } else { // Linear Interpolation = 1
993 y = yArray[0] + ( yArray[1]-yArray[0] ) * ( x-xArray[0] ) / ( xArray[1] - xArray[0] ) ;
1002 void AliTPCCorrection::Search( Int_t n, const Double_t xArray[], Double_t x, Int_t &low ) {
1004 // Search an ordered table by starting at the most recently used point
1007 Long_t middle, high ;
1008 Int_t ascend = 0, increment = 1 ;
1010 if ( xArray[n-1] >= xArray[0] ) ascend = 1 ; // Ascending ordered table if true
1012 if ( low < 0 || low > n-1 ) {
1013 low = -1 ; high = n ;
1014 } else { // Ordered Search phase
1015 if ( (Int_t)( x >= xArray[low] ) == ascend ) {
1016 if ( low == n-1 ) return ;
1018 while ( (Int_t)( x >= xArray[high] ) == ascend ) {
1021 high = low + increment ;
1022 if ( high > n-1 ) { high = n ; break ; }
1025 if ( low == 0 ) { low = -1 ; return ; }
1027 while ( (Int_t)( x < xArray[low] ) == ascend ) {
1030 if ( increment >= high ) { low = -1 ; break ; }
1031 else low = high - increment ;
1036 while ( (high-low) != 1 ) { // Binary Search Phase
1037 middle = ( high + low ) / 2 ;
1038 if ( (Int_t)( x >= xArray[middle] ) == ascend )
1044 if ( x == xArray[n-1] ) low = n-2 ;
1045 if ( x == xArray[0] ) low = 0 ;
1049 void AliTPCCorrection::InitLookUpfulcrums() {
1051 // Initialization of interpolation points - for main look up table
1052 // (course grid in the middle, fine grid on the borders)
1055 AliTPCROC * roc = AliTPCROC::Instance();
1056 const Double_t rLow = TMath::Floor(roc->GetPadRowRadii(0,0))-1; // first padRow plus some margin
1060 for (Int_t i = 1; i<kNR; i++) {
1061 fgkRList[i] = fgkRList[i-1] + 3.5; // 3.5 cm spacing
1062 if (fgkRList[i]<90 ||fgkRList[i]>245)
1063 fgkRList[i] = fgkRList[i-1] + 0.5; // 0.5 cm spacing
1064 else if (fgkRList[i]<100 || fgkRList[i]>235)
1065 fgkRList[i] = fgkRList[i-1] + 1.5; // 1.5 cm spacing
1066 else if (fgkRList[i]<120 || fgkRList[i]>225)
1067 fgkRList[i] = fgkRList[i-1] + 2.5; // 2.5 cm spacing
1071 fgkZList[0] = -249.5;
1072 fgkZList[kNZ-1] = 249.5;
1073 for (Int_t j = 1; j<kNZ/2; j++) {
1074 fgkZList[j] = fgkZList[j-1];
1075 if (TMath::Abs(fgkZList[j])< 0.15)
1076 fgkZList[j] = fgkZList[j-1] + 0.09; // 0.09 cm spacing
1077 else if(TMath::Abs(fgkZList[j])< 0.6)
1078 fgkZList[j] = fgkZList[j-1] + 0.4; // 0.4 cm spacing
1079 else if (TMath::Abs(fgkZList[j])< 2.5 || TMath::Abs(fgkZList[j])>248)
1080 fgkZList[j] = fgkZList[j-1] + 0.5; // 0.5 cm spacing
1081 else if (TMath::Abs(fgkZList[j])<10 || TMath::Abs(fgkZList[j])>235)
1082 fgkZList[j] = fgkZList[j-1] + 1.5; // 1.5 cm spacing
1083 else if (TMath::Abs(fgkZList[j])<25 || TMath::Abs(fgkZList[j])>225)
1084 fgkZList[j] = fgkZList[j-1] + 2.5; // 2.5 cm spacing
1086 fgkZList[j] = fgkZList[j-1] + 4; // 4 cm spacing
1088 fgkZList[kNZ-j-1] = -fgkZList[j];
1092 for (Int_t k = 0; k<kNPhi; k++)
1093 fgkPhiList[k] = TMath::TwoPi()*k/(kNPhi-1);
1099 void AliTPCCorrection::PoissonRelaxation2D(TMatrixD &arrayV, TMatrixD &chargeDensity,
1100 TMatrixD &arrayErOverEz, TMatrixD &arrayDeltaEz,
1101 Int_t rows, Int_t columns, Int_t iterations,
1102 Bool_t rocDisplacement ) {
1104 // Solve Poisson's Equation by Relaxation Technique in 2D (assuming cylindrical symmetry)
1106 // Solve Poissons equation in a cylindrical coordinate system. The arrayV matrix must be filled with the
1107 // boundary conditions on the first and last rows, and the first and last columns. The remainder of the
1108 // array can be blank or contain a preliminary guess at the solution. The Charge density matrix contains
1109 // the enclosed spacecharge density at each point. The charge density matrix can be full of zero's if
1110 // you wish to solve Laplaces equation however it should not contain random numbers or you will get
1111 // random numbers back as a solution.
1112 // Poisson's equation is solved by iteratively relaxing the matrix to the final solution. In order to
1113 // speed up the convergence to the best solution, this algorithm does a binary expansion of the solution
1114 // space. First it solves the problem on a very sparse grid by skipping rows and columns in the original
1115 // matrix. Then it doubles the number of points and solves the problem again. Then it doubles the
1116 // number of points and solves the problem again. This happens several times until the maximum number
1117 // of points has been included in the array.
1119 // NOTE: In order for this algorithmto work, the number of rows and columns must be a power of 2 plus one.
1120 // So rows == 2**M + 1 and columns == 2**N + 1. The number of rows and columns can be different.
1122 // NOTE: rocDisplacement is used to include (or ignore) the ROC misalignment in the dz calculation
1124 // Original code by Jim Thomas (STAR TPC Collaboration)
1127 Double_t ezField = (fgkCathodeV-fgkGG)/fgkTPCZ0; // = ALICE Electric Field (V/cm) Magnitude ~ -400 V/cm;
1129 const Float_t gridSizeR = (fgkOFCRadius-fgkIFCRadius) / (rows-1) ;
1130 const Float_t gridSizeZ = fgkTPCZ0 / (columns-1) ;
1131 const Float_t ratio = gridSizeR*gridSizeR / (gridSizeZ*gridSizeZ) ;
1133 TMatrixD arrayEr(rows,columns) ;
1134 TMatrixD arrayEz(rows,columns) ;
1136 //Check that number of rows and columns is suitable for a binary expansion
1138 if ( !IsPowerOfTwo(rows-1) ) {
1139 AliError("PoissonRelaxation - Error in the number of rows. Must be 2**M - 1");
1142 if ( !IsPowerOfTwo(columns-1) ) {
1143 AliError("PoissonRelaxation - Error in the number of columns. Must be 2**N - 1");
1147 // Solve Poisson's equation in cylindrical coordinates by relaxation technique
1148 // Allow for different size grid spacing in R and Z directions
1149 // Use a binary expansion of the size of the matrix to speed up the solution of the problem
1151 Int_t iOne = (rows-1)/4 ;
1152 Int_t jOne = (columns-1)/4 ;
1153 // Solve for N in 2**N, add one.
1154 Int_t loops = 1 + (int) ( 0.5 + TMath::Log2( (double) TMath::Max(iOne,jOne) ) ) ;
1156 for ( Int_t count = 0 ; count < loops ; count++ ) {
1157 // Loop while the matrix expands & the resolution increases.
1159 Float_t tempGridSizeR = gridSizeR * iOne ;
1160 Float_t tempRatio = ratio * iOne * iOne / ( jOne * jOne ) ;
1161 Float_t tempFourth = 1.0 / (2.0 + 2.0*tempRatio) ;
1163 // Do this the standard C++ way to avoid gcc extensions for Float_t coef1[rows]
1164 std::vector<float> coef1(rows) ;
1165 std::vector<float> coef2(rows) ;
1167 for ( Int_t i = iOne ; i < rows-1 ; i+=iOne ) {
1168 Float_t radius = fgkIFCRadius + i*gridSizeR ;
1169 coef1[i] = 1.0 + tempGridSizeR/(2*radius);
1170 coef2[i] = 1.0 - tempGridSizeR/(2*radius);
1173 TMatrixD sumChargeDensity(rows,columns) ;
1175 for ( Int_t i = iOne ; i < rows-1 ; i += iOne ) {
1176 Float_t radius = fgkIFCRadius + iOne*gridSizeR ;
1177 for ( Int_t j = jOne ; j < columns-1 ; j += jOne ) {
1178 if ( iOne == 1 && jOne == 1 ) sumChargeDensity(i,j) = chargeDensity(i,j) ;
1180 // Add up all enclosed charge density contributions within 1/2 unit in all directions
1181 Float_t weight = 0.0 ;
1183 sumChargeDensity(i,j) = 0.0 ;
1184 for ( Int_t ii = i-iOne/2 ; ii <= i+iOne/2 ; ii++ ) {
1185 for ( Int_t jj = j-jOne/2 ; jj <= j+jOne/2 ; jj++ ) {
1186 if ( ii == i-iOne/2 || ii == i+iOne/2 || jj == j-jOne/2 || jj == j+jOne/2 ) weight = 0.5 ;
1189 // Note that this is cylindrical geometry
1190 sumChargeDensity(i,j) += chargeDensity(ii,jj)*weight*radius ;
1191 sum += weight*radius ;
1194 sumChargeDensity(i,j) /= sum ;
1196 sumChargeDensity(i,j) *= tempGridSizeR*tempGridSizeR; // just saving a step later on
1200 for ( Int_t k = 1 ; k <= iterations; k++ ) {
1201 // Solve Poisson's Equation
1202 // Over-relaxation index, must be >= 1 but < 2. Arrange for it to evolve from 2 => 1
1203 // as interations increase.
1204 Float_t overRelax = 1.0 + TMath::Sqrt( TMath::Cos( (k*TMath::PiOver2())/iterations ) ) ;
1205 Float_t overRelaxM1 = overRelax - 1.0 ;
1206 Float_t overRelaxtempFourth, overRelaxcoef5 ;
1207 overRelaxtempFourth = overRelax * tempFourth ;
1208 overRelaxcoef5 = overRelaxM1 / overRelaxtempFourth ;
1210 for ( Int_t i = iOne ; i < rows-1 ; i += iOne ) {
1211 for ( Int_t j = jOne ; j < columns-1 ; j += jOne ) {
1213 arrayV(i,j) = ( coef2[i] * arrayV(i-iOne,j)
1214 + tempRatio * ( arrayV(i,j-jOne) + arrayV(i,j+jOne) )
1215 - overRelaxcoef5 * arrayV(i,j)
1216 + coef1[i] * arrayV(i+iOne,j)
1217 + sumChargeDensity(i,j)
1218 ) * overRelaxtempFourth;
1222 if ( k == iterations ) {
1223 // After full solution is achieved, copy low resolution solution into higher res array
1224 for ( Int_t i = iOne ; i < rows-1 ; i += iOne ) {
1225 for ( Int_t j = jOne ; j < columns-1 ; j += jOne ) {
1228 arrayV(i+iOne/2,j) = ( arrayV(i+iOne,j) + arrayV(i,j) ) / 2 ;
1229 if ( i == iOne ) arrayV(i-iOne/2,j) = ( arrayV(0,j) + arrayV(iOne,j) ) / 2 ;
1232 arrayV(i,j+jOne/2) = ( arrayV(i,j+jOne) + arrayV(i,j) ) / 2 ;
1233 if ( j == jOne ) arrayV(i,j-jOne/2) = ( arrayV(i,0) + arrayV(i,jOne) ) / 2 ;
1235 if ( iOne > 1 && jOne > 1 ) {
1236 arrayV(i+iOne/2,j+jOne/2) = ( arrayV(i+iOne,j+jOne) + arrayV(i,j) ) / 2 ;
1237 if ( i == iOne ) arrayV(i-iOne/2,j-jOne/2) = ( arrayV(0,j-jOne) + arrayV(iOne,j) ) / 2 ;
1238 if ( j == jOne ) arrayV(i-iOne/2,j-jOne/2) = ( arrayV(i-iOne,0) + arrayV(i,jOne) ) / 2 ;
1239 // Note that this leaves a point at the upper left and lower right corners uninitialized.
1240 // -> Not a big deal.
1249 iOne = iOne / 2 ; if ( iOne < 1 ) iOne = 1 ;
1250 jOne = jOne / 2 ; if ( jOne < 1 ) jOne = 1 ;
1252 sumChargeDensity.Clear();
1255 // Differentiate V(r) and solve for E(r) using special equations for the first and last rows
1256 for ( Int_t j = 0 ; j < columns ; j++ ) {
1257 for ( Int_t i = 1 ; i < rows-1 ; i++ ) arrayEr(i,j) = -1 * ( arrayV(i+1,j) - arrayV(i-1,j) ) / (2*gridSizeR) ;
1258 arrayEr(0,j) = -1 * ( -0.5*arrayV(2,j) + 2.0*arrayV(1,j) - 1.5*arrayV(0,j) ) / gridSizeR ;
1259 arrayEr(rows-1,j) = -1 * ( 1.5*arrayV(rows-1,j) - 2.0*arrayV(rows-2,j) + 0.5*arrayV(rows-3,j) ) / gridSizeR ;
1262 // Differentiate V(z) and solve for E(z) using special equations for the first and last columns
1263 for ( Int_t i = 0 ; i < rows ; i++) {
1264 for ( Int_t j = 1 ; j < columns-1 ; j++ ) arrayEz(i,j) = -1 * ( arrayV(i,j+1) - arrayV(i,j-1) ) / (2*gridSizeZ) ;
1265 arrayEz(i,0) = -1 * ( -0.5*arrayV(i,2) + 2.0*arrayV(i,1) - 1.5*arrayV(i,0) ) / gridSizeZ ;
1266 arrayEz(i,columns-1) = -1 * ( 1.5*arrayV(i,columns-1) - 2.0*arrayV(i,columns-2) + 0.5*arrayV(i,columns-3) ) / gridSizeZ ;
1269 for ( Int_t i = 0 ; i < rows ; i++) {
1270 // Note: go back and compare to old version of this code. See notes below.
1271 // JT Test ... attempt to divide by real Ez not Ez to first order
1272 for ( Int_t j = 0 ; j < columns ; j++ ) {
1273 arrayEz(i,j) += ezField;
1274 // This adds back the overall Z gradient of the field (main E field component)
1276 // Warning: (-=) assumes you are using an error potetial without the overall Field included
1279 // Integrate Er/Ez from Z to zero
1280 for ( Int_t j = 0 ; j < columns ; j++ ) {
1281 for ( Int_t i = 0 ; i < rows ; i++ ) {
1283 Int_t index = 1 ; // Simpsons rule if N=odd. If N!=odd then add extra point by trapezoidal rule.
1284 arrayErOverEz(i,j) = 0.0 ;
1285 arrayDeltaEz(i,j) = 0.0 ;
1287 for ( Int_t k = j ; k < columns ; k++ ) {
1288 arrayErOverEz(i,j) += index*(gridSizeZ/3.0)*arrayEr(i,k)/arrayEz(i,k) ;
1289 arrayDeltaEz(i,j) += index*(gridSizeZ/3.0)*(arrayEz(i,k)-ezField) ;
1290 if ( index != 4 ) index = 4; else index = 2 ;
1293 arrayErOverEz(i,j) -= (gridSizeZ/3.0)*arrayEr(i,columns-1)/arrayEz(i,columns-1) ;
1294 arrayDeltaEz(i,j) -= (gridSizeZ/3.0)*(arrayEz(i,columns-1)-ezField) ;
1297 arrayErOverEz(i,j) += (gridSizeZ/3.0) * ( 0.5*arrayEr(i,columns-2)/arrayEz(i,columns-2)
1298 -2.5*arrayEr(i,columns-1)/arrayEz(i,columns-1));
1299 arrayDeltaEz(i,j) += (gridSizeZ/3.0) * ( 0.5*(arrayEz(i,columns-2)-ezField)
1300 -2.5*(arrayEz(i,columns-1)-ezField));
1302 if ( j == columns-2 ) {
1303 arrayErOverEz(i,j) = (gridSizeZ/3.0) * ( 1.5*arrayEr(i,columns-2)/arrayEz(i,columns-2)
1304 +1.5*arrayEr(i,columns-1)/arrayEz(i,columns-1) ) ;
1305 arrayDeltaEz(i,j) = (gridSizeZ/3.0) * ( 1.5*(arrayEz(i,columns-2)-ezField)
1306 +1.5*(arrayEz(i,columns-1)-ezField) ) ;
1308 if ( j == columns-1 ) {
1309 arrayErOverEz(i,j) = 0.0 ;
1310 arrayDeltaEz(i,j) = 0.0 ;
1315 // calculate z distortion from the integrated Delta Ez residuals
1316 // and include the aquivalence (Volt to cm) of the ROC shift !!
1318 for ( Int_t j = 0 ; j < columns ; j++ ) {
1319 for ( Int_t i = 0 ; i < rows ; i++ ) {
1321 // Scale the Ez distortions with the drift velocity pertubation -> delivers cm
1322 arrayDeltaEz(i,j) = arrayDeltaEz(i,j)*fgkdvdE;
1324 // ROC Potential in cm aquivalent
1325 Double_t dzROCShift = arrayV(i, columns -1)/ezField;
1326 if ( rocDisplacement ) arrayDeltaEz(i,j) = arrayDeltaEz(i,j) + dzROCShift; // add the ROC misaligment
1336 void AliTPCCorrection::PoissonRelaxation3D( TMatrixD**arrayofArrayV, TMatrixD**arrayofChargeDensities,
1337 TMatrixD**arrayofEroverEz, TMatrixD**arrayofEPhioverEz, TMatrixD**arrayofDeltaEz,
1338 Int_t rows, Int_t columns, Int_t phislices,
1339 Float_t deltaphi, Int_t iterations, Int_t symmetry,
1340 Bool_t rocDisplacement ) {
1342 // 3D - Solve Poisson's Equation in 3D by Relaxation Technique
1344 // NOTE: In order for this algorith to work, the number of rows and columns must be a power of 2 plus one.
1345 // The number of rows and COLUMNS can be different.
1348 // COLUMNS == 2**N + 1
1349 // PHISLICES == Arbitrary but greater than 3
1351 // DeltaPhi in Radians
1353 // SYMMETRY = 0 if no phi symmetries, and no phi boundary conditions
1354 // = 1 if we have reflection symmetry at the boundaries (eg. sector symmetry or half sector symmetries).
1356 // NOTE: rocDisplacement is used to include (or ignore) the ROC misalignment in the dz calculation
1358 const Double_t ezField = (fgkCathodeV-fgkGG)/fgkTPCZ0; // = ALICE Electric Field (V/cm) Magnitude ~ -400 V/cm;
1360 const Float_t gridSizeR = (fgkOFCRadius-fgkIFCRadius) / (rows-1) ;
1361 const Float_t gridSizePhi = deltaphi ;
1362 const Float_t gridSizeZ = fgkTPCZ0 / (columns-1) ;
1363 const Float_t ratioPhi = gridSizeR*gridSizeR / (gridSizePhi*gridSizePhi) ;
1364 const Float_t ratioZ = gridSizeR*gridSizeR / (gridSizeZ*gridSizeZ) ;
1366 TMatrixD arrayE(rows,columns) ;
1368 // Check that the number of rows and columns is suitable for a binary expansion
1369 if ( !IsPowerOfTwo((rows-1)) ) {
1370 AliError("Poisson3DRelaxation - Error in the number of rows. Must be 2**M - 1");
1372 if ( !IsPowerOfTwo((columns-1)) ) {
1373 AliError("Poisson3DRelaxation - Error in the number of columns. Must be 2**N - 1");
1375 if ( phislices <= 3 ) {
1376 AliError("Poisson3DRelaxation - Error in the number of phislices. Must be larger than 3");
1378 if ( phislices > 1000 ) {
1379 AliError("Poisson3D phislices > 1000 is not allowed (nor wise) ");
1382 // Solve Poisson's equation in cylindrical coordinates by relaxation technique
1383 // Allow for different size grid spacing in R and Z directions
1384 // Use a binary expansion of the matrix to speed up the solution of the problem
1386 Int_t loops, mplus, mminus, signplus, signminus ;
1387 Int_t ione = (rows-1)/4 ;
1388 Int_t jone = (columns-1)/4 ;
1389 loops = TMath::Max(ione, jone) ; // Calculate the number of loops for the binary expansion
1390 loops = 1 + (int) ( 0.5 + TMath::Log2((double)loops) ) ; // Solve for N in 2**N
1392 TMatrixD* arrayofSumChargeDensities[1000] ; // Create temporary arrays to store low resolution charge arrays
1394 for ( Int_t i = 0 ; i < phislices ; i++ ) { arrayofSumChargeDensities[i] = new TMatrixD(rows,columns) ; }
1395 AliSysInfo::AddStamp("3DInit", 10,0,0);
1397 for ( Int_t count = 0 ; count < loops ; count++ ) { // START the master loop and do the binary expansion
1398 AliSysInfo::AddStamp("3Diter", 20,count,0);
1400 Float_t tempgridSizeR = gridSizeR * ione ;
1401 Float_t tempratioPhi = ratioPhi * ione * ione ; // Used tobe divided by ( m_one * m_one ) when m_one was != 1
1402 Float_t tempratioZ = ratioZ * ione * ione / ( jone * jone ) ;
1404 std::vector<float> coef1(rows) ; // Do this the standard C++ way to avoid gcc extensions for Float_t coef1[rows]
1405 std::vector<float> coef2(rows) ; // Do this the standard C++ way to avoid gcc extensions for Float_t coef1[rows]
1406 std::vector<float> coef3(rows) ; // Do this the standard C++ way to avoid gcc extensions for Float_t coef1[rows]
1407 std::vector<float> coef4(rows) ; // Do this the standard C++ way to avoid gcc extensions for Float_t coef1[rows]
1409 for ( Int_t i = ione ; i < rows-1 ; i+=ione ) {
1410 Float_t radius = fgkIFCRadius + i*gridSizeR ;
1411 coef1[i] = 1.0 + tempgridSizeR/(2*radius);
1412 coef2[i] = 1.0 - tempgridSizeR/(2*radius);
1413 coef3[i] = tempratioPhi/(radius*radius);
1414 coef4[i] = 0.5 / (1.0 + tempratioZ + coef3[i]);
1417 for ( Int_t m = 0 ; m < phislices ; m++ ) {
1418 TMatrixD &chargeDensity = *arrayofChargeDensities[m] ;
1419 TMatrixD &sumChargeDensity = *arrayofSumChargeDensities[m] ;
1420 for ( Int_t i = ione ; i < rows-1 ; i += ione ) {
1421 Float_t radius = fgkIFCRadius + i*gridSizeR ;
1422 for ( Int_t j = jone ; j < columns-1 ; j += jone ) {
1423 if ( ione == 1 && jone == 1 ) sumChargeDensity(i,j) = chargeDensity(i,j) ;
1424 else { // Add up all enclosed charge density contributions within 1/2 unit in all directions
1425 Float_t weight = 0.0 ;
1427 sumChargeDensity(i,j) = 0.0 ;
1428 for ( Int_t ii = i-ione/2 ; ii <= i+ione/2 ; ii++ ) {
1429 for ( Int_t jj = j-jone/2 ; jj <= j+jone/2 ; jj++ ) {
1430 if ( ii == i-ione/2 || ii == i+ione/2 || jj == j-jone/2 || jj == j+jone/2 ) weight = 0.5 ;
1433 sumChargeDensity(i,j) += chargeDensity(ii,jj)*weight*radius ;
1434 sum += weight*radius ;
1437 sumChargeDensity(i,j) /= sum ;
1439 sumChargeDensity(i,j) *= tempgridSizeR*tempgridSizeR; // just saving a step later on
1444 for ( Int_t k = 1 ; k <= iterations; k++ ) {
1446 // over-relaxation index, >= 1 but < 2
1447 Float_t overRelax = 1.0 + TMath::Sqrt( TMath::Cos( (k*TMath::PiOver2())/iterations ) ) ;
1448 Float_t overRelaxM1 = overRelax - 1.0 ;
1450 std::vector<float> overRelaxcoef4(rows) ; // Do this the standard C++ way to avoid gcc extensions
1451 std::vector<float> overRelaxcoef5(rows) ; // Do this the standard C++ way to avoid gcc extensions
1453 for ( Int_t i = ione ; i < rows-1 ; i+=ione ) {
1454 overRelaxcoef4[i] = overRelax * coef4[i] ;
1455 overRelaxcoef5[i] = overRelaxM1 / overRelaxcoef4[i] ;
1458 for ( Int_t m = 0 ; m < phislices ; m++ ) {
1460 mplus = m + 1; signplus = 1 ;
1461 mminus = m - 1 ; signminus = 1 ;
1462 if (symmetry==1) { // Reflection symmetry in phi (e.g. symmetry at sector boundaries, or half sectors, etc.)
1463 if ( mplus > phislices-1 ) mplus = phislices - 2 ;
1464 if ( mminus < 0 ) mminus = 1 ;
1466 else if (symmetry==-1) { // Anti-symmetry in phi
1467 if ( mplus > phislices-1 ) { mplus = phislices - 2 ; signplus = -1 ; }
1468 if ( mminus < 0 ) { mminus = 1 ; signminus = -1 ; }
1470 else { // No Symmetries in phi, no boundaries, the calculation is continuous across all phi
1471 if ( mplus > phislices-1 ) mplus = m + 1 - phislices ;
1472 if ( mminus < 0 ) mminus = m - 1 + phislices ;
1474 TMatrixD& arrayV = *arrayofArrayV[m] ;
1475 TMatrixD& arrayVP = *arrayofArrayV[mplus] ;
1476 TMatrixD& arrayVM = *arrayofArrayV[mminus] ;
1477 TMatrixD& sumChargeDensity = *arrayofSumChargeDensities[m] ;
1478 Double_t *arrayVfast = arrayV.GetMatrixArray();
1479 Double_t *arrayVPfast = arrayVP.GetMatrixArray();
1480 Double_t *arrayVMfast = arrayVM.GetMatrixArray();
1481 Double_t *sumChargeDensityFast=sumChargeDensity.GetMatrixArray();
1484 // slow implementation
1485 for ( Int_t i = ione ; i < rows-1 ; i+=ione ) {
1486 for ( Int_t j = jone ; j < columns-1 ; j+=jone ) {
1488 arrayV(i,j) = ( coef2[i] * arrayV(i-ione,j)
1489 + tempratioZ * ( arrayV(i,j-jone) + arrayV(i,j+jone) )
1490 - overRelaxcoef5[i] * arrayV(i,j)
1491 + coef1[i] * arrayV(i+ione,j)
1492 + coef3[i] * ( signplus*arrayVP(i,j) + signminus*arrayVM(i,j) )
1493 + sumChargeDensity(i,j)
1494 ) * overRelaxcoef4[i] ;
1495 // Note: over-relax the solution at each step. This speeds up the convergance.
1499 for ( Int_t i = ione ; i < rows-1 ; i+=ione ) {
1500 Double_t *arrayVfastI = &(arrayVfast[i*columns]);
1501 Double_t *arrayVPfastI = &(arrayVPfast[i*columns]);
1502 Double_t *arrayVMfastI = &(arrayVMfast[i*columns]);
1503 Double_t *sumChargeDensityFastI=&(sumChargeDensityFast[i*columns]);
1504 for ( Int_t j = jone ; j < columns-1 ; j+=jone ) {
1505 Double_t /*resSlow*/resFast;
1506 // resSlow = ( coef2[i] * arrayV(i-ione,j)
1507 // + tempratioZ * ( arrayV(i,j-jone) + arrayV(i,j+jone) )
1508 // - overRelaxcoef5[i] * arrayV(i,j)
1509 // + coef1[i] * arrayV(i+ione,j)
1510 // + coef3[i] * ( signplus*arrayVP(i,j) + signminus*arrayVM(i,j) )
1511 // + sumChargeDensity(i,j)
1512 // ) * overRelaxcoef4[i] ;
1513 resFast = ( coef2[i] * arrayVfastI[j-columns*ione]
1514 + tempratioZ * ( arrayVfastI[j-jone] + arrayVfastI[j+jone] )
1515 - overRelaxcoef5[i] * arrayVfastI[j]
1516 + coef1[i] * arrayVfastI[j+columns*ione]
1517 + coef3[i] * ( signplus* arrayVPfastI[j] + signminus*arrayVMfastI[j])
1518 + sumChargeDensityFastI[j]
1519 ) * overRelaxcoef4[i] ;
1520 // if (resSlow!=resFast){
1521 // printf("problem\t%d\t%d\t%f\t%f\t%f\n",i,j,resFast,resSlow,resFast-resSlow);
1523 arrayVfastI[j]=resFast;
1524 // Note: over-relax the solution at each step. This speeds up the convergance.
1529 if ( k == iterations ) { // After full solution is achieved, copy low resolution solution into higher res array
1530 for ( Int_t i = ione ; i < rows-1 ; i+=ione ) {
1531 for ( Int_t j = jone ; j < columns-1 ; j+=jone ) {
1534 arrayV(i+ione/2,j) = ( arrayV(i+ione,j) + arrayV(i,j) ) / 2 ;
1535 if ( i == ione ) arrayV(i-ione/2,j) = ( arrayV(0,j) + arrayV(ione,j) ) / 2 ;
1538 arrayV(i,j+jone/2) = ( arrayV(i,j+jone) + arrayV(i,j) ) / 2 ;
1539 if ( j == jone ) arrayV(i,j-jone/2) = ( arrayV(i,0) + arrayV(i,jone) ) / 2 ;
1541 if ( ione > 1 && jone > 1 ) {
1542 arrayV(i+ione/2,j+jone/2) = ( arrayV(i+ione,j+jone) + arrayV(i,j) ) / 2 ;
1543 if ( i == ione ) arrayV(i-ione/2,j-jone/2) = ( arrayV(0,j-jone) + arrayV(ione,j) ) / 2 ;
1544 if ( j == jone ) arrayV(i-ione/2,j-jone/2) = ( arrayV(i-ione,0) + arrayV(i,jone) ) / 2 ;
1545 // Note that this leaves a point at the upper left and lower right corners uninitialized. Not a big deal.
1554 ione = ione / 2 ; if ( ione < 1 ) ione = 1 ;
1555 jone = jone / 2 ; if ( jone < 1 ) jone = 1 ;
1559 //Differentiate V(r) and solve for E(r) using special equations for the first and last row
1560 //Integrate E(r)/E(z) from point of origin to pad plane
1561 AliSysInfo::AddStamp("CalcField", 100,0,0);
1563 for ( Int_t m = 0 ; m < phislices ; m++ ) {
1564 TMatrixD& arrayV = *arrayofArrayV[m] ;
1565 TMatrixD& eroverEz = *arrayofEroverEz[m] ;
1567 for ( Int_t j = columns-1 ; j >= 0 ; j-- ) { // Count backwards to facilitate integration over Z
1569 // Differentiate in R
1570 for ( Int_t i = 1 ; i < rows-1 ; i++ ) arrayE(i,j) = -1 * ( arrayV(i+1,j) - arrayV(i-1,j) ) / (2*gridSizeR) ;
1571 arrayE(0,j) = -1 * ( -0.5*arrayV(2,j) + 2.0*arrayV(1,j) - 1.5*arrayV(0,j) ) / gridSizeR ;
1572 arrayE(rows-1,j) = -1 * ( 1.5*arrayV(rows-1,j) - 2.0*arrayV(rows-2,j) + 0.5*arrayV(rows-3,j) ) / gridSizeR ;
1574 for ( Int_t i = 0 ; i < rows ; i++ ) {
1575 Int_t index = 1 ; // Simpsons rule if N=odd. If N!=odd then add extra point by trapezoidal rule.
1576 eroverEz(i,j) = 0.0 ;
1577 for ( Int_t k = j ; k < columns ; k++ ) {
1579 eroverEz(i,j) += index*(gridSizeZ/3.0)*arrayE(i,k)/(-1*ezField) ;
1580 if ( index != 4 ) index = 4; else index = 2 ;
1582 if ( index == 4 ) eroverEz(i,j) -= (gridSizeZ/3.0)*arrayE(i,columns-1)/ (-1*ezField) ;
1583 if ( index == 2 ) eroverEz(i,j) +=
1584 (gridSizeZ/3.0)*(0.5*arrayE(i,columns-2)-2.5*arrayE(i,columns-1))/(-1*ezField) ;
1585 if ( j == columns-2 ) eroverEz(i,j) =
1586 (gridSizeZ/3.0)*(1.5*arrayE(i,columns-2)+1.5*arrayE(i,columns-1))/(-1*ezField) ;
1587 if ( j == columns-1 ) eroverEz(i,j) = 0.0 ;
1590 // if ( m == 0 ) { TCanvas* c1 = new TCanvas("erOverEz","erOverEz",50,50,840,600) ; c1 -> cd() ;
1591 // eroverEz.Draw("surf") ; } // JT test
1593 AliSysInfo::AddStamp("IntegrateEr", 120,0,0);
1595 //Differentiate V(r) and solve for E(phi)
1596 //Integrate E(phi)/E(z) from point of origin to pad plane
1598 for ( Int_t m = 0 ; m < phislices ; m++ ) {
1600 mplus = m + 1; signplus = 1 ;
1601 mminus = m - 1 ; signminus = 1 ;
1602 if (symmetry==1) { // Reflection symmetry in phi (e.g. symmetry at sector boundaries, or half sectors, etc.)
1603 if ( mplus > phislices-1 ) mplus = phislices - 2 ;
1604 if ( mminus < 0 ) mminus = 1 ;
1606 else if (symmetry==-1) { // Anti-symmetry in phi
1607 if ( mplus > phislices-1 ) { mplus = phislices - 2 ; signplus = -1 ; }
1608 if ( mminus < 0 ) { mminus = 1 ; signminus = -1 ; }
1610 else { // No Symmetries in phi, no boundaries, the calculations is continuous across all phi
1611 if ( mplus > phislices-1 ) mplus = m + 1 - phislices ;
1612 if ( mminus < 0 ) mminus = m - 1 + phislices ;
1614 TMatrixD &arrayVP = *arrayofArrayV[mplus] ;
1615 TMatrixD &arrayVM = *arrayofArrayV[mminus] ;
1616 TMatrixD &ePhioverEz = *arrayofEPhioverEz[m] ;
1617 for ( Int_t j = columns-1 ; j >= 0 ; j-- ) { // Count backwards to facilitate integration over Z
1618 // Differentiate in Phi
1619 for ( Int_t i = 0 ; i < rows ; i++ ) {
1620 Float_t radius = fgkIFCRadius + i*gridSizeR ;
1621 arrayE(i,j) = -1 * (signplus * arrayVP(i,j) - signminus * arrayVM(i,j) ) / (2*radius*gridSizePhi) ;
1624 for ( Int_t i = 0 ; i < rows ; i++ ) {
1625 Int_t index = 1 ; // Simpsons rule if N=odd. If N!=odd then add extra point by trapezoidal rule.
1626 ePhioverEz(i,j) = 0.0 ;
1627 for ( Int_t k = j ; k < columns ; k++ ) {
1629 ePhioverEz(i,j) += index*(gridSizeZ/3.0)*arrayE(i,k)/(-1*ezField) ;
1630 if ( index != 4 ) index = 4; else index = 2 ;
1632 if ( index == 4 ) ePhioverEz(i,j) -= (gridSizeZ/3.0)*arrayE(i,columns-1)/ (-1*ezField) ;
1633 if ( index == 2 ) ePhioverEz(i,j) +=
1634 (gridSizeZ/3.0)*(0.5*arrayE(i,columns-2)-2.5*arrayE(i,columns-1))/(-1*ezField) ;
1635 if ( j == columns-2 ) ePhioverEz(i,j) =
1636 (gridSizeZ/3.0)*(1.5*arrayE(i,columns-2)+1.5*arrayE(i,columns-1))/(-1*ezField) ;
1637 if ( j == columns-1 ) ePhioverEz(i,j) = 0.0 ;
1640 // if ( m == 5 ) { TCanvas* c2 = new TCanvas("arrayE","arrayE",50,50,840,600) ; c2 -> cd() ;
1641 // arrayE.Draw("surf") ; } // JT test
1643 AliSysInfo::AddStamp("IntegrateEphi", 130,0,0);
1646 // Differentiate V(r) and solve for E(z) using special equations for the first and last row
1647 // Integrate (E(z)-Ezstd) from point of origin to pad plane
1649 for ( Int_t m = 0 ; m < phislices ; m++ ) {
1650 TMatrixD& arrayV = *arrayofArrayV[m] ;
1651 TMatrixD& deltaEz = *arrayofDeltaEz[m] ;
1653 // Differentiate V(z) and solve for E(z) using special equations for the first and last columns
1654 for ( Int_t i = 0 ; i < rows ; i++) {
1655 for ( Int_t j = 1 ; j < columns-1 ; j++ ) arrayE(i,j) = -1 * ( arrayV(i,j+1) - arrayV(i,j-1) ) / (2*gridSizeZ) ;
1656 arrayE(i,0) = -1 * ( -0.5*arrayV(i,2) + 2.0*arrayV(i,1) - 1.5*arrayV(i,0) ) / gridSizeZ ;
1657 arrayE(i,columns-1) = -1 * ( 1.5*arrayV(i,columns-1) - 2.0*arrayV(i,columns-2) + 0.5*arrayV(i,columns-3) ) / gridSizeZ ;
1660 for ( Int_t j = columns-1 ; j >= 0 ; j-- ) { // Count backwards to facilitate integration over Z
1662 for ( Int_t i = 0 ; i < rows ; i++ ) {
1663 Int_t index = 1 ; // Simpsons rule if N=odd. If N!=odd then add extra point by trapezoidal rule.
1664 deltaEz(i,j) = 0.0 ;
1665 for ( Int_t k = j ; k < columns ; k++ ) {
1666 deltaEz(i,j) += index*(gridSizeZ/3.0)*arrayE(i,k) ;
1667 if ( index != 4 ) index = 4; else index = 2 ;
1669 if ( index == 4 ) deltaEz(i,j) -= (gridSizeZ/3.0)*arrayE(i,columns-1) ;
1670 if ( index == 2 ) deltaEz(i,j) +=
1671 (gridSizeZ/3.0)*(0.5*arrayE(i,columns-2)-2.5*arrayE(i,columns-1)) ;
1672 if ( j == columns-2 ) deltaEz(i,j) =
1673 (gridSizeZ/3.0)*(1.5*arrayE(i,columns-2)+1.5*arrayE(i,columns-1)) ;
1674 if ( j == columns-1 ) deltaEz(i,j) = 0.0 ;
1678 // if ( m == 0 ) { TCanvas* c1 = new TCanvas("erOverEz","erOverEz",50,50,840,600) ; c1 -> cd() ;
1679 // eroverEz.Draw("surf") ; } // JT test
1681 // calculate z distortion from the integrated Delta Ez residuals
1682 // and include the aquivalence (Volt to cm) of the ROC shift !!
1684 for ( Int_t j = 0 ; j < columns ; j++ ) {
1685 for ( Int_t i = 0 ; i < rows ; i++ ) {
1687 // Scale the Ez distortions with the drift velocity pertubation -> delivers cm
1688 deltaEz(i,j) = deltaEz(i,j)*fgkdvdE;
1690 // ROC Potential in cm aquivalent
1691 Double_t dzROCShift = arrayV(i, columns -1)/ezField;
1692 if ( rocDisplacement ) deltaEz(i,j) = deltaEz(i,j) + dzROCShift; // add the ROC misaligment
1697 } // end loop over phi
1698 AliSysInfo::AddStamp("IntegrateEz", 140,0,0);
1701 for ( Int_t k = 0 ; k < phislices ; k++ )
1703 arrayofSumChargeDensities[k]->Delete() ;
1712 Int_t AliTPCCorrection::IsPowerOfTwo(Int_t i) const {
1714 // Helperfunction: Check if integer is a power of 2
1717 while( i > 0 ) { j += (i&1) ; i = (i>>1) ; }
1718 if ( j == 1 ) return(1) ; // True
1719 return(0) ; // False
1723 AliExternalTrackParam * AliTPCCorrection::FitDistortedTrack(AliExternalTrackParam & trackIn, Double_t refX, Int_t dir, TTreeSRedirector * const pcstream){
1725 // Fit the track parameters - without and with distortion
1726 // 1. Space points in the TPC are simulated along the trajectory
1727 // 2. Space points distorted
1728 // 3. Fits the non distorted and distroted track to the reference plane at refX
1729 // 4. For visualization and debugging purposes the space points and tracks can be stored in the tree - using the TTreeSRedirector functionality
1731 // trackIn - input track parameters
1732 // refX - reference X to fit the track
1733 // dir - direction - out=1 or in=-1
1734 // pcstream - debug streamer to check the results
1736 // see AliExternalTrackParam.h documentation:
1737 // track1.fP[0] - local y (rphi)
1739 // track1.fP[2] - sinus of local inclination angle
1740 // track1.fP[3] - tangent of deep angle
1741 // track1.fP[4] - 1/pt
1743 AliTPCROC * roc = AliTPCROC::Instance();
1744 const Int_t npoints0=roc->GetNRows(0)+roc->GetNRows(36);
1745 const Double_t kRTPC0 =roc->GetPadRowRadii(0,0);
1746 const Double_t kRTPC1 =roc->GetPadRowRadii(36,roc->GetNRows(36)-1);
1747 const Double_t kMaxSnp = 0.85;
1748 const Double_t kSigmaY=0.1;
1749 const Double_t kSigmaZ=0.1;
1750 const Double_t kMaxR=500;
1751 const Double_t kMaxZ=500;
1753 const Double_t kMaxZ0=220;
1754 const Double_t kZcut=3;
1755 const Double_t kMass = TDatabasePDG::Instance()->GetParticle("pi+")->Mass();
1759 AliExternalTrackParam track(trackIn); //
1761 AliTrackPointArray pointArray0(npoints0);
1762 AliTrackPointArray pointArray1(npoints0);
1764 if (!AliTrackerBase::PropagateTrackTo(&track,kRTPC0,kMass,5,kTRUE,kMaxSnp)) return 0;
1766 // simulate the track
1768 Float_t covPoint[6]={0,0,0, static_cast<Float_t>(kSigmaY*kSigmaY),0,static_cast<Float_t>(kSigmaZ*kSigmaZ)}; //covariance at the local frame
1769 for (Double_t radius=kRTPC0; radius<kRTPC1; radius++){
1770 if (!AliTrackerBase::PropagateTrackTo(&track,radius,kMass,5,kTRUE,kMaxSnp)) return 0;
1772 xyz[0]+=gRandom->Gaus(0,0.000005);
1773 xyz[1]+=gRandom->Gaus(0,0.000005);
1774 xyz[2]+=gRandom->Gaus(0,0.000005);
1775 if (TMath::Abs(track.GetZ())>kMaxZ0) continue;
1776 if (TMath::Abs(track.GetX())<kRTPC0) continue;
1777 if (TMath::Abs(track.GetX())>kRTPC1) continue;
1778 AliTrackPoint pIn0; // space point
1780 Int_t sector= (xyz[2]>0)? 0:18;
1781 pointArray0.GetPoint(pIn0,npoints);
1782 pointArray1.GetPoint(pIn1,npoints);
1783 Double_t alpha = TMath::ATan2(xyz[1],xyz[0]);
1784 Float_t distPoint[3]={static_cast<Float_t>(xyz[0]),static_cast<Float_t>(xyz[1]),static_cast<Float_t>(xyz[2])};
1785 DistortPoint(distPoint, sector);
1786 pIn0.SetXYZ(xyz[0], xyz[1],xyz[2]);
1787 pIn1.SetXYZ(distPoint[0], distPoint[1],distPoint[2]);
1789 track.Rotate(alpha);
1790 AliTrackPoint prot0 = pIn0.Rotate(alpha); // rotate to the local frame - non distoted point
1791 AliTrackPoint prot1 = pIn1.Rotate(alpha); // rotate to the local frame - distorted point
1792 prot0.SetXYZ(prot0.GetX(),prot0.GetY(), prot0.GetZ(),covPoint);
1793 prot1.SetXYZ(prot1.GetX(),prot1.GetY(), prot1.GetZ(),covPoint);
1794 pIn0=prot0.Rotate(-alpha); // rotate back to global frame
1795 pIn1=prot1.Rotate(-alpha); // rotate back to global frame
1796 pointArray0.AddPoint(npoints, &pIn0);
1797 pointArray1.AddPoint(npoints, &pIn1);
1799 if (npoints>=npoints0) break;
1801 if (npoints<npoints0/4.) return 0;
1805 AliExternalTrackParam *track0=0;
1806 AliExternalTrackParam *track1=0;
1807 AliTrackPoint point1,point2,point3;
1808 if (dir==1) { //make seed inner
1809 pointArray0.GetPoint(point1,1);
1810 pointArray0.GetPoint(point2,11);
1811 pointArray0.GetPoint(point3,21);
1813 if (dir==-1){ //make seed outer
1814 pointArray0.GetPoint(point1,npoints-21);
1815 pointArray0.GetPoint(point2,npoints-11);
1816 pointArray0.GetPoint(point3,npoints-1);
1818 if ((TMath::Abs(point1.GetX()-point3.GetX())+TMath::Abs(point1.GetY()-point3.GetY()))<10){
1819 printf("fit points not properly initialized\n");
1822 track0 = AliTrackerBase::MakeSeed(point1, point2, point3);
1823 track1 = AliTrackerBase::MakeSeed(point1, point2, point3);
1824 track0->ResetCovariance(10);
1825 track1->ResetCovariance(10);
1826 if (TMath::Abs(AliTrackerBase::GetBz())<0.01){
1827 ((Double_t*)track0->GetParameter())[4]= trackIn.GetParameter()[4];
1828 ((Double_t*)track1->GetParameter())[4]= trackIn.GetParameter()[4];
1830 for (Int_t jpoint=0; jpoint<npoints; jpoint++){
1831 Int_t ipoint= (dir>0) ? jpoint: npoints-1-jpoint;
1835 pointArray0.GetPoint(pIn0,ipoint);
1836 pointArray1.GetPoint(pIn1,ipoint);
1837 AliTrackPoint prot0 = pIn0.Rotate(track0->GetAlpha()); // rotate to the local frame - non distoted point
1838 AliTrackPoint prot1 = pIn1.Rotate(track1->GetAlpha()); // rotate to the local frame - distorted point
1839 if (TMath::Abs(prot0.GetX())<kRTPC0) continue;
1840 if (TMath::Abs(prot0.GetX())>kRTPC1) continue;
1842 if (!AliTrackerBase::PropagateTrackTo(track0,prot0.GetX(),kMass,5,kFALSE,kMaxSnp)) break;
1843 if (!AliTrackerBase::PropagateTrackTo(track1,prot0.GetX(),kMass,5,kFALSE,kMaxSnp)) break;
1844 if (TMath::Abs(track0->GetZ())>kMaxZ) break;
1845 if (TMath::Abs(track0->GetX())>kMaxR) break;
1846 if (TMath::Abs(track1->GetZ())>kMaxZ) break;
1847 if (TMath::Abs(track1->GetX())>kMaxR) break;
1848 if (dir>0 && track1->GetX()>refX) continue;
1849 if (dir<0 && track1->GetX()<refX) continue;
1850 if (TMath::Abs(track1->GetZ())<kZcut)continue;
1851 track.GetXYZ(xyz); // distorted track also propagated to the same reference radius
1853 Double_t pointPos[2]={0,0};
1854 Double_t pointCov[3]={0,0,0};
1855 pointPos[0]=prot0.GetY();//local y
1856 pointPos[1]=prot0.GetZ();//local z
1857 pointCov[0]=prot0.GetCov()[3];//simay^2
1858 pointCov[1]=prot0.GetCov()[4];//sigmayz
1859 pointCov[2]=prot0.GetCov()[5];//sigmaz^2
1860 if (!track0->Update(pointPos,pointCov)) break;
1862 Double_t deltaX=prot1.GetX()-prot0.GetX(); // delta X
1863 Double_t deltaYX=deltaX*TMath::Tan(TMath::ASin(track1->GetSnp())); // deltaY due delta X
1864 Double_t deltaZX=deltaX*track1->GetTgl(); // deltaZ due delta X
1866 pointPos[0]=prot1.GetY()-deltaYX;//local y is sign correct? should be minus
1867 pointPos[1]=prot1.GetZ()-deltaZX;//local z is sign correct? should be minus
1868 pointCov[0]=prot1.GetCov()[3];//simay^2
1869 pointCov[1]=prot1.GetCov()[4];//sigmayz
1870 pointCov[2]=prot1.GetCov()[5];//sigmaz^2
1871 if (!track1->Update(pointPos,pointCov)) break;
1875 if (npoints2<npoints/4.) return 0;
1876 AliTrackerBase::PropagateTrackTo(track0,refX,kMass,5.,kTRUE,kMaxSnp);
1877 AliTrackerBase::PropagateTrackTo(track0,refX,kMass,1.,kTRUE,kMaxSnp);
1878 track1->Rotate(track0->GetAlpha());
1879 AliTrackerBase::PropagateTrackTo(track1,track0->GetX(),kMass,5.,kFALSE,kMaxSnp);
1881 if (pcstream) (*pcstream)<<Form("fitDistort%s",GetName())<<
1882 "point0.="<<&pointArray0<< // points
1883 "point1.="<<&pointArray1<< // distorted points
1884 "trackIn.="<<&track<< // original track
1885 "track0.="<<track0<< // fitted track
1886 "track1.="<<track1<< // fitted distorted track
1888 new(&trackIn) AliExternalTrackParam(*track0);
1897 TTree* AliTPCCorrection::CreateDistortionTree(Double_t step){
1899 // create the distortion tree on a mesh with granularity given by step
1900 // return the tree with distortions at given position
1901 // Map is created on the mesh with given step size
1903 TTreeSRedirector *pcstream = new TTreeSRedirector(Form("correction%s.root",GetName()));
1905 for (Double_t x= -250; x<250; x+=step){
1906 for (Double_t y= -250; y<250; y+=step){
1907 Double_t r = TMath::Sqrt(x*x+y*y);
1909 if (r>250) continue;
1910 for (Double_t z= -250; z<250; z+=step){
1911 Int_t roc=(z>0)?0:18;
1915 Double_t phi = TMath::ATan2(y,x);
1916 DistortPoint(xyz,roc);
1917 Double_t r1 = TMath::Sqrt(xyz[0]*xyz[0]+xyz[1]*xyz[1]);
1918 Double_t phi1 = TMath::ATan2(xyz[1],xyz[0]);
1919 if ((phi1-phi)>TMath::Pi()) phi1-=TMath::Pi();
1920 if ((phi1-phi)<-TMath::Pi()) phi1+=TMath::Pi();
1921 Double_t dx = xyz[0]-x;
1922 Double_t dy = xyz[1]-y;
1923 Double_t dz = xyz[2]-z;
1925 Double_t drphi=(phi1-phi)*r;
1926 (*pcstream)<<"distortion"<<
1927 "x="<<x<< // original position
1932 "x1="<<xyz[0]<< // distorted position
1938 "dx="<<dx<< // delta position
1948 TFile f(Form("correction%s.root",GetName()));
1949 TTree * tree = (TTree*)f.Get("distortion");
1950 TTree * tree2= tree->CopyTree("1");
1951 tree2->SetName(Form("dist%s",GetName()));
1952 tree2->SetDirectory(0);
1960 void AliTPCCorrection::MakeTrackDistortionTree(TTree *tinput, Int_t dtype, Int_t ptype, const TObjArray * corrArray, Int_t step, Int_t offset, Bool_t debug ){
1963 // For each partial correction (specified in array) and given track topology (phi, theta, snp, refX)
1964 // calculates partial distortions
1965 // Partial distortion is stored in the resulting tree
1966 // Output is storred in the file distortion_<dettype>_<partype>.root
1967 // Partial distortion is stored with the name given by correction name
1970 // Parameters of function:
1971 // input - input tree
1972 // dtype - distortion type 0 - ITSTPC, 1 -TPCTRD, 2 - TPCvertex , 3 - TPC-TOF, 4 - TPCTPC track crossing
1973 // ppype - parameter type
1974 // corrArray - array with partial corrections
1975 // step - skipe entries - if 1 all entries processed - it is slow
1976 // debug 0 if debug on also space points dumped - it is slow
1978 const Double_t kMaxSnp = 0.85;
1979 const Double_t kcutSnp=0.25;
1980 const Double_t kcutTheta=1.;
1981 const Double_t kRadiusTPC=85;
1982 // AliTPCROC *tpcRoc =AliTPCROC::Instance();
1984 const Double_t kMass = TDatabasePDG::Instance()->GetParticle("pi+")->Mass();
1985 // const Double_t kB2C=-0.299792458e-3;
1986 const Int_t kMinEntries=20;
1987 Double_t phi,theta, snp, mean,rms, entries,sector,dsec;
1990 tinput->SetBranchAddress("run",&run);
1991 tinput->SetBranchAddress("theta",&theta);
1992 tinput->SetBranchAddress("phi", &phi);
1993 tinput->SetBranchAddress("snp",&snp);
1994 tinput->SetBranchAddress("mean",&mean);
1995 tinput->SetBranchAddress("rms",&rms);
1996 tinput->SetBranchAddress("entries",&entries);
1997 tinput->SetBranchAddress("sector",§or);
1998 tinput->SetBranchAddress("dsec",&dsec);
1999 tinput->SetBranchAddress("refX",&refX);
2000 TTreeSRedirector *pcstream = new TTreeSRedirector(Form("distortion%d_%d_%d.root",dtype,ptype,offset));
2002 Int_t nentries=tinput->GetEntries();
2003 Int_t ncorr=corrArray->GetEntries();
2004 Double_t corrections[100]={0}; //
2006 Double_t cov[15]={0,0,0,0,0,0,0,0,0,0,0,0,0,0};
2008 if (dtype==5 || dtype==6) dtype=4;
2009 if (dtype==0) { dir=-1;}
2010 if (dtype==1) { dir=1;}
2011 if (dtype==2) { dir=-1;}
2012 if (dtype==3) { dir=1;}
2013 if (dtype==4) { dir=-1;}
2015 for (Int_t ientry=offset; ientry<nentries; ientry+=step){
2016 tinput->GetEntry(ientry);
2017 if (TMath::Abs(snp)>kMaxSnp) continue;
2020 if (dtype==2) tPar[1]=theta*kRadiusTPC;
2023 tPar[4]=(gRandom->Rndm()-0.5)*0.02; // should be calculated - non equal to 0
2025 // tracks crossing CE
2026 tPar[1]=0; // track at the CE
2027 //if (TMath::Abs(theta) <0.05) continue; // deep cross
2030 if (TMath::Abs(snp) >kcutSnp) continue;
2031 if (TMath::Abs(theta) >kcutTheta) continue;
2032 printf("%f\t%f\t%f\t%f\t%f\t%f\n",entries, sector,theta,snp, mean,rms);
2033 Double_t bz=AliTrackerBase::GetBz();
2034 if (dtype !=4) { //exclude TPC - for TPC mainly non primary tracks
2035 if (dtype!=2 && TMath::Abs(bz)>0.1 ) tPar[4]=snp/(refX*bz*kB2C*2);
2037 if (dtype==2 && TMath::Abs(bz)>0.1 ) {
2038 tPar[4]=snp/(kRadiusTPC*bz*kB2C*2);//
2039 // snp at the TPC inner radius in case the vertex match used
2043 tPar[4]+=(gRandom->Rndm()-0.5)*0.02;
2044 AliExternalTrackParam track(refX,phi,tPar,cov);
2048 Double_t pt=1./tPar[4];
2049 Double_t dRrec=0; // dummy value - needed for points - e.g for laser
2050 //if (ptype==4 &&bz<0) mean*=-1; // interpret as curvature -- COMMENTED out - in lookup signed 1/pt used
2051 Double_t refXD=refX;
2052 (*pcstream)<<"fit"<<
2053 "run="<<run<< // run number
2054 "bz="<<bz<< // magnetic filed used
2055 "dtype="<<dtype<< // detector match type
2056 "ptype="<<ptype<< // parameter type
2057 "theta="<<theta<< // theta
2058 "phi="<<phi<< // phi
2059 "snp="<<snp<< // snp
2060 "mean="<<mean<< // mean dist value
2061 "rms="<<rms<< // rms
2064 "refX="<<refXD<< // referece X as double
2065 "gx="<<xyz[0]<< // global position at reference
2066 "gy="<<xyz[1]<< // global position at reference
2067 "gz="<<xyz[2]<< // global position at reference
2068 "dRrec="<<dRrec<< // delta Radius in reconstruction
2070 "id="<<id<< // track id
2071 "entries="<<entries;// number of entries in bin
2074 if (entries<kMinEntries) isOK=kFALSE;
2076 if (dtype!=4) for (Int_t icorr=0; icorr<ncorr; icorr++) {
2077 AliTPCCorrection *corr = (AliTPCCorrection*)corrArray->At(icorr);
2078 corrections[icorr]=0;
2079 if (entries>kMinEntries){
2080 AliExternalTrackParam trackIn(refX,phi,tPar,cov);
2081 AliExternalTrackParam *trackOut = 0;
2082 if (debug) trackOut=corr->FitDistortedTrack(trackIn, refX, dir,pcstream);
2083 if (!debug) trackOut=corr->FitDistortedTrack(trackIn, refX, dir,0);
2084 if (dtype==0) {dir= -1;}
2085 if (dtype==1) {dir= 1;}
2086 if (dtype==2) {dir= -1;}
2087 if (dtype==3) {dir= 1;}
2090 if (!AliTrackerBase::PropagateTrackTo(&trackIn,refX,kMass,5,kTRUE,kMaxSnp)) isOK=kFALSE;
2091 if (!trackOut->Rotate(trackIn.GetAlpha())) isOK=kFALSE;
2092 if (!AliTrackerBase::PropagateTrackTo(trackOut,trackIn.GetX(),kMass,5,kFALSE,kMaxSnp)) isOK=kFALSE;
2093 // trackOut->PropagateTo(trackIn.GetX(),AliTrackerBase::GetBz());
2095 corrections[icorr]= trackOut->GetParameter()[ptype]-trackIn.GetParameter()[ptype];
2098 corrections[icorr]=0;
2101 //if (ptype==4 &&bz<0) corrections[icorr]*=-1; // interpret as curvature - commented out
2103 (*pcstream)<<"fit"<<
2104 Form("%s=",corr->GetName())<<corrections[icorr]; // dump correction value
2107 if (dtype==4) for (Int_t icorr=0; icorr<ncorr; icorr++) {
2109 // special case of the TPC tracks crossing the CE
2111 AliTPCCorrection *corr = (AliTPCCorrection*)corrArray->At(icorr);
2112 corrections[icorr]=0;
2113 if (entries>kMinEntries){
2114 AliExternalTrackParam trackIn0(refX,phi,tPar,cov); //Outer - direction to vertex
2115 AliExternalTrackParam trackIn1(refX,phi,tPar,cov); //Inner - direction magnet
2116 AliExternalTrackParam *trackOut0 = 0;
2117 AliExternalTrackParam *trackOut1 = 0;
2119 if (debug) trackOut0=corr->FitDistortedTrack(trackIn0, refX, dir,pcstream);
2120 if (!debug) trackOut0=corr->FitDistortedTrack(trackIn0, refX, dir,0);
2121 if (debug) trackOut1=corr->FitDistortedTrack(trackIn1, refX, -dir,pcstream);
2122 if (!debug) trackOut1=corr->FitDistortedTrack(trackIn1, refX, -dir,0);
2124 if (trackOut0 && trackOut1){
2125 if (!AliTrackerBase::PropagateTrackTo(&trackIn0,refX,kMass,5,kTRUE,kMaxSnp)) isOK=kFALSE;
2126 if (!AliTrackerBase::PropagateTrackTo(&trackIn0,refX,kMass,1,kFALSE,kMaxSnp)) isOK=kFALSE;
2127 if (!trackOut0->Rotate(trackIn0.GetAlpha())) isOK=kFALSE;
2128 if (!AliTrackerBase::PropagateTrackTo(trackOut0,trackIn0.GetX(),kMass,5,kFALSE,kMaxSnp)) isOK=kFALSE;
2130 if (!AliTrackerBase::PropagateTrackTo(&trackIn1,refX,kMass,5,kTRUE,kMaxSnp)) isOK=kFALSE;
2131 if (!trackIn1.Rotate(trackIn0.GetAlpha())) isOK=kFALSE;
2132 if (!AliTrackerBase::PropagateTrackTo(&trackIn1,trackIn0.GetX(),kMass,1,kFALSE,kMaxSnp)) isOK=kFALSE;
2133 if (!trackOut1->Rotate(trackIn1.GetAlpha())) isOK=kFALSE;
2134 if (!AliTrackerBase::PropagateTrackTo(trackOut1,trackIn1.GetX(),kMass,5,kFALSE,kMaxSnp)) isOK=kFALSE;
2136 corrections[icorr] = (trackOut0->GetParameter()[ptype]-trackIn0.GetParameter()[ptype]);
2137 corrections[icorr]-= (trackOut1->GetParameter()[ptype]-trackIn1.GetParameter()[ptype]);
2139 if ((TMath::Abs(trackOut0->GetX()-trackOut1->GetX())>0.1)||
2140 (TMath::Abs(trackOut0->GetX()-trackIn1.GetX())>0.1)||
2141 (TMath::Abs(trackOut0->GetAlpha()-trackOut1->GetAlpha())>0.00001)||
2142 (TMath::Abs(trackOut0->GetAlpha()-trackIn1.GetAlpha())>0.00001)||
2143 (TMath::Abs(trackIn0.GetTgl()-trackIn1.GetTgl())>0.0001)||
2144 (TMath::Abs(trackIn0.GetSnp()-trackIn1.GetSnp())>0.0001)
2151 corrections[icorr]=0;
2155 //if (ptype==4 &&bz<0) corrections[icorr]*=-1; // interpret as curvature - commented out no in lookup
2157 (*pcstream)<<"fit"<<
2158 Form("%s=",corr->GetName())<<corrections[icorr]; // dump correction value
2161 (*pcstream)<<"fit"<<"isOK="<<isOK<<"\n";
2170 void AliTPCCorrection::MakeSectorDistortionTree(TTree *tinput, Int_t dtype, Int_t ptype, const TObjArray * corrArray, Int_t step, Int_t offset, Bool_t debug ){
2173 // For each partial correction (specified in array) and given track topology (phi, theta, snp, refX)
2174 // calculates partial distortions
2175 // Partial distortion is stored in the resulting tree
2176 // Output is storred in the file distortion_<dettype>_<partype>.root
2177 // Partial distortion is stored with the name given by correction name
2180 // Parameters of function:
2181 // input - input tree
2182 // dtype - distortion type 10 - IROC-OROC
2183 // ppype - parameter type
2184 // corrArray - array with partial corrections
2185 // step - skipe entries - if 1 all entries processed - it is slow
2186 // debug 0 if debug on also space points dumped - it is slow
2188 const Double_t kMaxSnp = 0.8;
2189 const Int_t kMinEntries=200;
2190 // AliTPCROC *tpcRoc =AliTPCROC::Instance();
2192 const Double_t kMass = TDatabasePDG::Instance()->GetParticle("pi+")->Mass();
2193 // const Double_t kB2C=-0.299792458e-3;
2194 Double_t phi,theta, snp, mean,rms, entries,sector,dsec,globalZ;
2199 tinput->SetBranchAddress("run",&run);
2200 tinput->SetBranchAddress("theta",&theta);
2201 tinput->SetBranchAddress("phi", &phi);
2202 tinput->SetBranchAddress("snp",&snp);
2203 tinput->SetBranchAddress("mean",&mean);
2204 tinput->SetBranchAddress("rms",&rms);
2205 tinput->SetBranchAddress("entries",&entries);
2206 tinput->SetBranchAddress("sector",§or);
2207 tinput->SetBranchAddress("dsec",&dsec);
2208 tinput->SetBranchAddress("refX",&refXD);
2209 tinput->SetBranchAddress("z",&globalZ);
2210 tinput->SetBranchAddress("isec0",&isec0);
2211 tinput->SetBranchAddress("isec1",&isec1);
2212 TTreeSRedirector *pcstream = new TTreeSRedirector(Form("distortionSector%d_%d_%d.root",dtype,ptype,offset));
2214 Int_t nentries=tinput->GetEntries();
2215 Int_t ncorr=corrArray->GetEntries();
2216 Double_t corrections[100]={0}; //
2218 Double_t cov[15]={0,0,0,0,0,0,0,0,0,0,0,0,0,0};
2221 for (Int_t ientry=offset; ientry<nentries; ientry+=step){
2222 tinput->GetEntry(ientry);
2225 if (TMath::Abs(TMath::Abs(isec0%18)-TMath::Abs(isec1%18))==0) id=1; // IROC-OROC - opposite side
2226 if (TMath::Abs(TMath::Abs(isec0%36)-TMath::Abs(isec1%36))==0) id=2; // IROC-OROC - same side
2227 if (dtype==10 && id==-1) continue;
2234 tPar[4]=(gRandom->Rndm()-0.1)*0.2; //
2235 Double_t pt=1./tPar[4];
2237 printf("%f\t%f\t%f\t%f\t%f\t%f\n",entries, sector,theta,snp, mean,rms);
2238 Double_t bz=AliTrackerBase::GetBz();
2239 AliExternalTrackParam track(refX,phi,tPar,cov);
2240 Double_t xyz[3],xyzIn[3],xyzOut[3];
2242 track.GetXYZAt(85,bz,xyzIn);
2243 track.GetXYZAt(245,bz,xyzOut);
2244 Double_t phiIn = TMath::ATan2(xyzIn[1],xyzIn[0]);
2245 Double_t phiOut = TMath::ATan2(xyzOut[1],xyzOut[0]);
2246 Double_t phiRef = TMath::ATan2(xyz[1],xyz[0]);
2247 Int_t sectorRef = TMath::Nint(9.*phiRef/TMath::Pi()-0.5);
2248 Int_t sectorIn = TMath::Nint(9.*phiIn/TMath::Pi()-0.5);
2249 Int_t sectorOut = TMath::Nint(9.*phiOut/TMath::Pi()-0.5);
2252 if (sectorIn!=sectorOut) isOK=kFALSE; // requironment - cluster in the same sector
2253 if (sectorIn!=sectorRef) isOK=kFALSE; // requironment - cluster in the same sector
2254 if (entries<kMinEntries/(1+TMath::Abs(globalZ/100.))) isOK=kFALSE; // requironment - minimal amount of tracks in bin
2256 if (TMath::Abs(theta)>1) isOK=kFALSE;
2258 Double_t dRrec=0; // dummy value - needed for points - e.g for laser
2260 (*pcstream)<<"fit"<<
2262 "bz="<<bz<< // magnetic filed used
2263 "dtype="<<dtype<< // detector match type
2264 "ptype="<<ptype<< // parameter type
2265 "theta="<<theta<< // theta
2266 "phi="<<phi<< // phi
2267 "snp="<<snp<< // snp
2268 "mean="<<mean<< // mean dist value
2269 "rms="<<rms<< // rms
2272 "refX="<<refXD<< // referece X
2273 "gx="<<xyz[0]<< // global position at reference
2274 "gy="<<xyz[1]<< // global position at reference
2275 "gz="<<xyz[2]<< // global position at reference
2276 "dRrec="<<dRrec<< // delta Radius in reconstruction
2278 "id="<<id<< // track id
2279 "entries="<<entries;// number of entries in bin
2281 AliExternalTrackParam *trackOut0 = 0;
2282 AliExternalTrackParam *trackOut1 = 0;
2283 AliExternalTrackParam *ptrackIn0 = 0;
2284 AliExternalTrackParam *ptrackIn1 = 0;
2286 for (Int_t icorr=0; icorr<ncorr; icorr++) {
2288 // special case of the TPC tracks crossing the CE
2290 AliTPCCorrection *corr = (AliTPCCorrection*)corrArray->At(icorr);
2291 corrections[icorr]=0;
2292 if (entries>kMinEntries &&isOK){
2293 AliExternalTrackParam trackIn0(refX,phi,tPar,cov);
2294 AliExternalTrackParam trackIn1(refX,phi,tPar,cov);
2295 ptrackIn1=&trackIn0;
2296 ptrackIn0=&trackIn1;
2298 if (debug) trackOut0=corr->FitDistortedTrack(trackIn0, refX, dir,pcstream);
2299 if (!debug) trackOut0=corr->FitDistortedTrack(trackIn0, refX, dir,0);
2300 if (debug) trackOut1=corr->FitDistortedTrack(trackIn1, refX, -dir,pcstream);
2301 if (!debug) trackOut1=corr->FitDistortedTrack(trackIn1, refX, -dir,0);
2303 if (trackOut0 && trackOut1){
2305 if (!AliTrackerBase::PropagateTrackTo(&trackIn0,refX,kMass,1,kTRUE,kMaxSnp)) isOK=kFALSE;
2306 if (!AliTrackerBase::PropagateTrackTo(&trackIn0,refX,kMass,1,kFALSE,kMaxSnp)) isOK=kFALSE;
2307 // rotate all tracks to the same frame
2308 if (!trackOut0->Rotate(trackIn0.GetAlpha())) isOK=kFALSE;
2309 if (!trackIn1.Rotate(trackIn0.GetAlpha())) isOK=kFALSE;
2310 if (!trackOut1->Rotate(trackIn0.GetAlpha())) isOK=kFALSE;
2312 if (!AliTrackerBase::PropagateTrackTo(trackOut0,refX,kMass,1,kFALSE,kMaxSnp)) isOK=kFALSE;
2313 if (!AliTrackerBase::PropagateTrackTo(&trackIn1,refX,kMass,1,kFALSE,kMaxSnp)) isOK=kFALSE;
2314 if (!AliTrackerBase::PropagateTrackTo(trackOut1,refX,kMass,1,kFALSE,kMaxSnp)) isOK=kFALSE;
2316 corrections[icorr] = (trackOut0->GetParameter()[ptype]-trackIn0.GetParameter()[ptype]);
2317 corrections[icorr]-= (trackOut1->GetParameter()[ptype]-trackIn1.GetParameter()[ptype]);
2318 (*pcstream)<<"fitDebug"<< // just to debug the correction
2320 "pIn0.="<<ptrackIn0<<
2321 "pIn1.="<<ptrackIn1<<
2322 "pOut0.="<<trackOut0<<
2323 "pOut1.="<<trackOut1<<
2329 corrections[icorr]=0;
2333 (*pcstream)<<"fit"<<
2334 Form("%s=",corr->GetName())<<corrections[icorr]; // dump correction value
2337 (*pcstream)<<"fit"<<"isOK="<<isOK<<"\n";
2344 void AliTPCCorrection::MakeLaserDistortionTreeOld(TTree* tree, TObjArray *corrArray, Int_t itype){
2346 // Make a laser fit tree for global minimization
2348 const Double_t cutErrY=0.1;
2349 const Double_t cutErrZ=0.1;
2350 const Double_t kEpsilon=0.00000001;
2351 const Double_t kMaxDist=1.; // max distance - space correction
2352 const Double_t kMaxRMS=0.05; // max distance -between point and local mean
2357 AliTPCLaserTrack *ltr=0;
2358 AliTPCLaserTrack::LoadTracks();
2359 tree->SetBranchAddress("dY.",&vecdY);
2360 tree->SetBranchAddress("dZ.",&vecdZ);
2361 tree->SetBranchAddress("eY.",&veceY);
2362 tree->SetBranchAddress("eZ.",&veceZ);
2363 tree->SetBranchAddress("LTr.",<r);
2364 Int_t entries= tree->GetEntries();
2365 TTreeSRedirector *pcstream= new TTreeSRedirector("distortionLaser_0.root");
2366 Double_t bz=AliTrackerBase::GetBz();
2369 for (Int_t ientry=0; ientry<entries; ientry++){
2370 tree->GetEntry(ientry);
2371 if (!ltr->GetVecGX()){
2372 ltr->UpdatePoints();
2374 TVectorD * delta= (itype==0)? vecdY:vecdZ;
2375 TVectorD * err= (itype==0)? veceY:veceZ;
2376 TLinearFitter fitter(2,"pol1");
2377 for (Int_t iter=0; iter<2; iter++){
2378 Double_t kfit0=0, kfit1=0;
2379 Int_t npoints=fitter.GetNpoints();
2382 kfit0=fitter.GetParameter(0);
2383 kfit1=fitter.GetParameter(1);
2385 for (Int_t irow=0; irow<159; irow++){
2388 Int_t nentries = 1000;
2389 if (veceY->GetMatrixArray()[irow]>cutErrY||veceZ->GetMatrixArray()[irow]>cutErrZ) nentries=0;
2390 if (veceY->GetMatrixArray()[irow]<kEpsilon||veceZ->GetMatrixArray()[irow]<kEpsilon) nentries=0;
2393 Int_t first3=TMath::Max(irow-3,0);
2394 Int_t last3 =TMath::Min(irow+3,159);
2396 if ((*ltr->GetVecSec())[irow]>=0 && err) {
2397 for (Int_t jrow=first3; jrow<=last3; jrow++){
2398 if ((*ltr->GetVecSec())[irow]!= (*ltr->GetVecSec())[jrow]) continue;
2399 if ((*err)[jrow]<kEpsilon) continue;
2400 array[counter]=(*delta)[jrow];
2407 rms3 = TMath::RMS(counter,array);
2408 mean3 = TMath::Mean(counter,array);
2412 Double_t phi =(*ltr->GetVecPhi())[irow];
2413 Double_t theta =ltr->GetTgl();
2414 Double_t mean=delta->GetMatrixArray()[irow];
2415 Double_t gx=0,gy=0,gz=0;
2416 Double_t snp = (*ltr->GetVecP2())[irow];
2418 // Double_t rms = err->GetMatrixArray()[irow];
2420 gx = (*ltr->GetVecGX())[irow];
2421 gy = (*ltr->GetVecGY())[irow];
2422 gz = (*ltr->GetVecGZ())[irow];
2424 // get delta R used in reconstruction
2425 AliTPCcalibDB* calib=AliTPCcalibDB::Instance();
2426 AliTPCCorrection * correction = calib->GetTPCComposedCorrection(AliTrackerBase::GetBz());
2427 // const AliTPCRecoParam * recoParam = calib->GetTransform()->GetCurrentRecoParam();
2428 //Double_t xyz0[3]={gx,gy,gz};
2429 Double_t oldR=TMath::Sqrt(gx*gx+gy*gy);
2430 Double_t fphi = TMath::ATan2(gy,gx);
2431 Double_t fsector = 9.*fphi/TMath::Pi();
2432 if (fsector<0) fsector+=18;
2433 Double_t dsec = fsector-Int_t(fsector)-0.5;
2435 Int_t id= ltr->GetId();
2439 Float_t xyz1[3]={static_cast<Float_t>(gx),static_cast<Float_t>(gy),static_cast<Float_t>(gz)};
2440 Int_t sector=(gz>0)?0:18;
2441 correction->CorrectPoint(xyz1, sector);
2442 refX=TMath::Sqrt(xyz1[0]*xyz1[0]+xyz1[1]*xyz1[1]);
2445 if (TMath::Abs(rms3)>kMaxRMS) isOK=kFALSE;
2446 if (TMath::Abs(mean-mean3)>kMaxRMS) isOK=kFALSE;
2447 if (counter<4) isOK=kFALSE;
2448 if (npoints<90) isOK=kFALSE;
2450 fitter.AddPoint(&refX,mean);
2452 Double_t deltaF=kfit0+kfit1*refX;
2454 (*pcstream)<<"fitFull"<< // dumpe also intermediate results
2455 "bz="<<bz<< // magnetic filed used
2456 "dtype="<<dtype<< // detector match type
2457 "ptype="<<itype<< // parameter type
2458 "theta="<<theta<< // theta
2459 "phi="<<phi<< // phi
2460 "snp="<<snp<< // snp
2461 "mean="<<mean3<< // mean dist value
2462 "rms="<<rms3<< // rms
2464 "npoints="<<npoints<< //number of points
2465 "mean3="<<mean3<< // mean dist value
2466 "rms3="<<rms3<< // rms
2467 "counter="<<counter<<
2468 "sector="<<fsector<<
2471 "refX="<<refX<< // reference radius
2472 "gx="<<gx<< // global position
2473 "gy="<<gy<< // global position
2474 "gz="<<gz<< // global position
2475 "dRrec="<<dRrec<< // delta Radius in reconstruction
2476 "id="<<id<< //bundle
2477 "entries="<<nentries<<// number of entries in bin
2480 if (iter==1) (*pcstream)<<"fit"<< // dump valus for fit
2481 "bz="<<bz<< // magnetic filed used
2482 "dtype="<<dtype<< // detector match type
2483 "ptype="<<itype<< // parameter type
2484 "theta="<<theta<< // theta
2485 "phi="<<phi<< // phi
2486 "snp="<<snp<< // snp
2487 "mean="<<mean3<< // mean dist value
2488 "rms="<<rms3<< // rms
2489 "sector="<<fsector<<
2492 "refX="<<refX<< // reference radius
2493 "gx="<<gx<< // global position
2494 "gy="<<gy<< // global position
2495 "gz="<<gz<< // global position
2496 "dRrec="<<dRrec<< // delta Radius in reconstruction
2498 "id="<<id<< //bundle
2499 "entries="<<nentries;// number of entries in bin
2502 Double_t ky = TMath::Tan(TMath::ASin(snp));
2503 Int_t ncorr = corrArray->GetEntries();
2504 Double_t r0 = TMath::Sqrt(gx*gx+gy*gy);
2505 Double_t phi0 = TMath::ATan2(gy,gx);
2506 Double_t distortions[1000]={0};
2507 Double_t distortionsR[1000]={0};
2509 for (Int_t icorr=0; icorr<ncorr; icorr++) {
2510 AliTPCCorrection *corr = (AliTPCCorrection*)corrArray->At(icorr);
2511 Float_t distPoint[3]={static_cast<Float_t>(gx),static_cast<Float_t>(gy),static_cast<Float_t>(gz)};
2512 Int_t sector= (gz>0)? 0:18;
2514 corr->DistortPoint(distPoint, sector);
2516 // Double_t value=distPoint[2]-gz;
2517 if (itype==0 && r0>1){
2518 Double_t r1 = TMath::Sqrt(distPoint[0]*distPoint[0]+distPoint[1]*distPoint[1]);
2519 Double_t phi1 = TMath::ATan2(distPoint[1],distPoint[0]);
2520 Double_t drphi= r0*(phi1-phi0);
2521 Double_t dr = r1-r0;
2522 distortions[icorr] = drphi-ky*dr;
2523 distortionsR[icorr] = dr;
2525 if (TMath::Abs(distortions[icorr])>kMaxDist) {isOKF=icorr+1; isOK=kFALSE; }
2526 if (TMath::Abs(distortionsR[icorr])>kMaxDist) {isOKF=icorr+1; isOK=kFALSE;}
2527 (*pcstream)<<"fit"<<
2528 Form("%s=",corr->GetName())<<distortions[icorr]; // dump correction value
2530 (*pcstream)<<"fit"<<"isOK="<<isOK<<"\n";
2540 void AliTPCCorrection::MakeDistortionMap(THnSparse * his0, TTreeSRedirector * const pcstream, const char* hname, Int_t run, Float_t refX, Int_t type, Int_t integ){
2542 // make a distortion map out ou fthe residual histogram
2543 // Results are written to the debug streamer - pcstream
2545 // his0 - input (4D) residual histogram
2546 // pcstream - file to write the tree
2548 // refX - track matching reference X
2549 // type - 0- y 1-z,2 -snp, 3-theta, 4=1/pt
2551 // OBJ: TAxis #Delta #Delta
2552 // OBJ: TAxis tanTheta tan(#Theta)
2553 // OBJ: TAxis phi #phi
2554 // OBJ: TAxis snp snp
2556 // marian.ivanov@cern.ch
2557 const Int_t kMinEntries=10;
2558 Double_t bz=AliTrackerBase::GetBz();
2559 Int_t idim[4]={0,1,2,3};
2563 Int_t nbins3=his0->GetAxis(3)->GetNbins();
2564 Int_t first3=his0->GetAxis(3)->GetFirst();
2565 Int_t last3 =his0->GetAxis(3)->GetLast();
2567 for (Int_t ibin3=first3; ibin3<last3; ibin3+=1){ // axis 3 - local angle
2568 his0->GetAxis(3)->SetRange(TMath::Max(ibin3-integ,1),TMath::Min(ibin3+integ,nbins3));
2569 Double_t x3= his0->GetAxis(3)->GetBinCenter(ibin3);
2570 THnSparse * his3= his0->Projection(3,idim); //projected histogram according selection 3
2572 Int_t nbins2 = his3->GetAxis(2)->GetNbins();
2573 Int_t first2 = his3->GetAxis(2)->GetFirst();
2574 Int_t last2 = his3->GetAxis(2)->GetLast();
2576 for (Int_t ibin2=first2; ibin2<last2; ibin2+=1){ // axis 2 - phi
2577 his3->GetAxis(2)->SetRange(TMath::Max(ibin2-integ,1),TMath::Min(ibin2+integ,nbins2));
2578 Double_t x2= his3->GetAxis(2)->GetBinCenter(ibin2);
2579 THnSparse * his2= his3->Projection(2,idim); //projected histogram according selection 2
2580 Int_t nbins1 = his2->GetAxis(1)->GetNbins();
2581 Int_t first1 = his2->GetAxis(1)->GetFirst();
2582 Int_t last1 = his2->GetAxis(1)->GetLast();
2583 for (Int_t ibin1=first1; ibin1<last1; ibin1++){ //axis 1 - theta
2585 Double_t x1= his2->GetAxis(1)->GetBinCenter(ibin1);
2586 his2->GetAxis(1)->SetRange(TMath::Max(ibin1-1,1),TMath::Min(ibin1+1,nbins1));
2587 if (TMath::Abs(x1)<0.1){
2588 if (x1<0) his2->GetAxis(1)->SetRange(TMath::Max(ibin1-1,1),TMath::Min(ibin1,nbins1));
2589 if (x1>0) his2->GetAxis(1)->SetRange(TMath::Max(ibin1,1),TMath::Min(ibin1+1,nbins1));
2591 if (TMath::Abs(x1)<0.06){
2592 his2->GetAxis(1)->SetRange(TMath::Max(ibin1,1),TMath::Min(ibin1,nbins1));
2594 TH1 * hisDelta = his2->Projection(0);
2596 Double_t entries = hisDelta->GetEntries();
2597 Double_t mean=0, rms=0;
2598 if (entries>kMinEntries){
2599 mean = hisDelta->GetMean();
2600 rms = hisDelta->GetRMS();
2602 Double_t sector = 9.*x2/TMath::Pi();
2603 if (sector<0) sector+=18;
2604 Double_t dsec = sector-Int_t(sector)-0.5;
2606 (*pcstream)<<hname<<
2611 "z="<<z<< // dummy z
2613 "entries="<<entries<<
2616 "refX="<<refX<< // track matching refernce plane
2622 //printf("%f\t%f\t%f\t%f\t%f\n",x3,x2,x1, entries,mean);
2633 void AliTPCCorrection::MakeDistortionMapCosmic(THnSparse * hisInput, TTreeSRedirector * const pcstream, const char* hname, Int_t run, Float_t refX, Int_t type){
2635 // make a distortion map out ou fthe residual histogram
2636 // Results are written to the debug streamer - pcstream
2638 // his0 - input (4D) residual histogram
2639 // pcstream - file to write the tree
2641 // refX - track matching reference X
2642 // type - 0- y 1-z,2 -snp, 3-theta, 4=1/pt
2643 // marian.ivanov@cern.ch
2646 // Collection name='TObjArray', class='TObjArray', size=16
2647 // 0. OBJ: TAxis #Delta #Delta
2648 // 1. OBJ: TAxis N_{cl} N_{cl}
2649 // 2. OBJ: TAxis dca_{r} (cm) dca_{r} (cm)
2650 // 3. OBJ: TAxis z (cm) z (cm)
2651 // 4. OBJ: TAxis sin(#phi) sin(#phi)
2652 // 5. OBJ: TAxis tan(#theta) tan(#theta)
2653 // 6. OBJ: TAxis 1/pt (1/GeV) 1/pt (1/GeV)
2654 // 7. OBJ: TAxis pt (GeV) pt (GeV)
2655 // 8. OBJ: TAxis alpha alpha
2656 const Int_t kMinEntries=10;
2658 // 1. make default selections
2661 Int_t idim0[4]={0 , 5, 8, 3}; // delta, theta, alpha, z
2662 hisInput->GetAxis(1)->SetRangeUser(110,190); //long tracks
2663 hisInput->GetAxis(2)->SetRangeUser(-10,35); //tracks close to beam pipe
2664 hisInput->GetAxis(4)->SetRangeUser(-0.3,0.3); //small snp at TPC entrance
2665 hisInput->GetAxis(7)->SetRangeUser(3,100); //"high pt tracks"
2666 hisDelta= hisInput->Projection(0);
2667 hisInput->GetAxis(0)->SetRangeUser(-6.*hisDelta->GetRMS(), +6.*hisDelta->GetRMS());
2669 THnSparse *his0= hisInput->Projection(4,idim0);
2671 // 2. Get mean in diferent bins
2673 Int_t nbins1=his0->GetAxis(1)->GetNbins();
2674 Int_t first1=his0->GetAxis(1)->GetFirst();
2675 Int_t last1 =his0->GetAxis(1)->GetLast();
2677 Double_t bz=AliTrackerBase::GetBz();
2678 Int_t idim[4]={0,1, 2, 3}; // delta, theta,alpha,z
2680 for (Int_t ibin1=first1; ibin1<=last1; ibin1++){ //axis 1 - theta
2682 Double_t x1= his0->GetAxis(1)->GetBinCenter(ibin1);
2683 his0->GetAxis(1)->SetRange(TMath::Max(ibin1-1,1),TMath::Min(ibin1+1,nbins1));
2685 THnSparse * his1 = his0->Projection(4,idim); // projected histogram according range1
2686 Int_t nbins3 = his1->GetAxis(3)->GetNbins();
2687 Int_t first3 = his1->GetAxis(3)->GetFirst();
2688 Int_t last3 = his1->GetAxis(3)->GetLast();
2690 for (Int_t ibin3=first3-1; ibin3<=last3; ibin3+=1){ // axis 3 - z at "vertex"
2691 his1->GetAxis(3)->SetRange(TMath::Max(ibin3-1,1),TMath::Min(ibin3+1,nbins3));
2692 Double_t x3= his1->GetAxis(3)->GetBinCenter(ibin3);
2694 his1->GetAxis(3)->SetRangeUser(-1,1);
2697 THnSparse * his3= his1->Projection(4,idim); //projected histogram according selection 3
2698 Int_t nbins2 = his3->GetAxis(2)->GetNbins();
2699 Int_t first2 = his3->GetAxis(2)->GetFirst();
2700 Int_t last2 = his3->GetAxis(2)->GetLast();
2702 for (Int_t ibin2=first2; ibin2<=last2; ibin2+=1){
2703 his3->GetAxis(2)->SetRange(TMath::Max(ibin2-1,1),TMath::Min(ibin2+1,nbins2));
2704 Double_t x2= his3->GetAxis(2)->GetBinCenter(ibin2);
2705 hisDelta = his3->Projection(0);
2707 Double_t entries = hisDelta->GetEntries();
2708 Double_t mean=0, rms=0;
2709 if (entries>kMinEntries){
2710 mean = hisDelta->GetMean();
2711 rms = hisDelta->GetRMS();
2713 Double_t sector = 9.*x2/TMath::Pi();
2714 if (sector<0) sector+=18;
2715 Double_t dsec = sector-Int_t(sector)-0.5;
2716 Double_t snp=0; // dummy snp - equal 0
2717 (*pcstream)<<hname<<
2719 "bz="<<bz<< // magnetic field
2720 "theta="<<x1<< // theta
2721 "phi="<<x2<< // phi (alpha)
2722 "z="<<x3<< // z at "vertex"
2723 "snp="<<snp<< // dummy snp
2724 "entries="<<entries<< // entries in bin
2725 "mean="<<mean<< // mean
2727 "refX="<<refX<< // track matching refernce plane
2728 "type="<<type<< // parameter type
2729 "sector="<<sector<< // sector
2730 "dsec="<<dsec<< // dummy delta sector
2733 printf("%f\t%f\t%f\t%f\t%f\n",x1,x3,x2, entries,mean);
2744 void AliTPCCorrection::MakeDistortionMapSector(THnSparse * hisInput, TTreeSRedirector * const pcstream, const char* hname, Int_t run, Int_t type){
2746 // make a distortion map out of the residual histogram
2747 // Results are written to the debug streamer - pcstream
2749 // his0 - input (4D) residual histogram
2750 // pcstream - file to write the tree
2752 // type - 0- y 1-z,2 -snp, 3-theta
2753 // marian.ivanov@cern.ch
2755 //Collection name='TObjArray', class='TObjArray', size=16
2756 //0 OBJ: TAxis delta delta
2757 //1 OBJ: TAxis phi phi
2758 //2 OBJ: TAxis localX localX
2759 //3 OBJ: TAxis kY kY
2760 //4 OBJ: TAxis kZ kZ
2761 //5 OBJ: TAxis is1 is1
2762 //6 OBJ: TAxis is0 is0
2764 //8. OBJ: TAxis IsPrimary IsPrimary
2766 const Int_t kMinEntries=10;
2767 THnSparse * hisSector0=0;
2768 TH1 * htemp=0; // histogram to calculate mean value of parameter
2769 Double_t bz=AliTrackerBase::GetBz();
2772 // Loop over pair of sector:
2783 for (Int_t isec0=0; isec0<72; isec0++){
2784 Int_t index0[9]={0, 4, 3, 7, 1, 2, 5, 6,8}; //regroup indeces
2786 //hisInput->GetAxis(8)->SetRangeUser(-0.1,0.4); // select secondaries only ? - get out later ?
2787 hisInput->GetAxis(6)->SetRangeUser(isec0-0.1,isec0+0.1);
2788 hisSector0=hisInput->Projection(7,index0);
2791 for (Int_t isec1=isec0+1; isec1<72; isec1++){
2792 //if (isec1!=isec0+36) continue;
2793 if ( TMath::Abs((isec0%18)-(isec1%18))>1.5 && TMath::Abs((isec0%18)-(isec1%18))<16.5) continue;
2794 printf("Sectors %d\t%d\n",isec1,isec0);
2795 hisSector0->GetAxis(6)->SetRangeUser(isec1-0.1,isec1+0.1);
2796 TH1 * hisX=hisSector0->Projection(5);
2797 Double_t refX= hisX->GetMean();
2799 TH1 *hisDelta=hisSector0->Projection(0);
2800 Double_t dmean = hisDelta->GetMean();
2801 Double_t drms = hisDelta->GetRMS();
2802 hisSector0->GetAxis(0)->SetRangeUser(dmean-5.*drms, dmean+5.*drms);
2805 // 1. make default selections
2807 Int_t idim0[5]={0 , 1, 2, 3, 4}; // {delta, theta, snp, z, phi }
2808 THnSparse *hisSector1= hisSector0->Projection(5,idim0);
2810 // 2. Get mean in diferent bins
2812 Int_t idim[5]={0, 1, 2, 3, 4}; // {delta, theta-1,snp-2 ,z-3, phi-4}
2814 // Int_t nbinsPhi=hisSector1->GetAxis(4)->GetNbins();
2815 Int_t firstPhi=hisSector1->GetAxis(4)->GetFirst();
2816 Int_t lastPhi =hisSector1->GetAxis(4)->GetLast();
2818 for (Int_t ibinPhi=firstPhi; ibinPhi<=lastPhi; ibinPhi+=1){ //axis 4 - phi
2822 Double_t xPhi= hisSector1->GetAxis(4)->GetBinCenter(ibinPhi);
2823 Double_t psec = (9*xPhi/TMath::Pi());
2824 if (psec<0) psec+=18;
2825 Bool_t isOK0=kFALSE;
2826 Bool_t isOK1=kFALSE;
2827 if (TMath::Abs(psec-isec0%18-0.5)<1. || TMath::Abs(psec-isec0%18-17.5)<1.) isOK0=kTRUE;
2828 if (TMath::Abs(psec-isec1%18-0.5)<1. || TMath::Abs(psec-isec1%18-17.5)<1.) isOK1=kTRUE;
2829 if (!isOK0) continue;
2830 if (!isOK1) continue;
2832 hisSector1->GetAxis(4)->SetRange(TMath::Max(ibinPhi-2,firstPhi),TMath::Min(ibinPhi+2,lastPhi));
2833 if (isec1!=isec0+36) {
2834 hisSector1->GetAxis(4)->SetRange(TMath::Max(ibinPhi-3,firstPhi),TMath::Min(ibinPhi+3,lastPhi));
2837 htemp = hisSector1->Projection(4);
2838 xPhi=htemp->GetMean();
2840 THnSparse * hisPhi = hisSector1->Projection(4,idim);
2841 //Int_t nbinsZ = hisPhi->GetAxis(3)->GetNbins();
2842 Int_t firstZ = hisPhi->GetAxis(3)->GetFirst();
2843 Int_t lastZ = hisPhi->GetAxis(3)->GetLast();
2845 for (Int_t ibinZ=firstZ; ibinZ<=lastZ; ibinZ+=1){ // axis 3 - z
2849 hisPhi->GetAxis(3)->SetRange(TMath::Max(ibinZ,firstZ),TMath::Min(ibinZ,lastZ));
2850 if (isec1!=isec0+36) {
2851 hisPhi->GetAxis(3)->SetRange(TMath::Max(ibinZ-1,firstZ),TMath::Min(ibinZ-1,lastZ));
2853 htemp = hisPhi->Projection(3);
2854 Double_t xZ= htemp->GetMean();
2856 THnSparse * hisZ= hisPhi->Projection(3,idim);
2857 //projected histogram according selection 3 -z
2860 //Int_t nbinsSnp = hisZ->GetAxis(2)->GetNbins();
2861 Int_t firstSnp = hisZ->GetAxis(2)->GetFirst();
2862 Int_t lastSnp = hisZ->GetAxis(2)->GetLast();
2863 for (Int_t ibinSnp=firstSnp; ibinSnp<=lastSnp; ibinSnp+=2){ // axis 2 - snp
2867 hisZ->GetAxis(2)->SetRange(TMath::Max(ibinSnp-1,firstSnp),TMath::Min(ibinSnp+1,lastSnp));
2868 if (isec1!=isec0+36) {
2869 hisZ->GetAxis(2)->SetRange(TMath::Max(ibinSnp-2,firstSnp),TMath::Min(ibinSnp+2,lastSnp));
2871 htemp = hisZ->Projection(2);
2872 Double_t xSnp= htemp->GetMean();
2874 THnSparse * hisSnp= hisZ->Projection(2,idim);
2875 //projected histogram according selection 2 - snp
2877 //Int_t nbinsTheta = hisSnp->GetAxis(1)->GetNbins();
2878 Int_t firstTheta = hisSnp->GetAxis(1)->GetFirst();
2879 Int_t lastTheta = hisSnp->GetAxis(1)->GetLast();
2881 for (Int_t ibinTheta=firstTheta; ibinTheta<=lastTheta; ibinTheta+=2){ // axis1 theta
2884 hisSnp->GetAxis(1)->SetRange(TMath::Max(ibinTheta-2,firstTheta),TMath::Min(ibinTheta+2,lastTheta));
2885 if (isec1!=isec0+36) {
2886 hisSnp->GetAxis(1)->SetRange(TMath::Max(ibinTheta-3,firstTheta),TMath::Min(ibinTheta+3,lastTheta));
2888 htemp = hisSnp->Projection(1);
2889 Double_t xTheta=htemp->GetMean();
2891 hisDelta = hisSnp->Projection(0);
2893 Double_t entries = hisDelta->GetEntries();
2894 Double_t mean=0, rms=0;
2895 if (entries>kMinEntries){
2896 mean = hisDelta->GetMean();
2897 rms = hisDelta->GetRMS();
2899 Double_t sector = 9.*xPhi/TMath::Pi();
2900 if (sector<0) sector+=18;
2901 Double_t dsec = sector-Int_t(sector)-0.5;
2902 Int_t dtype=1; // TPC alignment type
2903 (*pcstream)<<hname<<
2905 "bz="<<bz<< // magnetic field
2906 "ptype="<<type<< // parameter type
2907 "dtype="<<dtype<< // parameter type
2908 "isec0="<<isec0<< // sector 0
2909 "isec1="<<isec1<< // sector 1
2910 "sector="<<sector<< // sector as float
2911 "dsec="<<dsec<< // delta sector
2913 "theta="<<xTheta<< // theta
2914 "phi="<<xPhi<< // phi (alpha)
2916 "snp="<<xSnp<< // snp
2918 "entries="<<entries<< // entries in bin
2919 "mean="<<mean<< // mean
2920 "rms="<<rms<< // rms
2921 "refX="<<refX<< // track matching reference plane
2924 printf("%d\t%d\t%f\t%f\t%f\t%f\t%f\t%f\n",isec0, isec1, xPhi,xZ,xSnp, xTheta, entries,mean);
2945 void AliTPCCorrection::StoreInOCDB(Int_t startRun, Int_t endRun, const char *comment){
2947 // Store object in the OCDB
2948 // By default the object is stored in the current directory
2949 // default comment consit of user name and the date
2951 TString ocdbStorage="";
2952 ocdbStorage+="local://"+gSystem->GetFromPipe("pwd")+"/OCDB";
2953 AliCDBMetaData *metaData= new AliCDBMetaData();
2954 metaData->SetObjectClassName("AliTPCCorrection");
2955 metaData->SetResponsible("Marian Ivanov");
2956 metaData->SetBeamPeriod(1);
2957 metaData->SetAliRootVersion("05-25-01"); //root version
2958 TString userName=gSystem->GetFromPipe("echo $USER");
2959 TString date=gSystem->GetFromPipe("date");
2961 if (!comment) metaData->SetComment(Form("Space point distortion calibration\n User: %s\n Data%s",userName.Data(),date.Data()));
2962 if (comment) metaData->SetComment(comment);
2964 id1=new AliCDBId("TPC/Calib/Correction", startRun, endRun);
2965 AliCDBStorage* gStorage = AliCDBManager::Instance()->GetStorage(ocdbStorage);
2966 gStorage->Put(this, (*id1), metaData);
2970 void AliTPCCorrection::FastSimDistortedVertex(Double_t orgVertex[3], Int_t nTracks, AliESDVertex &aV, AliESDVertex &avOrg, AliESDVertex &cV, AliESDVertex &cvOrg, TTreeSRedirector * const pcstream, Double_t etaCuts){
2972 // Fast method to simulate the influence of the given distortion on the vertex reconstruction
2975 AliMagF* magF= (AliMagF*)TGeoGlobalMagField::Instance()->GetField();
2976 if (!magF) AliError("Magneticd field - not initialized");
2977 Double_t bz = magF->SolenoidField(); //field in kGauss
2978 printf("bz: %f\n",bz);
2979 AliVertexerTracks *vertexer = new AliVertexerTracks(bz); // bz in kGauss
2981 TObjArray aTrk; // Original Track array of Aside
2982 TObjArray daTrk; // Distorted Track array of A side
2983 UShort_t *aId = new UShort_t[nTracks]; // A side Track ID
2986 UShort_t *cId = new UShort_t [nTracks];
2988 Double_t mass = TDatabasePDG::Instance()->GetParticle("pi+")->Mass();
2989 TF1 fpt("fpt",Form("x*(1+(sqrt(x*x+%f^2)-%f)/([0]*[1]))^(-[0])",mass,mass),0.4,10);
2990 fpt.SetParameters(7.24,0.120);
2992 for(Int_t nt=0; nt<nTracks; nt++){
2993 Double_t phi = gRandom->Uniform(0.0, 2*TMath::Pi());
2994 Double_t eta = gRandom->Uniform(-etaCuts, etaCuts);
2995 Double_t pt = fpt.GetRandom(); // momentum for f1
2996 // printf("phi %lf eta %lf pt %lf\n",phi,eta,pt);
2998 if(gRandom->Rndm() < 0.5){
3004 Double_t theta = 2*TMath::ATan(TMath::Exp(-eta))-TMath::Pi()/2.;
3006 pxyz[0]=pt*TMath::Cos(phi);
3007 pxyz[1]=pt*TMath::Sin(phi);
3008 pxyz[2]=pt*TMath::Tan(theta);
3009 Double_t cv[21]={0};
3010 AliExternalTrackParam *t= new AliExternalTrackParam(orgVertex, pxyz, cv, sign);
3014 AliExternalTrackParam *td = FitDistortedTrack(*t, refX, dir, NULL);
3016 if (pcstream) (*pcstream)<<"track"<<
3022 if(( eta>0.07 )&&( eta<etaCuts )) { // - log(tan(0.5*theta)), theta = 0.5*pi - ATan(5.0/80.0)
3026 Int_t nn=aTrk.GetEntriesFast();
3029 }else if(( eta<-0.07 )&&( eta>-etaCuts )){
3033 Int_t nn=cTrk.GetEntriesFast();
3038 }// end of track loop
3040 vertexer->SetTPCMode();
3041 vertexer->SetConstraintOff();
3043 aV = *((AliESDVertex*)vertexer->FindPrimaryVertex(&daTrk,aId));
3044 avOrg = *((AliESDVertex*)vertexer->FindPrimaryVertex(&aTrk,aId));
3045 cV = *((AliESDVertex*)vertexer->FindPrimaryVertex(&dcTrk,cId));
3046 cvOrg = *((AliESDVertex*)vertexer->FindPrimaryVertex(&cTrk,cId));
3047 if (pcstream) (*pcstream)<<"vertex"<<
3048 "x="<<orgVertex[0]<<
3049 "y="<<orgVertex[1]<<
3050 "z="<<orgVertex[2]<<
3051 "av.="<<&aV<< // distorted vertex A side
3052 "cv.="<<&cV<< // distroted vertex C side
3053 "avO.="<<&avOrg<< // original vertex A side
3060 void AliTPCCorrection::AddVisualCorrection(AliTPCCorrection* corr, Int_t position){
3062 // make correction available for visualization using
3063 // TFormula, TFX and TTree::Draw
3064 // important in order to check corrections and also compute dervied variables
3065 // e.g correction partial derivatives
3067 // NOTE - class is not owner of correction
3069 if (!fgVisualCorrection) fgVisualCorrection=new TObjArray(10000);
3070 if (position>=fgVisualCorrection->GetEntriesFast())
3071 fgVisualCorrection->Expand((position+10)*2);
3072 fgVisualCorrection->AddAt(corr, position);
3075 AliTPCCorrection* AliTPCCorrection::GetVisualCorrection(Int_t position) {
3077 // Get visula correction registered at index=position
3079 return fgVisualCorrection? (AliTPCCorrection*)fgVisualCorrection->At(position):0;
3084 Double_t AliTPCCorrection::GetCorrSector(Double_t sector, Double_t r, Double_t kZ, Int_t axisType, Int_t corrType){
3086 // calculate the correction at given position - check the geffCorr
3088 // corrType return values
3094 if (!fgVisualCorrection) return 0;
3095 AliTPCCorrection *corr = (AliTPCCorrection*)fgVisualCorrection->At(corrType);
3096 if (!corr) return 0;
3098 Double_t phi=sector*TMath::Pi()/9.;
3099 Double_t gx = r*TMath::Cos(phi);
3100 Double_t gy = r*TMath::Sin(phi);
3102 Int_t nsector=(gz>=0) ? 0:18;
3106 Float_t distPoint[3]={static_cast<Float_t>(gx),static_cast<Float_t>(gy),static_cast<Float_t>(gz)};
3107 corr->DistortPoint(distPoint, nsector);
3108 Double_t r0=TMath::Sqrt(gx*gx+gy*gy);
3109 Double_t r1=TMath::Sqrt(distPoint[0]*distPoint[0]+distPoint[1]*distPoint[1]);
3110 Double_t phi0=TMath::ATan2(gy,gx);
3111 Double_t phi1=TMath::ATan2(distPoint[1],distPoint[0]);
3112 if (axisType==0) return r1-r0;
3113 if (axisType==1) return (phi1-phi0)*r0;
3114 if (axisType==2) return distPoint[2]-gz;
3115 if (axisType==3) return (TMath::Cos(phi)*(distPoint[0]-gx)+ TMath::Cos(phi)*(distPoint[1]-gy));
3119 Double_t AliTPCCorrection::GetCorrXYZ(Double_t gx, Double_t gy, Double_t gz, Int_t axisType, Int_t corrType){
3121 // return correction at given x,y,z
3123 if (!fgVisualCorrection) return 0;
3124 AliTPCCorrection *corr = (AliTPCCorrection*)fgVisualCorrection->At(corrType);
3125 if (!corr) return 0;
3126 Double_t phi0= TMath::ATan2(gy,gx);
3127 Int_t nsector=(gz>=0) ? 0:18;
3128 Float_t distPoint[3]={static_cast<Float_t>(gx),static_cast<Float_t>(gy),static_cast<Float_t>(gz)};
3129 corr->CorrectPoint(distPoint, nsector);
3130 Double_t r0=TMath::Sqrt(gx*gx+gy*gy);
3131 Double_t r1=TMath::Sqrt(distPoint[0]*distPoint[0]+distPoint[1]*distPoint[1]);
3132 Double_t phi1=TMath::ATan2(distPoint[1],distPoint[0]);
3133 if (axisType==0) return r1-r0;
3134 if (axisType==1) return (phi1-phi0)*r0;
3135 if (axisType==2) return distPoint[2]-gz;
3139 Double_t AliTPCCorrection::GetCorrXYZDz(Double_t gx, Double_t gy, Double_t gz, Int_t axisType, Int_t corrType,Double_t delta){
3141 // return correction at given x,y,z
3143 if (!fgVisualCorrection) return 0;
3144 AliTPCCorrection *corr = (AliTPCCorrection*)fgVisualCorrection->At(corrType);
3145 if (!corr) return 0;
3146 Double_t phi0= TMath::ATan2(gy,gx);
3147 Int_t nsector=(gz>=0) ? 0:18;
3148 Float_t distPoint[3]={static_cast<Float_t>(gx),static_cast<Float_t>(gy),static_cast<Float_t>(gz)};
3149 Float_t dxyz[3]={static_cast<Float_t>(gx),static_cast<Float_t>(gy),static_cast<Float_t>(gz)};
3151 corr->GetCorrectionDz(distPoint, nsector,dxyz,delta);
3152 distPoint[0]+=dxyz[0];
3153 distPoint[1]+=dxyz[1];
3154 distPoint[2]+=dxyz[2];
3155 Double_t r0=TMath::Sqrt(gx*gx+gy*gy);
3156 Double_t r1=TMath::Sqrt(distPoint[0]*distPoint[0]+distPoint[1]*distPoint[1]);
3157 Double_t phi1=TMath::ATan2(distPoint[1],distPoint[0]);
3158 if (axisType==0) return r1-r0;
3159 if (axisType==1) return (phi1-phi0)*r0;
3160 if (axisType==2) return distPoint[2]-gz;
3164 Double_t AliTPCCorrection::GetCorrXYZIntegrateZ(Double_t gx, Double_t gy, Double_t gz, Int_t axisType, Int_t corrType,Double_t delta){
3166 // return correction at given x,y,z
3168 if (!fgVisualCorrection) return 0;
3169 AliTPCCorrection *corr = (AliTPCCorrection*)fgVisualCorrection->At(corrType);
3170 if (!corr) return 0;
3171 Double_t phi0= TMath::ATan2(gy,gx);
3172 Int_t nsector=(gz>=0) ? 0:18;
3173 Float_t distPoint[3]={static_cast<Float_t>(gx),static_cast<Float_t>(gy),static_cast<Float_t>(gz)};
3174 Float_t dxyz[3]={static_cast<Float_t>(gx),static_cast<Float_t>(gy),static_cast<Float_t>(gz)};
3176 corr->GetCorrectionIntegralDz(distPoint, nsector,dxyz,delta);
3177 distPoint[0]+=dxyz[0];
3178 distPoint[1]+=dxyz[1];
3179 distPoint[2]+=dxyz[2];
3180 Double_t r0=TMath::Sqrt(gx*gx+gy*gy);
3181 Double_t r1=TMath::Sqrt(distPoint[0]*distPoint[0]+distPoint[1]*distPoint[1]);
3182 Double_t phi1=TMath::ATan2(distPoint[1],distPoint[0]);
3183 if (axisType==0) return r1-r0;
3184 if (axisType==1) return (phi1-phi0)*r0;
3185 if (axisType==2) return distPoint[2]-gz;
3190 Double_t AliTPCCorrection::GetDistXYZ(Double_t gx, Double_t gy, Double_t gz, Int_t axisType, Int_t corrType){
3192 // return correction at given x,y,z
3194 if (!fgVisualCorrection) return 0;
3195 AliTPCCorrection *corr = (AliTPCCorrection*)fgVisualCorrection->At(corrType);
3196 if (!corr) return 0;
3197 Double_t phi0= TMath::ATan2(gy,gx);
3198 Int_t nsector=(gz>=0) ? 0:18;
3199 Float_t distPoint[3]={static_cast<Float_t>(gx),static_cast<Float_t>(gy),static_cast<Float_t>(gz)};
3200 corr->DistortPoint(distPoint, nsector);
3201 Double_t r0=TMath::Sqrt(gx*gx+gy*gy);
3202 Double_t r1=TMath::Sqrt(distPoint[0]*distPoint[0]+distPoint[1]*distPoint[1]);
3203 Double_t phi1=TMath::ATan2(distPoint[1],distPoint[0]);
3204 if (axisType==0) return r1-r0;
3205 if (axisType==1) return (phi1-phi0)*r0;
3206 if (axisType==2) return distPoint[2]-gz;
3210 Double_t AliTPCCorrection::GetDistXYZDz(Double_t gx, Double_t gy, Double_t gz, Int_t axisType, Int_t corrType,Double_t delta){
3212 // return correction at given x,y,z
3214 if (!fgVisualCorrection) return 0;
3215 AliTPCCorrection *corr = (AliTPCCorrection*)fgVisualCorrection->At(corrType);
3216 if (!corr) return 0;
3217 Double_t phi0= TMath::ATan2(gy,gx);
3218 Int_t nsector=(gz>=0) ? 0:18;
3219 Float_t distPoint[3]={static_cast<Float_t>(gx),static_cast<Float_t>(gy),static_cast<Float_t>(gz)};
3220 Float_t dxyz[3]={static_cast<Float_t>(gx),static_cast<Float_t>(gy),static_cast<Float_t>(gz)};
3222 corr->GetDistortionDz(distPoint, nsector,dxyz,delta);
3223 distPoint[0]+=dxyz[0];
3224 distPoint[1]+=dxyz[1];
3225 distPoint[2]+=dxyz[2];
3226 Double_t r0=TMath::Sqrt(gx*gx+gy*gy);
3227 Double_t r1=TMath::Sqrt(distPoint[0]*distPoint[0]+distPoint[1]*distPoint[1]);
3228 Double_t phi1=TMath::ATan2(distPoint[1],distPoint[0]);
3229 if (axisType==0) return r1-r0;
3230 if (axisType==1) return (phi1-phi0)*r0;
3231 if (axisType==2) return distPoint[2]-gz;
3235 Double_t AliTPCCorrection::GetDistXYZIntegrateZ(Double_t gx, Double_t gy, Double_t gz, Int_t axisType, Int_t corrType,Double_t delta){
3237 // return correction at given x,y,z
3239 if (!fgVisualCorrection) return 0;
3240 AliTPCCorrection *corr = (AliTPCCorrection*)fgVisualCorrection->At(corrType);
3241 if (!corr) return 0;
3242 Double_t phi0= TMath::ATan2(gy,gx);
3243 Int_t nsector=(gz>=0) ? 0:18;
3244 Float_t distPoint[3]={static_cast<Float_t>(gx),static_cast<Float_t>(gy),static_cast<Float_t>(gz)};
3245 Float_t dxyz[3]={static_cast<Float_t>(gx),static_cast<Float_t>(gy),static_cast<Float_t>(gz)};
3247 corr->GetDistortionIntegralDz(distPoint, nsector,dxyz,delta);
3248 distPoint[0]+=dxyz[0];
3249 distPoint[1]+=dxyz[1];
3250 distPoint[2]+=dxyz[2];
3251 Double_t r0=TMath::Sqrt(gx*gx+gy*gy);
3252 Double_t r1=TMath::Sqrt(distPoint[0]*distPoint[0]+distPoint[1]*distPoint[1]);
3253 Double_t phi1=TMath::ATan2(distPoint[1],distPoint[0]);
3254 if (axisType==0) return r1-r0;
3255 if (axisType==1) return (phi1-phi0)*r0;
3256 if (axisType==2) return distPoint[2]-gz;
3262 void AliTPCCorrection::MakeLaserDistortionTree(TTree* tree, TObjArray */*corrArray*/, Int_t /*itype*/){
3264 // Make a laser fit tree for global minimization
3266 AliTPCcalibDB* calib=AliTPCcalibDB::Instance();
3267 AliTPCCorrection * correction = calib->GetTPCComposedCorrection();
3268 if (!correction) correction = calib->GetTPCComposedCorrection(AliTrackerBase::GetBz());
3269 correction->AddVisualCorrection(correction,0); //register correction
3271 // AliTPCTransform *transform = AliTPCcalibDB::Instance()->GetTransform() ;
3272 //AliTPCParam *param = AliTPCcalibDB::Instance()->GetParameters();
3274 const Double_t cutErrY=0.05;
3275 const Double_t kSigmaCut=4;
3276 // const Double_t cutErrZ=0.03;
3277 const Double_t kEpsilon=0.00000001;
3278 // const Double_t kMaxDist=1.; // max distance - space correction
3283 AliTPCLaserTrack *ltr=0;
3284 AliTPCLaserTrack::LoadTracks();
3285 tree->SetBranchAddress("dY.",&vecdY);
3286 tree->SetBranchAddress("dZ.",&vecdZ);
3287 tree->SetBranchAddress("eY.",&veceY);
3288 tree->SetBranchAddress("eZ.",&veceZ);
3289 tree->SetBranchAddress("LTr.",<r);
3290 Int_t entries= tree->GetEntries();
3291 TTreeSRedirector *pcstream= new TTreeSRedirector("distortionLaser_0.root");
3292 Double_t bz=AliTrackerBase::GetBz();
3294 // Double_t globalXYZ[3];
3295 //Double_t globalXYZCorr[3];
3296 for (Int_t ientry=0; ientry<entries; ientry++){
3297 tree->GetEntry(ientry);
3298 if (!ltr->GetVecGX()){
3299 ltr->UpdatePoints();
3304 printf("Entry\t%d\n",ientry);
3305 for (Int_t irow0=0; irow0<158; irow0+=1){
3307 TLinearFitter fitter10(4,"hyp3");
3308 TLinearFitter fitter5(2,"hyp1");
3309 Int_t sector= (Int_t)(*ltr->GetVecSec())[irow0];
3310 if (sector<0) continue;
3311 //if (TMath::Abs(vecdY->GetMatrixArray()[irow0])<kEpsilon) continue;
3313 Double_t refX= (*ltr->GetVecLX())[irow0];
3314 Int_t firstRow1 = TMath::Max(irow0-10,0);
3315 Int_t lastRow1 = TMath::Min(irow0+10,158);
3316 Double_t padWidth=(irow0<64)?0.4:0.6;
3317 // make long range fit
3318 for (Int_t irow1=firstRow1; irow1<=lastRow1; irow1++){
3319 if (TMath::Abs((*ltr->GetVecSec())[irow1]-sector)>kEpsilon) continue;
3320 if (veceY->GetMatrixArray()[irow1]>cutErrY) continue;
3321 if (TMath::Abs(vecdY->GetMatrixArray()[irow1])<kEpsilon) continue;
3322 Double_t idealX= (*ltr->GetVecLX())[irow1];
3323 Double_t idealY= (*ltr->GetVecLY())[irow1];
3324 // Double_t idealZ= (*ltr->GetVecLZ())[irow1];
3325 Double_t gx= (*ltr->GetVecGX())[irow1];
3326 Double_t gy= (*ltr->GetVecGY())[irow1];
3327 Double_t gz= (*ltr->GetVecGZ())[irow1];
3328 Double_t measY=(*vecdY)[irow1]+idealY;
3329 Double_t deltaR = GetCorrXYZ(gx, gy, gz, 0,0);
3330 // deltaR = R distorted -R ideal
3331 Double_t xxx[4]={idealX+deltaR-refX,TMath::Cos(idealY/padWidth), TMath::Sin(idealY/padWidth)};
3332 fitter10.AddPoint(xxx,measY,1);
3335 Double_t rms10=0;//TMath::Sqrt(fitter10.GetChisquare()/(fitter10.GetNpoints()-4));
3336 Double_t mean10 =0;// fitter10.GetParameter(0);
3337 Double_t slope10 =0;// fitter10.GetParameter(0);
3338 Double_t cosPart10 = 0;// fitter10.GetParameter(2);
3339 Double_t sinPart10 =0;// fitter10.GetParameter(3);
3341 if (fitter10.GetNpoints()>10){
3343 rms10=TMath::Sqrt(fitter10.GetChisquare()/(fitter10.GetNpoints()-4));
3344 mean10 = fitter10.GetParameter(0);
3345 slope10 = fitter10.GetParameter(1);
3346 cosPart10 = fitter10.GetParameter(2);
3347 sinPart10 = fitter10.GetParameter(3);
3349 // make short range fit
3351 for (Int_t irow1=firstRow1+5; irow1<=lastRow1-5; irow1++){
3352 if (TMath::Abs((*ltr->GetVecSec())[irow1]-sector)>kEpsilon) continue;
3353 if (veceY->GetMatrixArray()[irow1]>cutErrY) continue;
3354 if (TMath::Abs(vecdY->GetMatrixArray()[irow1])<kEpsilon) continue;
3355 Double_t idealX= (*ltr->GetVecLX())[irow1];
3356 Double_t idealY= (*ltr->GetVecLY())[irow1];
3357 // Double_t idealZ= (*ltr->GetVecLZ())[irow1];
3358 Double_t gx= (*ltr->GetVecGX())[irow1];
3359 Double_t gy= (*ltr->GetVecGY())[irow1];
3360 Double_t gz= (*ltr->GetVecGZ())[irow1];
3361 Double_t measY=(*vecdY)[irow1]+idealY;
3362 Double_t deltaR = GetCorrXYZ(gx, gy, gz, 0,0);
3363 // deltaR = R distorted -R ideal
3364 Double_t expY= mean10+slope10*(idealX+deltaR-refX);
3365 if (TMath::Abs(measY-expY)>kSigmaCut*rms10) continue;
3367 Double_t corr=cosPart10*TMath::Cos(idealY/padWidth)+sinPart10*TMath::Sin(idealY/padWidth);
3368 Double_t xxx[4]={idealX+deltaR-refX,TMath::Cos(idealY/padWidth), TMath::Sin(idealY/padWidth)};
3369 fitter5.AddPoint(xxx,measY-corr,1);
3374 if (fitter5.GetNpoints()<8) isOK=kFALSE;
3376 Double_t rms5=0;//TMath::Sqrt(fitter5.GetChisquare()/(fitter5.GetNpoints()-4));
3377 Double_t offset5 =0;// fitter5.GetParameter(0);
3378 Double_t slope5 =0;// fitter5.GetParameter(0);
3381 rms5=TMath::Sqrt(fitter5.GetChisquare()/(fitter5.GetNpoints()-4));
3382 offset5 = fitter5.GetParameter(0);
3383 slope5 = fitter5.GetParameter(0);
3388 Double_t phi =(*ltr->GetVecPhi())[irow0];
3389 Double_t theta =ltr->GetTgl();
3390 Double_t mean=(vecdY)->GetMatrixArray()[irow0];
3391 Double_t gx=0,gy=0,gz=0;
3392 Double_t snp = (*ltr->GetVecP2())[irow0];
3393 Int_t bundle= ltr->GetBundle();
3394 Int_t id= ltr->GetId();
3395 // Double_t rms = err->GetMatrixArray()[irow];
3397 gx = (*ltr->GetVecGX())[irow0];
3398 gy = (*ltr->GetVecGY())[irow0];
3399 gz = (*ltr->GetVecGZ())[irow0];
3400 Double_t dRrec = GetCorrXYZ(gx, gy, gz, 0,0);
3401 fitter10.GetParameters(fit10);
3402 fitter5.GetParameters(fit5);
3403 Double_t idealY= (*ltr->GetVecLY())[irow0];
3404 Double_t measY=(*vecdY)[irow0]+idealY;
3405 Double_t corr=cosPart10*TMath::Cos(idealY/padWidth)+sinPart10*TMath::Sin(idealY/padWidth);
3406 if (TMath::Max(rms5,rms10)>0.06) isOK=kFALSE;
3408 (*pcstream)<<"fitFull"<< // dumpe also intermediate results
3409 "bz="<<bz<< // magnetic filed used
3410 "dtype="<<dtype<< // detector match type
3411 "ptype="<<ptype<< // parameter type
3412 "theta="<<theta<< // theta
3413 "phi="<<phi<< // phi
3414 "snp="<<snp<< // snp
3417 // // "dsec="<<dsec<<
3418 "refX="<<refX<< // reference radius
3419 "gx="<<gx<< // global position
3420 "gy="<<gy<< // global position
3421 "gz="<<gz<< // global position
3422 "dRrec="<<dRrec<< // delta Radius in reconstruction
3423 "id="<<id<< //bundle