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.SetOmegaTauT1T2(0,1,1); // B=0
64 // Float_t z0 = 1; // at +1 cm -> A side
65 // c2->cd(1); roc.CreateHistoDRinXY(1.,300,300)->Draw("cont4z");
66 // c2->cd(3);roc.CreateHistoDRPhiinXY(1.,300,300)->Draw("cont4z");
67 // c2->cd(5);roc.CreateHistoDZinXY(1.,300,300)->Draw("cont4z");
69 // c2->cd(2);roc.CreateHistoDRinZR(phi0)->Draw("surf2");
70 // c2->cd(4);roc.CreateHistoDRPhiinZR(phi0)->Draw("surf2");
71 // c2->cd(6);roc.CreateHistoDZinZR(phi0)->Draw("surf2");
78 // Date: 27/04/2010 <br>
79 // Authors: Magnus Mager, Stefan Rossegger, Jim Thomas
81 // _________________________________________________________________
84 #include "Riostream.h"
89 #include <TTreeStream.h>
92 #include <TTimeStamp.h>
93 #include <AliCDBStorage.h>
95 #include <AliCDBMetaData.h>
97 #include "AliTPCParamSR.h"
99 #include "AliTPCCorrection.h"
102 #include "AliExternalTrackParam.h"
103 #include "AliTrackPointArray.h"
104 #include "TDatabasePDG.h"
105 #include "AliTrackerBase.h"
106 #include "AliTPCROC.h"
107 #include "THnSparse.h"
109 #include "AliTPCLaserTrack.h"
110 #include "AliESDVertex.h"
111 #include "AliVertexerTracks.h"
112 #include "TDatabasePDG.h"
116 #include "TDatabasePDG.h"
118 #include "AliTPCTransform.h"
119 #include "AliTPCcalibDB.h"
120 #include "AliTPCExB.h"
122 #include "AliTPCRecoParam.h"
123 #include "TLinearFitter.h"
124 #include <AliSysInfo.h>
126 ClassImp(AliTPCCorrection)
129 TObjArray *AliTPCCorrection::fgVisualCorrection=0;
130 // instance of correction for visualization
133 // FIXME: the following values should come from the database
134 const Double_t AliTPCCorrection::fgkTPCZ0 = 249.7; // nominal gating grid position
135 const Double_t AliTPCCorrection::fgkIFCRadius= 83.5; // radius which renders the "18 rod manifold" best -> compare calc. of Jim Thomas
136 // compare gkIFCRadius= 83.05: Mean Radius of the Inner Field Cage ( 82.43 min, 83.70 max) (cm)
137 const Double_t AliTPCCorrection::fgkOFCRadius= 254.5; // Mean Radius of the Outer Field Cage (252.55 min, 256.45 max) (cm)
138 const Double_t AliTPCCorrection::fgkZOffSet = 0.2; // Offset from CE: calculate all distortions closer to CE as if at this point
139 const Double_t AliTPCCorrection::fgkCathodeV = -100000.0; // Cathode Voltage (volts)
140 const Double_t AliTPCCorrection::fgkGG = -70.0; // Gating Grid voltage (volts)
142 const Double_t AliTPCCorrection::fgkdvdE = 0.0024; // [cm/V] drift velocity dependency on the E field (from Magboltz for NeCO2N2 at standard environment)
144 const Double_t AliTPCCorrection::fgkEM = -1.602176487e-19/9.10938215e-31; // charge/mass in [C/kg]
145 const Double_t AliTPCCorrection::fgke0 = 8.854187817e-12; // vacuum permittivity [A·s/(V·m)]
148 AliTPCCorrection::AliTPCCorrection()
149 : TNamed("correction_unity","unity"),fILow(0),fJLow(0),fKLow(0), fT1(1), fT2(1), fIsLocal(kFALSE)
152 // default constructor
154 if (!fgVisualCorrection) fgVisualCorrection= new TObjArray;
156 InitLookUpfulcrums();
160 AliTPCCorrection::AliTPCCorrection(const char *name,const char *title)
161 : TNamed(name,title),fILow(0),fJLow(0),fKLow(0), fT1(1), fT2(1), fIsLocal(kFALSE)
164 // default constructor, that set the name and title
166 if (!fgVisualCorrection) fgVisualCorrection= new TObjArray;
168 InitLookUpfulcrums();
172 AliTPCCorrection::~AliTPCCorrection() {
174 // virtual destructor
178 void AliTPCCorrection::CorrectPoint(Float_t x[],const Short_t roc) {
180 // Corrects the initial coordinates x (cartesian coordinates)
181 // according to the given effect (inherited classes)
182 // roc represents the TPC read out chamber (offline numbering convention)
185 GetCorrection(x,roc,dx);
186 for (Int_t j=0;j<3;++j) x[j]+=dx[j];
189 void AliTPCCorrection::CorrectPoint(const Float_t x[],const Short_t roc,Float_t xp[]) {
191 // Corrects the initial coordinates x (cartesian coordinates) and stores the new
192 // (distorted) coordinates in xp. The distortion is set according to the given effect (inherited classes)
193 // roc represents the TPC read out chamber (offline numbering convention)
196 GetCorrection(x,roc,dx);
197 for (Int_t j=0;j<3;++j) xp[j]=x[j]+dx[j];
200 void AliTPCCorrection::DistortPoint(Float_t x[],const Short_t roc) {
202 // Distorts the initial coordinates x (cartesian coordinates)
203 // according to the given effect (inherited classes)
204 // roc represents the TPC read out chamber (offline numbering convention)
207 GetDistortion(x,roc,dx);
208 for (Int_t j=0;j<3;++j) x[j]+=dx[j];
211 void AliTPCCorrection::DistortPointLocal(Float_t x[],const Short_t roc) {
213 // Distorts the initial coordinates x (cartesian coordinates)
214 // according to the given effect (inherited classes)
215 // roc represents the TPC read out chamber (offline numbering convention)
217 Float_t gxyz[3]={0,0,0};
218 Double_t alpha = TMath::Pi()*(roc%18+0.5)/18;
219 Double_t ca=TMath::Cos(alpha), sa= TMath::Sin(alpha);
220 gxyz[0]= ca*x[0]+sa*x[1];
221 gxyz[1]= -sa*x[0]+ca*x[1];
223 DistortPoint(gxyz,roc);
224 x[0]= ca*gxyz[0]-sa*gxyz[1];
225 x[1]= +sa*gxyz[0]+ca*gxyz[1];
228 void AliTPCCorrection::CorrectPointLocal(Float_t x[],const Short_t roc) {
230 // Distorts the initial coordinates x (cartesian coordinates)
231 // according to the given effect (inherited classes)
232 // roc represents the TPC read out chamber (offline numbering convention)
234 Float_t gxyz[3]={0,0,0};
235 Double_t alpha = TMath::Pi()*(roc%18+0.5)/18;
236 Double_t ca=TMath::Cos(alpha), sa= TMath::Sin(alpha);
237 gxyz[0]= ca*x[0]+sa*x[1];
238 gxyz[1]= -sa*x[0]+ca*x[1];
240 CorrectPoint(gxyz,roc);
241 x[0]= ca*gxyz[0]-sa*gxyz[1];
242 x[1]= sa*gxyz[0]+ca*gxyz[1];
246 void AliTPCCorrection::DistortPoint(const Float_t x[],const Short_t roc,Float_t xp[]) {
248 // Distorts the initial coordinates x (cartesian coordinates) and stores the new
249 // (distorted) coordinates in xp. The distortion is set according to the given effect (inherited classes)
250 // roc represents the TPC read out chamber (offline numbering convention)
253 GetDistortion(x,roc,dx);
254 for (Int_t j=0;j<3;++j) xp[j]=x[j]+dx[j];
257 void AliTPCCorrection::GetCorrection(const Float_t /*x*/[],const Short_t /*roc*/,Float_t dx[]) {
259 // This function delivers the correction values dx in respect to the inital coordinates x
260 // roc represents the TPC read out chamber (offline numbering convention)
261 // Note: The dx is overwritten by the inherited effectice class ...
263 for (Int_t j=0;j<3;++j) { dx[j]=0.; }
266 void AliTPCCorrection::GetDistortion(const Float_t x[],const Short_t roc,Float_t dx[]) {
268 // This function delivers the distortion values dx in respect to the inital coordinates x
269 // roc represents the TPC read out chamber (offline numbering convention)
271 GetCorrection(x,roc,dx);
272 for (Int_t j=0;j<3;++j) dx[j]=-dx[j];
275 void AliTPCCorrection::GetCorrectionDz(const Float_t x[],const Short_t roc,Float_t dx[], Float_t delta) {
276 // author: marian.ivanov@cern.ch
278 // In this (virtual)function calculates the dx'/dz, dy'/dz and dz'/dz at given point (x,y,z)
279 // Generic implementation. Better precision can be acchieved knowing the internal structure
280 // of underlying trasnformation. Derived classes can reimplement it.
281 // To calculate correction is fitted in small neighberhood:
282 // (x+-delta,y+-delta,z+-delta) where delta is an argument
285 // x[] - space point corrdinate
286 // roc - readout chamber identifier (important e.g to do not miss the side of detector)
287 // delta - define the size of neighberhood
289 // dx[] - array {dx'/dz, dy'/dz , dz'/dz }
291 // if (fIsLocal){ //standard implemenation provides the correction/distortion integrated over full drift length
294 // GetCorrection(xyz,roc,dxyz);
296 static TLinearFitter fitx(2,"pol1");
297 static TLinearFitter fity(2,"pol1");
298 static TLinearFitter fitz(2,"pol1");
304 //adjust limits around CE to stay on one side
307 if ((x[2]+zmin*delta)<0){
319 if ((x[2]+zmax*delta)>0){
329 for (Int_t xdelta=-1; xdelta<=1; xdelta++)
330 for (Int_t ydelta=-1; ydelta<=1; ydelta++){
331 // for (Int_t zdelta=-1; zdelta<=1; zdelta++){
332 // for (Int_t xdelta=-2; xdelta<=0; xdelta++)
333 // for (Int_t ydelta=-2; ydelta<=0; ydelta++){
334 for (Int_t zdelta=zmin; zdelta<=zmax; zdelta++){
335 //TODO: what happens if x[2] is on the A-Side, but x[2]+zdelta*delta
336 // will be on the C-Side?
337 Float_t xyz[3]={x[0]+xdelta*delta, x[1]+ydelta*delta, x[2]+zdelta*delta};
339 GetCorrection(xyz,roc,dxyz);
340 Double_t adelta=zdelta*delta;
341 fitx.AddPoint(&adelta, dxyz[0]);
342 fity.AddPoint(&adelta, dxyz[1]);
343 fitz.AddPoint(&adelta, dxyz[2]);
349 dx[0] = fitx.GetParameter(1);
350 dx[1] = fity.GetParameter(1);
351 dx[2] = fitz.GetParameter(1);
354 void AliTPCCorrection::GetDistortionDz(const Float_t x[],const Short_t roc,Float_t dx[], Float_t delta) {
355 // author: marian.ivanov@cern.ch
357 // In this (virtual)function calculates the dx'/dz, dy'/dz and dz'/dz at given point (x,y,z)
358 // Generic implementation. Better precision can be acchieved knowing the internal structure
359 // of underlying trasnformation. Derived classes can reimplement it.
360 // To calculate distortion is fitted in small neighberhood:
361 // (x+-delta,y+-delta,z+-delta) where delta is an argument
364 // x[] - space point corrdinate
365 // roc - readout chamber identifier (important e.g to do not miss the side of detector)
366 // delta - define the size of neighberhood
368 // dx[] - array {dx'/dz, dy'/dz , dz'/dz }
370 static TLinearFitter fitx(2,"pol1");
371 static TLinearFitter fity(2,"pol1");
372 static TLinearFitter fitz(2,"pol1");
379 //adjust limits around CE to stay on one side
382 if ((x[2]+zmin*delta)<0){
388 if ((x[2]+zmax*delta)>0){
394 //TODO: in principle one shuld check that x[2]+zdelta*delta does not get 'out of' bounds,
395 // so close to the CE it doesn't change the sign, since then the corrections will be wrong ...
396 for (Int_t xdelta=-1; xdelta<=1; xdelta++)
397 for (Int_t ydelta=-1; ydelta<=1; ydelta++){
398 for (Int_t zdelta=zmin; zdelta<=zmax; zdelta++){
399 //TODO: what happens if x[2] is on the A-Side, but x[2]+zdelta*delta
400 // will be on the C-Side?
401 //TODO: For the C-Side, does this have the correct sign?
402 Float_t xyz[3]={x[0]+xdelta*delta, x[1]+ydelta*delta, x[2]+zdelta*delta};
404 GetDistortion(xyz,roc,dxyz);
405 Double_t adelta=zdelta*delta;
406 fitx.AddPoint(&adelta, dxyz[0]);
407 fity.AddPoint(&adelta, dxyz[1]);
408 fitz.AddPoint(&adelta, dxyz[2]);
414 dx[0] = fitx.GetParameter(1);
415 dx[1] = fity.GetParameter(1);
416 dx[2] = fitz.GetParameter(1);
419 void AliTPCCorrection::GetCorrectionIntegralDz(const Float_t x[],const Short_t roc,Float_t dx[], Float_t delta){
421 // Integrate 3D distortion along drift lines starting from the roc plane
422 // to the expected z position of the point, this assumes that dz is small
423 // and the error propagating to z' instead of the correct z is negligible
424 // To define the drift lines virtual function AliTPCCorrection::GetCorrectionDz is used
427 // x[] - space point corrdinate
428 // roc - readout chamber identifier (important e.g to do not miss the side of detector)
429 // delta - define the size of neighberhood
431 // dx[] - array { integral(dx'/dz), integral(dy'/dz) , integral(dz'/dz) }
433 Float_t zroc= ((roc%36)<18) ? fgkTPCZ0:-fgkTPCZ0;
434 Double_t zdrift = TMath::Abs(x[2]-zroc);
435 Int_t nsteps = Int_t(zdrift/delta)+1;
438 Float_t xyz[3]={x[0],x[1],zroc};
439 Float_t dxyz[3]={x[0],x[1],x[2]};
440 Short_t side=(roc/18)%2;
441 Float_t sign=1-2*side;
443 for (Int_t i=0;i<nsteps; i++){
444 //propagate backwards, therefore opposite signs
445 Float_t deltaZ=delta*(-sign);
446 // if (xyz[2]+deltaZ>fgkTPCZ0) deltaZ=TMath::Abs(xyz[2]-fgkTPCZ0);
447 // if (xyz[2]-deltaZ<-fgkTPCZ0) deltaZ=TMath::Abs(xyz[2]-fgkTPCZ0);
448 // protect again integrating through the CE
450 if (xyz[2]+deltaZ<0) deltaZ=-xyz[2]+1e-20;
452 if (xyz[2]+deltaZ>0) deltaZ=xyz[2]-+1e-20;
454 // since at larger drift (smaller z) the corrections are larger (absolute, but negative)
455 // the slopes will be positive.
456 // but since we chose deltaZ opposite sign the singn of the corretion should be fine
458 GetCorrectionDz(xyz,roc,dxyz,delta);
459 xyz[0]+=deltaZ*dxyz[0];
460 xyz[1]+=deltaZ*dxyz[1];
462 sumdz+=deltaZ*dxyz[2];
467 dx[2]= sumdz; //TODO: is sumdz correct?
470 void AliTPCCorrection::GetDistortionIntegralDz(const Float_t x[],const Short_t roc,Float_t dx[], Float_t delta){
472 // Integrate 3D distortion along drift lines
473 // To define the drift lines virtual function AliTPCCorrection::GetCorrectionDz is used
476 // x[] - space point corrdinate
477 // roc - readout chamber identifier (important e.g to do not miss the side of detector)
478 // delta - define the size of neighberhood
480 // dx[] - array { integral(dx'/dz), integral(dy'/dz) , integral(dz'/dz) }
482 Float_t zroc= ((roc%36)<18) ? fgkTPCZ0:-fgkTPCZ0;
483 Double_t zdrift = TMath::Abs(x[2]-zroc);
484 Int_t nsteps = Int_t(zdrift/delta)+1;
487 Float_t xyz[3]={x[0],x[1],x[2]};
488 Float_t dxyz[3]={x[0],x[1],x[2]};
489 Float_t sign=((roc%36)<18) ? 1.:-1.;
491 for (Int_t i=0;i<nsteps; i++){
492 Float_t deltaZ=delta;
493 if (xyz[2]+deltaZ>fgkTPCZ0) deltaZ=TMath::Abs(xyz[2]-zroc);
494 if (xyz[2]-deltaZ<-fgkTPCZ0) deltaZ=TMath::Abs(xyz[2]-zroc);
495 // since at larger drift (smaller z) the distortions are larger
496 // the slopes will be negative.
497 // and since we are moving towards the read-out plane the deltaZ for
498 // weighting the dK/dz should have the opposite sign
500 GetDistortionDz(xyz,roc,dxyz,delta);
501 xyz[0]+=-deltaZ*dxyz[0];
502 xyz[1]+=-deltaZ*dxyz[1];
503 xyz[2]+=deltaZ; //TODO: Should this also be corrected for the dxyz[2]
504 sumdz+=-deltaZ*dxyz[2];
509 dx[2]= sumdz; //TODO: is sumdz correct?
514 void AliTPCCorrection::Init() {
516 // Initialization funtion (not used at the moment)
520 void AliTPCCorrection::Update(const TTimeStamp &/*timeStamp*/) {
526 void AliTPCCorrection::Print(Option_t* /*option*/) const {
528 // Print function to check which correction classes are used
529 // option=="d" prints details regarding the setted magnitude
530 // option=="a" prints the C0 and C1 coefficents for calibration purposes
532 printf("TPC spacepoint correction: \"%s\"\n",GetTitle());
535 void AliTPCCorrection:: SetOmegaTauT1T2(Float_t /*omegaTau*/,Float_t t1,Float_t t2) {
537 // Virtual funtion to pass the wt values (might become event dependent) to the inherited classes
538 // t1 and t2 represent the "effective omegaTau" corrections and were measured in a dedicated
543 //SetOmegaTauT1T2(omegaTau, t1, t2);
546 TH2F* AliTPCCorrection::CreateHistoDRinXY(Float_t z,Int_t nx,Int_t ny) {
548 // Simple plot functionality.
549 // Returns a 2d hisogram which represents the corrections in radial direction (dr)
550 // in respect to position z within the XY plane.
551 // The histogramm has nx times ny entries.
553 AliTPCParam* tpcparam = new AliTPCParamSR;
555 TH2F *h=CreateTH2F("dr_xy",GetTitle(),"x [cm]","y [cm]","dr [cm]",
556 nx,-250.,250.,ny,-250.,250.);
559 Int_t roc=z>0.?0:18; // FIXME
560 for (Int_t iy=1;iy<=ny;++iy) {
561 x[1]=h->GetYaxis()->GetBinCenter(iy);
562 for (Int_t ix=1;ix<=nx;++ix) {
563 x[0]=h->GetXaxis()->GetBinCenter(ix);
564 GetCorrection(x,roc,dx);
565 Float_t r0=TMath::Sqrt((x[0] )*(x[0] )+(x[1] )*(x[1] ));
566 if (tpcparam->GetPadRowRadii(0,0)<=r0 && r0<=tpcparam->GetPadRowRadii(36,95)) {
567 Float_t r1=TMath::Sqrt((x[0]+dx[0])*(x[0]+dx[0])+(x[1]+dx[1])*(x[1]+dx[1]));
568 h->SetBinContent(ix,iy,r1-r0);
571 h->SetBinContent(ix,iy,0.);
578 TH2F* AliTPCCorrection::CreateHistoDRPhiinXY(Float_t z,Int_t nx,Int_t ny) {
580 // Simple plot functionality.
581 // Returns a 2d hisogram which represents the corrections in rphi direction (drphi)
582 // in respect to position z within the XY plane.
583 // The histogramm has nx times ny entries.
586 AliTPCParam* tpcparam = new AliTPCParamSR;
588 TH2F *h=CreateTH2F("drphi_xy",GetTitle(),"x [cm]","y [cm]","drphi [cm]",
589 nx,-250.,250.,ny,-250.,250.);
592 Int_t roc=z>0.?0:18; // FIXME
593 for (Int_t iy=1;iy<=ny;++iy) {
594 x[1]=h->GetYaxis()->GetBinCenter(iy);
595 for (Int_t ix=1;ix<=nx;++ix) {
596 x[0]=h->GetXaxis()->GetBinCenter(ix);
597 GetCorrection(x,roc,dx);
598 Float_t r0=TMath::Sqrt((x[0] )*(x[0] )+(x[1] )*(x[1] ));
599 if (tpcparam->GetPadRowRadii(0,0)<=r0 && r0<=tpcparam->GetPadRowRadii(36,95)) {
600 Float_t phi0=TMath::ATan2(x[1] ,x[0] );
601 Float_t phi1=TMath::ATan2(x[1]+dx[1],x[0]+dx[0]);
603 Float_t dphi=phi1-phi0;
604 if (dphi<TMath::Pi()) dphi+=TMath::TwoPi();
605 if (dphi>TMath::Pi()) dphi-=TMath::TwoPi();
607 h->SetBinContent(ix,iy,r0*dphi);
610 h->SetBinContent(ix,iy,0.);
617 TH2F* AliTPCCorrection::CreateHistoDZinXY(Float_t z,Int_t nx,Int_t ny) {
619 // Simple plot functionality.
620 // Returns a 2d hisogram which represents the corrections in longitudinal direction (dz)
621 // in respect to position z within the XY plane.
622 // The histogramm has nx times ny entries.
625 AliTPCParam* tpcparam = new AliTPCParamSR;
627 TH2F *h=CreateTH2F("dz_xy",GetTitle(),"x [cm]","y [cm]","dz [cm]",
628 nx,-250.,250.,ny,-250.,250.);
631 Int_t roc=z>0.?0:18; // FIXME
632 for (Int_t iy=1;iy<=ny;++iy) {
633 x[1]=h->GetYaxis()->GetBinCenter(iy);
634 for (Int_t ix=1;ix<=nx;++ix) {
635 x[0]=h->GetXaxis()->GetBinCenter(ix);
636 GetCorrection(x,roc,dx);
637 Float_t r0=TMath::Sqrt((x[0] )*(x[0] )+(x[1] )*(x[1] ));
638 if (tpcparam->GetPadRowRadii(0,0)<=r0 && r0<=tpcparam->GetPadRowRadii(36,95)) {
639 h->SetBinContent(ix,iy,dx[2]);
642 h->SetBinContent(ix,iy,0.);
649 TH2F* AliTPCCorrection::CreateHistoDRinZR(Float_t phi,Int_t nz,Int_t nr) {
651 // Simple plot functionality.
652 // Returns a 2d hisogram which represents the corrections in r direction (dr)
653 // in respect to angle phi within the ZR plane.
654 // The histogramm has nx times ny entries.
656 TH2F *h=CreateTH2F("dr_zr",GetTitle(),"z [cm]","r [cm]","dr [cm]",
657 nz,-250.,250.,nr,85.,250.);
659 for (Int_t ir=1;ir<=nr;++ir) {
660 Float_t radius=h->GetYaxis()->GetBinCenter(ir);
661 x[0]=radius*TMath::Cos(phi);
662 x[1]=radius*TMath::Sin(phi);
663 for (Int_t iz=1;iz<=nz;++iz) {
664 x[2]=h->GetXaxis()->GetBinCenter(iz);
665 Int_t roc=x[2]>0.?0:18; // FIXME
666 GetCorrection(x,roc,dx);
667 Float_t r0=TMath::Sqrt((x[0] )*(x[0] )+(x[1] )*(x[1] ));
668 Float_t r1=TMath::Sqrt((x[0]+dx[0])*(x[0]+dx[0])+(x[1]+dx[1])*(x[1]+dx[1]));
669 h->SetBinContent(iz,ir,r1-r0);
676 TH2F* AliTPCCorrection::CreateHistoDRPhiinZR(Float_t phi,Int_t nz,Int_t nr) {
678 // Simple plot functionality.
679 // Returns a 2d hisogram which represents the corrections in rphi direction (drphi)
680 // in respect to angle phi within the ZR plane.
681 // The histogramm has nx times ny entries.
683 TH2F *h=CreateTH2F("drphi_zr",GetTitle(),"z [cm]","r [cm]","drphi [cm]",
684 nz,-250.,250.,nr,85.,250.);
686 for (Int_t iz=1;iz<=nz;++iz) {
687 x[2]=h->GetXaxis()->GetBinCenter(iz);
688 Int_t roc=x[2]>0.?0:18; // FIXME
689 for (Int_t ir=1;ir<=nr;++ir) {
690 Float_t radius=h->GetYaxis()->GetBinCenter(ir);
691 x[0]=radius*TMath::Cos(phi);
692 x[1]=radius*TMath::Sin(phi);
693 GetCorrection(x,roc,dx);
694 Float_t r0=TMath::Sqrt((x[0] )*(x[0] )+(x[1] )*(x[1] ));
695 Float_t phi0=TMath::ATan2(x[1] ,x[0] );
696 Float_t phi1=TMath::ATan2(x[1]+dx[1],x[0]+dx[0]);
698 Float_t dphi=phi1-phi0;
699 if (dphi<TMath::Pi()) dphi+=TMath::TwoPi();
700 if (dphi>TMath::Pi()) dphi-=TMath::TwoPi();
702 h->SetBinContent(iz,ir,r0*dphi);
708 TH2F* AliTPCCorrection::CreateHistoDZinZR(Float_t phi,Int_t nz,Int_t nr) {
710 // Simple plot functionality.
711 // Returns a 2d hisogram which represents the corrections in longitudinal direction (dz)
712 // in respect to angle phi within the ZR plane.
713 // The histogramm has nx times ny entries.
715 TH2F *h=CreateTH2F("dz_zr",GetTitle(),"z [cm]","r [cm]","dz [cm]",
716 nz,-250.,250.,nr,85.,250.);
718 for (Int_t ir=1;ir<=nr;++ir) {
719 Float_t radius=h->GetYaxis()->GetBinCenter(ir);
720 x[0]=radius*TMath::Cos(phi);
721 x[1]=radius*TMath::Sin(phi);
722 for (Int_t iz=1;iz<=nz;++iz) {
723 x[2]=h->GetXaxis()->GetBinCenter(iz);
724 Int_t roc=x[2]>0.?0:18; // FIXME
725 GetCorrection(x,roc,dx);
726 h->SetBinContent(iz,ir,dx[2]);
734 TH2F* AliTPCCorrection::CreateTH2F(const char *name,const char *title,
735 const char *xlabel,const char *ylabel,const char *zlabel,
736 Int_t nbinsx,Double_t xlow,Double_t xup,
737 Int_t nbinsy,Double_t ylow,Double_t yup) {
739 // Helper function to create a 2d histogramm of given size
745 while (gDirectory->FindObject(hname.Data())) {
752 TH2F *h=new TH2F(hname.Data(),title,
755 h->GetXaxis()->SetTitle(xlabel);
756 h->GetYaxis()->SetTitle(ylabel);
757 h->GetZaxis()->SetTitle(zlabel);
762 // Simple Interpolation functions: e.g. with bi(tri)cubic interpolations (not yet in TH2 and TH3)
764 void AliTPCCorrection::Interpolate2DEdistortion( const Int_t order, const Double_t r, const Double_t z,
765 const Double_t er[kNZ][kNR], Double_t &erValue ) {
767 // Interpolate table - 2D interpolation
769 Double_t saveEr[5] = {0,0,0,0,0};
771 Search( kNZ, fgkZList, z, fJLow ) ;
772 Search( kNR, fgkRList, r, fKLow ) ;
773 if ( fJLow < 0 ) fJLow = 0 ; // check if out of range
774 if ( fKLow < 0 ) fKLow = 0 ;
775 if ( fJLow + order >= kNZ - 1 ) fJLow = kNZ - 1 - order ;
776 if ( fKLow + order >= kNR - 1 ) fKLow = kNR - 1 - order ;
778 for ( Int_t j = fJLow ; j < fJLow + order + 1 ; j++ ) {
779 saveEr[j-fJLow] = Interpolate( &fgkRList[fKLow], &er[j][fKLow], order, r ) ;
781 erValue = Interpolate( &fgkZList[fJLow], saveEr, order, z ) ;
785 void AliTPCCorrection::Interpolate3DEdistortion( const Int_t order, const Double_t r, const Float_t phi, const Double_t z,
786 const Double_t er[kNZ][kNPhi][kNR], const Double_t ephi[kNZ][kNPhi][kNR], const Double_t ez[kNZ][kNPhi][kNR],
787 Double_t &erValue, Double_t &ephiValue, Double_t &ezValue) {
789 // Interpolate table - 3D interpolation
792 Double_t saveEr[5]= {0,0,0,0,0};
793 Double_t savedEr[5]= {0,0,0,0,0} ;
795 Double_t saveEphi[5]= {0,0,0,0,0};
796 Double_t savedEphi[5]= {0,0,0,0,0} ;
798 Double_t saveEz[5]= {0,0,0,0,0};
799 Double_t savedEz[5]= {0,0,0,0,0} ;
801 Search( kNZ, fgkZList, z, fILow ) ;
802 Search( kNPhi, fgkPhiList, z, fJLow ) ;
803 Search( kNR, fgkRList, r, fKLow ) ;
805 if ( fILow < 0 ) fILow = 0 ; // check if out of range
806 if ( fJLow < 0 ) fJLow = 0 ;
807 if ( fKLow < 0 ) fKLow = 0 ;
809 if ( fILow + order >= kNZ - 1 ) fILow = kNZ - 1 - order ;
810 if ( fJLow + order >= kNPhi - 1 ) fJLow = kNPhi - 1 - order ;
811 if ( fKLow + order >= kNR - 1 ) fKLow = kNR - 1 - order ;
813 for ( Int_t i = fILow ; i < fILow + order + 1 ; i++ ) {
814 for ( Int_t j = fJLow ; j < fJLow + order + 1 ; j++ ) {
815 saveEr[j-fJLow] = Interpolate( &fgkRList[fKLow], &er[i][j][fKLow], order, r ) ;
816 saveEphi[j-fJLow] = Interpolate( &fgkRList[fKLow], &ephi[i][j][fKLow], order, r ) ;
817 saveEz[j-fJLow] = Interpolate( &fgkRList[fKLow], &ez[i][j][fKLow], order, r ) ;
819 savedEr[i-fILow] = Interpolate( &fgkPhiList[fJLow], saveEr, order, phi ) ;
820 savedEphi[i-fILow] = Interpolate( &fgkPhiList[fJLow], saveEphi, order, phi ) ;
821 savedEz[i-fILow] = Interpolate( &fgkPhiList[fJLow], saveEz, order, phi ) ;
823 erValue = Interpolate( &fgkZList[fILow], savedEr, order, z ) ;
824 ephiValue = Interpolate( &fgkZList[fILow], savedEphi, order, z ) ;
825 ezValue = Interpolate( &fgkZList[fILow], savedEz, order, z ) ;
829 Double_t AliTPCCorrection::Interpolate2DTable( const Int_t order, const Double_t x, const Double_t y,
830 const Int_t nx, const Int_t ny, const Double_t xv[], const Double_t yv[],
831 const TMatrixD &array ) {
833 // Interpolate table (TMatrix format) - 2D interpolation
836 static Int_t jlow = 0, klow = 0 ;
837 Double_t saveArray[5] = {0,0,0,0,0} ;
839 Search( nx, xv, x, jlow ) ;
840 Search( ny, yv, y, klow ) ;
841 if ( jlow < 0 ) jlow = 0 ; // check if out of range
842 if ( klow < 0 ) klow = 0 ;
843 if ( jlow + order >= nx - 1 ) jlow = nx - 1 - order ;
844 if ( klow + order >= ny - 1 ) klow = ny - 1 - order ;
846 for ( Int_t j = jlow ; j < jlow + order + 1 ; j++ )
848 Double_t *ajkl = &((TMatrixD&)array)(j,klow);
849 saveArray[j-jlow] = Interpolate( &yv[klow], ajkl , order, y ) ;
852 return( Interpolate( &xv[jlow], saveArray, order, x ) ) ;
856 Double_t AliTPCCorrection::Interpolate3DTable( const Int_t order, const Double_t x, const Double_t y, const Double_t z,
857 const Int_t nx, const Int_t ny, const Int_t nz,
858 const Double_t xv[], const Double_t yv[], const Double_t zv[],
859 TMatrixD **arrayofArrays ) {
861 // Interpolate table (TMatrix format) - 3D interpolation
864 static Int_t ilow = 0, jlow = 0, klow = 0 ;
865 Double_t saveArray[5]= {0,0,0,0,0};
866 Double_t savedArray[5]= {0,0,0,0,0} ;
868 Search( nx, xv, x, ilow ) ;
869 Search( ny, yv, y, jlow ) ;
870 Search( nz, zv, z, klow ) ;
872 if ( ilow < 0 ) ilow = 0 ; // check if out of range
873 if ( jlow < 0 ) jlow = 0 ;
874 if ( klow < 0 ) klow = 0 ;
876 if ( ilow + order >= nx - 1 ) ilow = nx - 1 - order ;
877 if ( jlow + order >= ny - 1 ) jlow = ny - 1 - order ;
878 if ( klow + order >= nz - 1 ) klow = nz - 1 - order ;
880 for ( Int_t k = klow ; k < klow + order + 1 ; k++ )
882 TMatrixD &table = *arrayofArrays[k] ;
883 for ( Int_t i = ilow ; i < ilow + order + 1 ; i++ )
885 saveArray[i-ilow] = Interpolate( &yv[jlow], &table(i,jlow), order, y ) ;
887 savedArray[k-klow] = Interpolate( &xv[ilow], saveArray, order, x ) ;
889 return( Interpolate( &zv[klow], savedArray, order, z ) ) ;
893 Double_t AliTPCCorrection::Interpolate( const Double_t xArray[], const Double_t yArray[],
894 const Int_t order, const Double_t x ) {
896 // Interpolate function Y(x) using linear (order=1) or quadratic (order=2) interpolation.
900 if ( order == 2 ) { // Quadratic Interpolation = 2
901 y = (x-xArray[1]) * (x-xArray[2]) * yArray[0] / ( (xArray[0]-xArray[1]) * (xArray[0]-xArray[2]) ) ;
902 y += (x-xArray[2]) * (x-xArray[0]) * yArray[1] / ( (xArray[1]-xArray[2]) * (xArray[1]-xArray[0]) ) ;
903 y += (x-xArray[0]) * (x-xArray[1]) * yArray[2] / ( (xArray[2]-xArray[0]) * (xArray[2]-xArray[1]) ) ;
904 } else { // Linear Interpolation = 1
905 y = yArray[0] + ( yArray[1]-yArray[0] ) * ( x-xArray[0] ) / ( xArray[1] - xArray[0] ) ;
912 Float_t AliTPCCorrection::Interpolate2DTable( const Int_t order, const Double_t x, const Double_t y,
913 const Int_t nx, const Int_t ny, const Double_t xv[], const Double_t yv[],
914 const TMatrixF &array ) {
916 // Interpolate table (TMatrix format) - 2D interpolation
917 // Float version (in order to decrease the OCDB size)
920 static Int_t jlow = 0, klow = 0 ;
921 Float_t saveArray[5] = {0.,0.,0.,0.,0.} ;
923 Search( nx, xv, x, jlow ) ;
924 Search( ny, yv, y, klow ) ;
925 if ( jlow < 0 ) jlow = 0 ; // check if out of range
926 if ( klow < 0 ) klow = 0 ;
927 if ( jlow + order >= nx - 1 ) jlow = nx - 1 - order ;
928 if ( klow + order >= ny - 1 ) klow = ny - 1 - order ;
930 for ( Int_t j = jlow ; j < jlow + order + 1 ; j++ )
932 Float_t *ajkl = &((TMatrixF&)array)(j,klow);
933 saveArray[j-jlow] = Interpolate( &yv[klow], ajkl , order, y ) ;
936 return( Interpolate( &xv[jlow], saveArray, order, x ) ) ;
940 Float_t AliTPCCorrection::Interpolate3DTable( const Int_t order, const Double_t x, const Double_t y, const Double_t z,
941 const Int_t nx, const Int_t ny, const Int_t nz,
942 const Double_t xv[], const Double_t yv[], const Double_t zv[],
943 TMatrixF **arrayofArrays ) {
945 // Interpolate table (TMatrix format) - 3D interpolation
946 // Float version (in order to decrease the OCDB size)
949 static Int_t ilow = 0, jlow = 0, klow = 0 ;
950 Float_t saveArray[5]= {0.,0.,0.,0.,0.};
951 Float_t savedArray[5]= {0.,0.,0.,0.,0.} ;
953 Search( nx, xv, x, ilow ) ;
954 Search( ny, yv, y, jlow ) ;
955 Search( nz, zv, z, klow ) ;
957 if ( ilow < 0 ) ilow = 0 ; // check if out of range
958 if ( jlow < 0 ) jlow = 0 ;
959 if ( klow < 0 ) klow = 0 ;
961 if ( ilow + order >= nx - 1 ) ilow = nx - 1 - order ;
962 if ( jlow + order >= ny - 1 ) jlow = ny - 1 - order ;
963 if ( klow + order >= nz - 1 ) klow = nz - 1 - order ;
965 for ( Int_t k = klow ; k < klow + order + 1 ; k++ )
967 TMatrixF &table = *arrayofArrays[k] ;
968 for ( Int_t i = ilow ; i < ilow + order + 1 ; i++ )
970 saveArray[i-ilow] = Interpolate( &yv[jlow], &table(i,jlow), order, y ) ;
972 savedArray[k-klow] = Interpolate( &xv[ilow], saveArray, order, x ) ;
974 return( Interpolate( &zv[klow], savedArray, order, z ) ) ;
977 Float_t AliTPCCorrection::Interpolate( const Double_t xArray[], const Float_t yArray[],
978 const Int_t order, const Double_t x ) {
980 // Interpolate function Y(x) using linear (order=1) or quadratic (order=2) interpolation.
981 // Float version (in order to decrease the OCDB size)
985 if ( order == 2 ) { // Quadratic Interpolation = 2
986 y = (x-xArray[1]) * (x-xArray[2]) * yArray[0] / ( (xArray[0]-xArray[1]) * (xArray[0]-xArray[2]) ) ;
987 y += (x-xArray[2]) * (x-xArray[0]) * yArray[1] / ( (xArray[1]-xArray[2]) * (xArray[1]-xArray[0]) ) ;
988 y += (x-xArray[0]) * (x-xArray[1]) * yArray[2] / ( (xArray[2]-xArray[0]) * (xArray[2]-xArray[1]) ) ;
989 } else { // Linear Interpolation = 1
990 y = yArray[0] + ( yArray[1]-yArray[0] ) * ( x-xArray[0] ) / ( xArray[1] - xArray[0] ) ;
999 void AliTPCCorrection::Search( const Int_t n, const Double_t xArray[], const Double_t x, Int_t &low ) {
1001 // Search an ordered table by starting at the most recently used point
1004 Long_t middle, high ;
1005 Int_t ascend = 0, increment = 1 ;
1007 if ( xArray[n-1] >= xArray[0] ) ascend = 1 ; // Ascending ordered table if true
1009 if ( low < 0 || low > n-1 ) {
1010 low = -1 ; high = n ;
1011 } else { // Ordered Search phase
1012 if ( (Int_t)( x >= xArray[low] ) == ascend ) {
1013 if ( low == n-1 ) return ;
1015 while ( (Int_t)( x >= xArray[high] ) == ascend ) {
1018 high = low + increment ;
1019 if ( high > n-1 ) { high = n ; break ; }
1022 if ( low == 0 ) { low = -1 ; return ; }
1024 while ( (Int_t)( x < xArray[low] ) == ascend ) {
1027 if ( increment >= high ) { low = -1 ; break ; }
1028 else low = high - increment ;
1033 while ( (high-low) != 1 ) { // Binary Search Phase
1034 middle = ( high + low ) / 2 ;
1035 if ( (Int_t)( x >= xArray[middle] ) == ascend )
1041 if ( x == xArray[n-1] ) low = n-2 ;
1042 if ( x == xArray[0] ) low = 0 ;
1046 void AliTPCCorrection::InitLookUpfulcrums() {
1048 // Initialization of interpolation points - for main look up table
1049 // (course grid in the middle, fine grid on the borders)
1052 AliTPCROC * roc = AliTPCROC::Instance();
1053 const Double_t rLow = TMath::Floor(roc->GetPadRowRadii(0,0))-1; // first padRow plus some margin
1057 for (Int_t i = 1; i<kNR; i++) {
1058 fgkRList[i] = fgkRList[i-1] + 3.5; // 3.5 cm spacing
1059 if (fgkRList[i]<90 ||fgkRList[i]>245)
1060 fgkRList[i] = fgkRList[i-1] + 0.5; // 0.5 cm spacing
1061 else if (fgkRList[i]<100 || fgkRList[i]>235)
1062 fgkRList[i] = fgkRList[i-1] + 1.5; // 1.5 cm spacing
1063 else if (fgkRList[i]<120 || fgkRList[i]>225)
1064 fgkRList[i] = fgkRList[i-1] + 2.5; // 2.5 cm spacing
1068 fgkZList[0] = -249.5;
1069 fgkZList[kNZ-1] = 249.5;
1070 for (Int_t j = 1; j<kNZ/2; j++) {
1071 fgkZList[j] = fgkZList[j-1];
1072 if (TMath::Abs(fgkZList[j])< 0.15)
1073 fgkZList[j] = fgkZList[j-1] + 0.09; // 0.09 cm spacing
1074 else if(TMath::Abs(fgkZList[j])< 0.6)
1075 fgkZList[j] = fgkZList[j-1] + 0.4; // 0.4 cm spacing
1076 else if (TMath::Abs(fgkZList[j])< 2.5 || TMath::Abs(fgkZList[j])>248)
1077 fgkZList[j] = fgkZList[j-1] + 0.5; // 0.5 cm spacing
1078 else if (TMath::Abs(fgkZList[j])<10 || TMath::Abs(fgkZList[j])>235)
1079 fgkZList[j] = fgkZList[j-1] + 1.5; // 1.5 cm spacing
1080 else if (TMath::Abs(fgkZList[j])<25 || TMath::Abs(fgkZList[j])>225)
1081 fgkZList[j] = fgkZList[j-1] + 2.5; // 2.5 cm spacing
1083 fgkZList[j] = fgkZList[j-1] + 4; // 4 cm spacing
1085 fgkZList[kNZ-j-1] = -fgkZList[j];
1089 for (Int_t k = 0; k<kNPhi; k++)
1090 fgkPhiList[k] = TMath::TwoPi()*k/(kNPhi-1);
1096 void AliTPCCorrection::PoissonRelaxation2D(TMatrixD &arrayV, TMatrixD &chargeDensity,
1097 TMatrixD &arrayErOverEz, TMatrixD &arrayDeltaEz,
1098 const Int_t rows, const Int_t columns, const Int_t iterations,
1099 const Bool_t rocDisplacement ) {
1101 // Solve Poisson's Equation by Relaxation Technique in 2D (assuming cylindrical symmetry)
1103 // Solve Poissons equation in a cylindrical coordinate system. The arrayV matrix must be filled with the
1104 // boundary conditions on the first and last rows, and the first and last columns. The remainder of the
1105 // array can be blank or contain a preliminary guess at the solution. The Charge density matrix contains
1106 // the enclosed spacecharge density at each point. The charge density matrix can be full of zero's if
1107 // you wish to solve Laplaces equation however it should not contain random numbers or you will get
1108 // random numbers back as a solution.
1109 // Poisson's equation is solved by iteratively relaxing the matrix to the final solution. In order to
1110 // speed up the convergence to the best solution, this algorithm does a binary expansion of the solution
1111 // space. First it solves the problem on a very sparse grid by skipping rows and columns in the original
1112 // matrix. Then it doubles the number of points and solves the problem again. Then it doubles the
1113 // number of points and solves the problem again. This happens several times until the maximum number
1114 // of points has been included in the array.
1116 // NOTE: In order for this algorithmto work, the number of rows and columns must be a power of 2 plus one.
1117 // So rows == 2**M + 1 and columns == 2**N + 1. The number of rows and columns can be different.
1119 // NOTE: rocDisplacement is used to include (or ignore) the ROC misalignment in the dz calculation
1121 // Original code by Jim Thomas (STAR TPC Collaboration)
1124 Double_t ezField = (fgkCathodeV-fgkGG)/fgkTPCZ0; // = ALICE Electric Field (V/cm) Magnitude ~ -400 V/cm;
1126 const Float_t gridSizeR = (fgkOFCRadius-fgkIFCRadius) / (rows-1) ;
1127 const Float_t gridSizeZ = fgkTPCZ0 / (columns-1) ;
1128 const Float_t ratio = gridSizeR*gridSizeR / (gridSizeZ*gridSizeZ) ;
1130 TMatrixD arrayEr(rows,columns) ;
1131 TMatrixD arrayEz(rows,columns) ;
1133 //Check that number of rows and columns is suitable for a binary expansion
1135 if ( !IsPowerOfTwo(rows-1) ) {
1136 AliError("PoissonRelaxation - Error in the number of rows. Must be 2**M - 1");
1139 if ( !IsPowerOfTwo(columns-1) ) {
1140 AliError("PoissonRelaxation - Error in the number of columns. Must be 2**N - 1");
1144 // Solve Poisson's equation in cylindrical coordinates by relaxation technique
1145 // Allow for different size grid spacing in R and Z directions
1146 // Use a binary expansion of the size of the matrix to speed up the solution of the problem
1148 Int_t iOne = (rows-1)/4 ;
1149 Int_t jOne = (columns-1)/4 ;
1150 // Solve for N in 2**N, add one.
1151 Int_t loops = 1 + (int) ( 0.5 + TMath::Log2( (double) TMath::Max(iOne,jOne) ) ) ;
1153 for ( Int_t count = 0 ; count < loops ; count++ ) {
1154 // Loop while the matrix expands & the resolution increases.
1156 Float_t tempGridSizeR = gridSizeR * iOne ;
1157 Float_t tempRatio = ratio * iOne * iOne / ( jOne * jOne ) ;
1158 Float_t tempFourth = 1.0 / (2.0 + 2.0*tempRatio) ;
1160 // Do this the standard C++ way to avoid gcc extensions for Float_t coef1[rows]
1161 std::vector<float> coef1(rows) ;
1162 std::vector<float> coef2(rows) ;
1164 for ( Int_t i = iOne ; i < rows-1 ; i+=iOne ) {
1165 Float_t radius = fgkIFCRadius + i*gridSizeR ;
1166 coef1[i] = 1.0 + tempGridSizeR/(2*radius);
1167 coef2[i] = 1.0 - tempGridSizeR/(2*radius);
1170 TMatrixD sumChargeDensity(rows,columns) ;
1172 for ( Int_t i = iOne ; i < rows-1 ; i += iOne ) {
1173 Float_t radius = fgkIFCRadius + iOne*gridSizeR ;
1174 for ( Int_t j = jOne ; j < columns-1 ; j += jOne ) {
1175 if ( iOne == 1 && jOne == 1 ) sumChargeDensity(i,j) = chargeDensity(i,j) ;
1177 // Add up all enclosed charge density contributions within 1/2 unit in all directions
1178 Float_t weight = 0.0 ;
1180 sumChargeDensity(i,j) = 0.0 ;
1181 for ( Int_t ii = i-iOne/2 ; ii <= i+iOne/2 ; ii++ ) {
1182 for ( Int_t jj = j-jOne/2 ; jj <= j+jOne/2 ; jj++ ) {
1183 if ( ii == i-iOne/2 || ii == i+iOne/2 || jj == j-jOne/2 || jj == j+jOne/2 ) weight = 0.5 ;
1186 // Note that this is cylindrical geometry
1187 sumChargeDensity(i,j) += chargeDensity(ii,jj)*weight*radius ;
1188 sum += weight*radius ;
1191 sumChargeDensity(i,j) /= sum ;
1193 sumChargeDensity(i,j) *= tempGridSizeR*tempGridSizeR; // just saving a step later on
1197 for ( Int_t k = 1 ; k <= iterations; k++ ) {
1198 // Solve Poisson's Equation
1199 // Over-relaxation index, must be >= 1 but < 2. Arrange for it to evolve from 2 => 1
1200 // as interations increase.
1201 Float_t overRelax = 1.0 + TMath::Sqrt( TMath::Cos( (k*TMath::PiOver2())/iterations ) ) ;
1202 Float_t overRelaxM1 = overRelax - 1.0 ;
1203 Float_t overRelaxtempFourth, overRelaxcoef5 ;
1204 overRelaxtempFourth = overRelax * tempFourth ;
1205 overRelaxcoef5 = overRelaxM1 / overRelaxtempFourth ;
1207 for ( Int_t i = iOne ; i < rows-1 ; i += iOne ) {
1208 for ( Int_t j = jOne ; j < columns-1 ; j += jOne ) {
1210 arrayV(i,j) = ( coef2[i] * arrayV(i-iOne,j)
1211 + tempRatio * ( arrayV(i,j-jOne) + arrayV(i,j+jOne) )
1212 - overRelaxcoef5 * arrayV(i,j)
1213 + coef1[i] * arrayV(i+iOne,j)
1214 + sumChargeDensity(i,j)
1215 ) * overRelaxtempFourth;
1219 if ( k == iterations ) {
1220 // After full solution is achieved, copy low resolution solution into higher res array
1221 for ( Int_t i = iOne ; i < rows-1 ; i += iOne ) {
1222 for ( Int_t j = jOne ; j < columns-1 ; j += jOne ) {
1225 arrayV(i+iOne/2,j) = ( arrayV(i+iOne,j) + arrayV(i,j) ) / 2 ;
1226 if ( i == iOne ) arrayV(i-iOne/2,j) = ( arrayV(0,j) + arrayV(iOne,j) ) / 2 ;
1229 arrayV(i,j+jOne/2) = ( arrayV(i,j+jOne) + arrayV(i,j) ) / 2 ;
1230 if ( j == jOne ) arrayV(i,j-jOne/2) = ( arrayV(i,0) + arrayV(i,jOne) ) / 2 ;
1232 if ( iOne > 1 && jOne > 1 ) {
1233 arrayV(i+iOne/2,j+jOne/2) = ( arrayV(i+iOne,j+jOne) + arrayV(i,j) ) / 2 ;
1234 if ( i == iOne ) arrayV(i-iOne/2,j-jOne/2) = ( arrayV(0,j-jOne) + arrayV(iOne,j) ) / 2 ;
1235 if ( j == jOne ) arrayV(i-iOne/2,j-jOne/2) = ( arrayV(i-iOne,0) + arrayV(i,jOne) ) / 2 ;
1236 // Note that this leaves a point at the upper left and lower right corners uninitialized.
1237 // -> Not a big deal.
1246 iOne = iOne / 2 ; if ( iOne < 1 ) iOne = 1 ;
1247 jOne = jOne / 2 ; if ( jOne < 1 ) jOne = 1 ;
1249 sumChargeDensity.Clear();
1252 // Differentiate V(r) and solve for E(r) using special equations for the first and last rows
1253 for ( Int_t j = 0 ; j < columns ; j++ ) {
1254 for ( Int_t i = 1 ; i < rows-1 ; i++ ) arrayEr(i,j) = -1 * ( arrayV(i+1,j) - arrayV(i-1,j) ) / (2*gridSizeR) ;
1255 arrayEr(0,j) = -1 * ( -0.5*arrayV(2,j) + 2.0*arrayV(1,j) - 1.5*arrayV(0,j) ) / gridSizeR ;
1256 arrayEr(rows-1,j) = -1 * ( 1.5*arrayV(rows-1,j) - 2.0*arrayV(rows-2,j) + 0.5*arrayV(rows-3,j) ) / gridSizeR ;
1259 // Differentiate V(z) and solve for E(z) using special equations for the first and last columns
1260 for ( Int_t i = 0 ; i < rows ; i++) {
1261 for ( Int_t j = 1 ; j < columns-1 ; j++ ) arrayEz(i,j) = -1 * ( arrayV(i,j+1) - arrayV(i,j-1) ) / (2*gridSizeZ) ;
1262 arrayEz(i,0) = -1 * ( -0.5*arrayV(i,2) + 2.0*arrayV(i,1) - 1.5*arrayV(i,0) ) / gridSizeZ ;
1263 arrayEz(i,columns-1) = -1 * ( 1.5*arrayV(i,columns-1) - 2.0*arrayV(i,columns-2) + 0.5*arrayV(i,columns-3) ) / gridSizeZ ;
1266 for ( Int_t i = 0 ; i < rows ; i++) {
1267 // Note: go back and compare to old version of this code. See notes below.
1268 // JT Test ... attempt to divide by real Ez not Ez to first order
1269 for ( Int_t j = 0 ; j < columns ; j++ ) {
1270 arrayEz(i,j) += ezField;
1271 // This adds back the overall Z gradient of the field (main E field component)
1273 // Warning: (-=) assumes you are using an error potetial without the overall Field included
1276 // Integrate Er/Ez from Z to zero
1277 for ( Int_t j = 0 ; j < columns ; j++ ) {
1278 for ( Int_t i = 0 ; i < rows ; i++ ) {
1280 Int_t index = 1 ; // Simpsons rule if N=odd. If N!=odd then add extra point by trapezoidal rule.
1281 arrayErOverEz(i,j) = 0.0 ;
1282 arrayDeltaEz(i,j) = 0.0 ;
1284 for ( Int_t k = j ; k < columns ; k++ ) {
1285 arrayErOverEz(i,j) += index*(gridSizeZ/3.0)*arrayEr(i,k)/arrayEz(i,k) ;
1286 arrayDeltaEz(i,j) += index*(gridSizeZ/3.0)*(arrayEz(i,k)-ezField) ;
1287 if ( index != 4 ) index = 4; else index = 2 ;
1290 arrayErOverEz(i,j) -= (gridSizeZ/3.0)*arrayEr(i,columns-1)/arrayEz(i,columns-1) ;
1291 arrayDeltaEz(i,j) -= (gridSizeZ/3.0)*(arrayEz(i,columns-1)-ezField) ;
1294 arrayErOverEz(i,j) += (gridSizeZ/3.0) * ( 0.5*arrayEr(i,columns-2)/arrayEz(i,columns-2)
1295 -2.5*arrayEr(i,columns-1)/arrayEz(i,columns-1));
1296 arrayDeltaEz(i,j) += (gridSizeZ/3.0) * ( 0.5*(arrayEz(i,columns-2)-ezField)
1297 -2.5*(arrayEz(i,columns-1)-ezField));
1299 if ( j == columns-2 ) {
1300 arrayErOverEz(i,j) = (gridSizeZ/3.0) * ( 1.5*arrayEr(i,columns-2)/arrayEz(i,columns-2)
1301 +1.5*arrayEr(i,columns-1)/arrayEz(i,columns-1) ) ;
1302 arrayDeltaEz(i,j) = (gridSizeZ/3.0) * ( 1.5*(arrayEz(i,columns-2)-ezField)
1303 +1.5*(arrayEz(i,columns-1)-ezField) ) ;
1305 if ( j == columns-1 ) {
1306 arrayErOverEz(i,j) = 0.0 ;
1307 arrayDeltaEz(i,j) = 0.0 ;
1312 // calculate z distortion from the integrated Delta Ez residuals
1313 // and include the aquivalence (Volt to cm) of the ROC shift !!
1315 for ( Int_t j = 0 ; j < columns ; j++ ) {
1316 for ( Int_t i = 0 ; i < rows ; i++ ) {
1318 // Scale the Ez distortions with the drift velocity pertubation -> delivers cm
1319 arrayDeltaEz(i,j) = arrayDeltaEz(i,j)*fgkdvdE;
1321 // ROC Potential in cm aquivalent
1322 Double_t dzROCShift = arrayV(i, columns -1)/ezField;
1323 if ( rocDisplacement ) arrayDeltaEz(i,j) = arrayDeltaEz(i,j) + dzROCShift; // add the ROC misaligment
1333 void AliTPCCorrection::PoissonRelaxation3D( TMatrixD**arrayofArrayV, TMatrixD**arrayofChargeDensities,
1334 TMatrixD**arrayofEroverEz, TMatrixD**arrayofEPhioverEz, TMatrixD**arrayofDeltaEz,
1335 const Int_t rows, const Int_t columns, const Int_t phislices,
1336 const Float_t deltaphi, const Int_t iterations, const Int_t symmetry,
1337 Bool_t rocDisplacement ) {
1339 // 3D - Solve Poisson's Equation in 3D by Relaxation Technique
1341 // NOTE: In order for this algorith to work, the number of rows and columns must be a power of 2 plus one.
1342 // The number of rows and COLUMNS can be different.
1345 // COLUMNS == 2**N + 1
1346 // PHISLICES == Arbitrary but greater than 3
1348 // DeltaPhi in Radians
1350 // SYMMETRY = 0 if no phi symmetries, and no phi boundary conditions
1351 // = 1 if we have reflection symmetry at the boundaries (eg. sector symmetry or half sector symmetries).
1353 // NOTE: rocDisplacement is used to include (or ignore) the ROC misalignment in the dz calculation
1355 const Double_t ezField = (fgkCathodeV-fgkGG)/fgkTPCZ0; // = ALICE Electric Field (V/cm) Magnitude ~ -400 V/cm;
1357 const Float_t gridSizeR = (fgkOFCRadius-fgkIFCRadius) / (rows-1) ;
1358 const Float_t gridSizePhi = deltaphi ;
1359 const Float_t gridSizeZ = fgkTPCZ0 / (columns-1) ;
1360 const Float_t ratioPhi = gridSizeR*gridSizeR / (gridSizePhi*gridSizePhi) ;
1361 const Float_t ratioZ = gridSizeR*gridSizeR / (gridSizeZ*gridSizeZ) ;
1363 TMatrixD arrayE(rows,columns) ;
1365 // Check that the number of rows and columns is suitable for a binary expansion
1366 if ( !IsPowerOfTwo((rows-1)) ) {
1367 AliError("Poisson3DRelaxation - Error in the number of rows. Must be 2**M - 1");
1369 if ( !IsPowerOfTwo((columns-1)) ) {
1370 AliError("Poisson3DRelaxation - Error in the number of columns. Must be 2**N - 1");
1372 if ( phislices <= 3 ) {
1373 AliError("Poisson3DRelaxation - Error in the number of phislices. Must be larger than 3");
1375 if ( phislices > 1000 ) {
1376 AliError("Poisson3D phislices > 1000 is not allowed (nor wise) ");
1379 // Solve Poisson's equation in cylindrical coordinates by relaxation technique
1380 // Allow for different size grid spacing in R and Z directions
1381 // Use a binary expansion of the matrix to speed up the solution of the problem
1383 Int_t loops, mplus, mminus, signplus, signminus ;
1384 Int_t ione = (rows-1)/4 ;
1385 Int_t jone = (columns-1)/4 ;
1386 loops = TMath::Max(ione, jone) ; // Calculate the number of loops for the binary expansion
1387 loops = 1 + (int) ( 0.5 + TMath::Log2((double)loops) ) ; // Solve for N in 2**N
1389 TMatrixD* arrayofSumChargeDensities[1000] ; // Create temporary arrays to store low resolution charge arrays
1391 for ( Int_t i = 0 ; i < phislices ; i++ ) { arrayofSumChargeDensities[i] = new TMatrixD(rows,columns) ; }
1392 AliSysInfo::AddStamp("3DInit", 10,0,0);
1394 for ( Int_t count = 0 ; count < loops ; count++ ) { // START the master loop and do the binary expansion
1395 AliSysInfo::AddStamp("3Diter", 20,count,0);
1397 Float_t tempgridSizeR = gridSizeR * ione ;
1398 Float_t tempratioPhi = ratioPhi * ione * ione ; // Used tobe divided by ( m_one * m_one ) when m_one was != 1
1399 Float_t tempratioZ = ratioZ * ione * ione / ( jone * jone ) ;
1401 std::vector<float> coef1(rows) ; // Do this the standard C++ way to avoid gcc extensions for Float_t coef1[rows]
1402 std::vector<float> coef2(rows) ; // Do this the standard C++ way to avoid gcc extensions for Float_t coef1[rows]
1403 std::vector<float> coef3(rows) ; // Do this the standard C++ way to avoid gcc extensions for Float_t coef1[rows]
1404 std::vector<float> coef4(rows) ; // Do this the standard C++ way to avoid gcc extensions for Float_t coef1[rows]
1406 for ( Int_t i = ione ; i < rows-1 ; i+=ione ) {
1407 Float_t radius = fgkIFCRadius + i*gridSizeR ;
1408 coef1[i] = 1.0 + tempgridSizeR/(2*radius);
1409 coef2[i] = 1.0 - tempgridSizeR/(2*radius);
1410 coef3[i] = tempratioPhi/(radius*radius);
1411 coef4[i] = 0.5 / (1.0 + tempratioZ + coef3[i]);
1414 for ( Int_t m = 0 ; m < phislices ; m++ ) {
1415 TMatrixD &chargeDensity = *arrayofChargeDensities[m] ;
1416 TMatrixD &sumChargeDensity = *arrayofSumChargeDensities[m] ;
1417 for ( Int_t i = ione ; i < rows-1 ; i += ione ) {
1418 Float_t radius = fgkIFCRadius + i*gridSizeR ;
1419 for ( Int_t j = jone ; j < columns-1 ; j += jone ) {
1420 if ( ione == 1 && jone == 1 ) sumChargeDensity(i,j) = chargeDensity(i,j) ;
1421 else { // Add up all enclosed charge density contributions within 1/2 unit in all directions
1422 Float_t weight = 0.0 ;
1424 sumChargeDensity(i,j) = 0.0 ;
1425 for ( Int_t ii = i-ione/2 ; ii <= i+ione/2 ; ii++ ) {
1426 for ( Int_t jj = j-jone/2 ; jj <= j+jone/2 ; jj++ ) {
1427 if ( ii == i-ione/2 || ii == i+ione/2 || jj == j-jone/2 || jj == j+jone/2 ) weight = 0.5 ;
1430 sumChargeDensity(i,j) += chargeDensity(ii,jj)*weight*radius ;
1431 sum += weight*radius ;
1434 sumChargeDensity(i,j) /= sum ;
1436 sumChargeDensity(i,j) *= tempgridSizeR*tempgridSizeR; // just saving a step later on
1441 for ( Int_t k = 1 ; k <= iterations; k++ ) {
1443 // over-relaxation index, >= 1 but < 2
1444 Float_t overRelax = 1.0 + TMath::Sqrt( TMath::Cos( (k*TMath::PiOver2())/iterations ) ) ;
1445 Float_t overRelaxM1 = overRelax - 1.0 ;
1447 std::vector<float> overRelaxcoef4(rows) ; // Do this the standard C++ way to avoid gcc extensions
1448 std::vector<float> overRelaxcoef5(rows) ; // Do this the standard C++ way to avoid gcc extensions
1450 for ( Int_t i = ione ; i < rows-1 ; i+=ione ) {
1451 overRelaxcoef4[i] = overRelax * coef4[i] ;
1452 overRelaxcoef5[i] = overRelaxM1 / overRelaxcoef4[i] ;
1455 for ( Int_t m = 0 ; m < phislices ; m++ ) {
1457 mplus = m + 1; signplus = 1 ;
1458 mminus = m - 1 ; signminus = 1 ;
1459 if (symmetry==1) { // Reflection symmetry in phi (e.g. symmetry at sector boundaries, or half sectors, etc.)
1460 if ( mplus > phislices-1 ) mplus = phislices - 2 ;
1461 if ( mminus < 0 ) mminus = 1 ;
1463 else if (symmetry==-1) { // Anti-symmetry in phi
1464 if ( mplus > phislices-1 ) { mplus = phislices - 2 ; signplus = -1 ; }
1465 if ( mminus < 0 ) { mminus = 1 ; signminus = -1 ; }
1467 else { // No Symmetries in phi, no boundaries, the calculation is continuous across all phi
1468 if ( mplus > phislices-1 ) mplus = m + 1 - phislices ;
1469 if ( mminus < 0 ) mminus = m - 1 + phislices ;
1471 TMatrixD& arrayV = *arrayofArrayV[m] ;
1472 TMatrixD& arrayVP = *arrayofArrayV[mplus] ;
1473 TMatrixD& arrayVM = *arrayofArrayV[mminus] ;
1474 TMatrixD& sumChargeDensity = *arrayofSumChargeDensities[m] ;
1475 Double_t *arrayVfast = arrayV.GetMatrixArray();
1476 Double_t *arrayVPfast = arrayVP.GetMatrixArray();
1477 Double_t *arrayVMfast = arrayVM.GetMatrixArray();
1478 Double_t *sumChargeDensityFast=sumChargeDensity.GetMatrixArray();
1481 // slow implementation
1482 for ( Int_t i = ione ; i < rows-1 ; i+=ione ) {
1483 for ( Int_t j = jone ; j < columns-1 ; j+=jone ) {
1485 arrayV(i,j) = ( coef2[i] * arrayV(i-ione,j)
1486 + tempratioZ * ( arrayV(i,j-jone) + arrayV(i,j+jone) )
1487 - overRelaxcoef5[i] * arrayV(i,j)
1488 + coef1[i] * arrayV(i+ione,j)
1489 + coef3[i] * ( signplus*arrayVP(i,j) + signminus*arrayVM(i,j) )
1490 + sumChargeDensity(i,j)
1491 ) * overRelaxcoef4[i] ;
1492 // Note: over-relax the solution at each step. This speeds up the convergance.
1496 for ( Int_t i = ione ; i < rows-1 ; i+=ione ) {
1497 Double_t *arrayVfastI = &(arrayVfast[i*columns]);
1498 Double_t *arrayVPfastI = &(arrayVPfast[i*columns]);
1499 Double_t *arrayVMfastI = &(arrayVMfast[i*columns]);
1500 Double_t *sumChargeDensityFastI=&(sumChargeDensityFast[i*columns]);
1501 for ( Int_t j = jone ; j < columns-1 ; j+=jone ) {
1502 Double_t resSlow,resFast;
1503 // resSlow = ( coef2[i] * arrayV(i-ione,j)
1504 // + tempratioZ * ( arrayV(i,j-jone) + arrayV(i,j+jone) )
1505 // - overRelaxcoef5[i] * arrayV(i,j)
1506 // + coef1[i] * arrayV(i+ione,j)
1507 // + coef3[i] * ( signplus*arrayVP(i,j) + signminus*arrayVM(i,j) )
1508 // + sumChargeDensity(i,j)
1509 // ) * overRelaxcoef4[i] ;
1510 resFast = ( coef2[i] * arrayVfastI[j-columns*ione]
1511 + tempratioZ * ( arrayVfastI[j-jone] + arrayVfastI[j+jone] )
1512 - overRelaxcoef5[i] * arrayVfastI[j]
1513 + coef1[i] * arrayVfastI[j+columns*ione]
1514 + coef3[i] * ( signplus* arrayVPfastI[j] + signminus*arrayVMfastI[j])
1515 + sumChargeDensityFastI[j]
1516 ) * overRelaxcoef4[i] ;
1517 // if (resSlow!=resFast){
1518 // printf("problem\t%d\t%d\t%f\t%f\t%f\n",i,j,resFast,resSlow,resFast-resSlow);
1520 arrayVfastI[j]=resFast;
1521 // Note: over-relax the solution at each step. This speeds up the convergance.
1526 if ( k == iterations ) { // After full solution is achieved, copy low resolution solution into higher res array
1527 for ( Int_t i = ione ; i < rows-1 ; i+=ione ) {
1528 for ( Int_t j = jone ; j < columns-1 ; j+=jone ) {
1531 arrayV(i+ione/2,j) = ( arrayV(i+ione,j) + arrayV(i,j) ) / 2 ;
1532 if ( i == ione ) arrayV(i-ione/2,j) = ( arrayV(0,j) + arrayV(ione,j) ) / 2 ;
1535 arrayV(i,j+jone/2) = ( arrayV(i,j+jone) + arrayV(i,j) ) / 2 ;
1536 if ( j == jone ) arrayV(i,j-jone/2) = ( arrayV(i,0) + arrayV(i,jone) ) / 2 ;
1538 if ( ione > 1 && jone > 1 ) {
1539 arrayV(i+ione/2,j+jone/2) = ( arrayV(i+ione,j+jone) + arrayV(i,j) ) / 2 ;
1540 if ( i == ione ) arrayV(i-ione/2,j-jone/2) = ( arrayV(0,j-jone) + arrayV(ione,j) ) / 2 ;
1541 if ( j == jone ) arrayV(i-ione/2,j-jone/2) = ( arrayV(i-ione,0) + arrayV(i,jone) ) / 2 ;
1542 // Note that this leaves a point at the upper left and lower right corners uninitialized. Not a big deal.
1551 ione = ione / 2 ; if ( ione < 1 ) ione = 1 ;
1552 jone = jone / 2 ; if ( jone < 1 ) jone = 1 ;
1556 //Differentiate V(r) and solve for E(r) using special equations for the first and last row
1557 //Integrate E(r)/E(z) from point of origin to pad plane
1558 AliSysInfo::AddStamp("CalcField", 100,0,0);
1560 for ( Int_t m = 0 ; m < phislices ; m++ ) {
1561 TMatrixD& arrayV = *arrayofArrayV[m] ;
1562 TMatrixD& eroverEz = *arrayofEroverEz[m] ;
1564 for ( Int_t j = columns-1 ; j >= 0 ; j-- ) { // Count backwards to facilitate integration over Z
1566 // Differentiate in R
1567 for ( Int_t i = 1 ; i < rows-1 ; i++ ) arrayE(i,j) = -1 * ( arrayV(i+1,j) - arrayV(i-1,j) ) / (2*gridSizeR) ;
1568 arrayE(0,j) = -1 * ( -0.5*arrayV(2,j) + 2.0*arrayV(1,j) - 1.5*arrayV(0,j) ) / gridSizeR ;
1569 arrayE(rows-1,j) = -1 * ( 1.5*arrayV(rows-1,j) - 2.0*arrayV(rows-2,j) + 0.5*arrayV(rows-3,j) ) / gridSizeR ;
1571 for ( Int_t i = 0 ; i < rows ; i++ ) {
1572 Int_t index = 1 ; // Simpsons rule if N=odd. If N!=odd then add extra point by trapezoidal rule.
1573 eroverEz(i,j) = 0.0 ;
1574 for ( Int_t k = j ; k < columns ; k++ ) {
1576 eroverEz(i,j) += index*(gridSizeZ/3.0)*arrayE(i,k)/(-1*ezField) ;
1577 if ( index != 4 ) index = 4; else index = 2 ;
1579 if ( index == 4 ) eroverEz(i,j) -= (gridSizeZ/3.0)*arrayE(i,columns-1)/ (-1*ezField) ;
1580 if ( index == 2 ) eroverEz(i,j) +=
1581 (gridSizeZ/3.0)*(0.5*arrayE(i,columns-2)-2.5*arrayE(i,columns-1))/(-1*ezField) ;
1582 if ( j == columns-2 ) eroverEz(i,j) =
1583 (gridSizeZ/3.0)*(1.5*arrayE(i,columns-2)+1.5*arrayE(i,columns-1))/(-1*ezField) ;
1584 if ( j == columns-1 ) eroverEz(i,j) = 0.0 ;
1587 // if ( m == 0 ) { TCanvas* c1 = new TCanvas("erOverEz","erOverEz",50,50,840,600) ; c1 -> cd() ;
1588 // eroverEz.Draw("surf") ; } // JT test
1590 AliSysInfo::AddStamp("IntegrateEr", 120,0,0);
1592 //Differentiate V(r) and solve for E(phi)
1593 //Integrate E(phi)/E(z) from point of origin to pad plane
1595 for ( Int_t m = 0 ; m < phislices ; m++ ) {
1597 mplus = m + 1; signplus = 1 ;
1598 mminus = m - 1 ; signminus = 1 ;
1599 if (symmetry==1) { // Reflection symmetry in phi (e.g. symmetry at sector boundaries, or half sectors, etc.)
1600 if ( mplus > phislices-1 ) mplus = phislices - 2 ;
1601 if ( mminus < 0 ) mminus = 1 ;
1603 else if (symmetry==-1) { // Anti-symmetry in phi
1604 if ( mplus > phislices-1 ) { mplus = phislices - 2 ; signplus = -1 ; }
1605 if ( mminus < 0 ) { mminus = 1 ; signminus = -1 ; }
1607 else { // No Symmetries in phi, no boundaries, the calculations is continuous across all phi
1608 if ( mplus > phislices-1 ) mplus = m + 1 - phislices ;
1609 if ( mminus < 0 ) mminus = m - 1 + phislices ;
1611 TMatrixD &arrayVP = *arrayofArrayV[mplus] ;
1612 TMatrixD &arrayVM = *arrayofArrayV[mminus] ;
1613 TMatrixD &ePhioverEz = *arrayofEPhioverEz[m] ;
1614 for ( Int_t j = columns-1 ; j >= 0 ; j-- ) { // Count backwards to facilitate integration over Z
1615 // Differentiate in Phi
1616 for ( Int_t i = 0 ; i < rows ; i++ ) {
1617 Float_t radius = fgkIFCRadius + i*gridSizeR ;
1618 arrayE(i,j) = -1 * (signplus * arrayVP(i,j) - signminus * arrayVM(i,j) ) / (2*radius*gridSizePhi) ;
1621 for ( Int_t i = 0 ; i < rows ; i++ ) {
1622 Int_t index = 1 ; // Simpsons rule if N=odd. If N!=odd then add extra point by trapezoidal rule.
1623 ePhioverEz(i,j) = 0.0 ;
1624 for ( Int_t k = j ; k < columns ; k++ ) {
1626 ePhioverEz(i,j) += index*(gridSizeZ/3.0)*arrayE(i,k)/(-1*ezField) ;
1627 if ( index != 4 ) index = 4; else index = 2 ;
1629 if ( index == 4 ) ePhioverEz(i,j) -= (gridSizeZ/3.0)*arrayE(i,columns-1)/ (-1*ezField) ;
1630 if ( index == 2 ) ePhioverEz(i,j) +=
1631 (gridSizeZ/3.0)*(0.5*arrayE(i,columns-2)-2.5*arrayE(i,columns-1))/(-1*ezField) ;
1632 if ( j == columns-2 ) ePhioverEz(i,j) =
1633 (gridSizeZ/3.0)*(1.5*arrayE(i,columns-2)+1.5*arrayE(i,columns-1))/(-1*ezField) ;
1634 if ( j == columns-1 ) ePhioverEz(i,j) = 0.0 ;
1637 // if ( m == 5 ) { TCanvas* c2 = new TCanvas("arrayE","arrayE",50,50,840,600) ; c2 -> cd() ;
1638 // arrayE.Draw("surf") ; } // JT test
1640 AliSysInfo::AddStamp("IntegrateEphi", 130,0,0);
1643 // Differentiate V(r) and solve for E(z) using special equations for the first and last row
1644 // Integrate (E(z)-Ezstd) from point of origin to pad plane
1646 for ( Int_t m = 0 ; m < phislices ; m++ ) {
1647 TMatrixD& arrayV = *arrayofArrayV[m] ;
1648 TMatrixD& deltaEz = *arrayofDeltaEz[m] ;
1650 // Differentiate V(z) and solve for E(z) using special equations for the first and last columns
1651 for ( Int_t i = 0 ; i < rows ; i++) {
1652 for ( Int_t j = 1 ; j < columns-1 ; j++ ) arrayE(i,j) = -1 * ( arrayV(i,j+1) - arrayV(i,j-1) ) / (2*gridSizeZ) ;
1653 arrayE(i,0) = -1 * ( -0.5*arrayV(i,2) + 2.0*arrayV(i,1) - 1.5*arrayV(i,0) ) / gridSizeZ ;
1654 arrayE(i,columns-1) = -1 * ( 1.5*arrayV(i,columns-1) - 2.0*arrayV(i,columns-2) + 0.5*arrayV(i,columns-3) ) / gridSizeZ ;
1657 for ( Int_t j = columns-1 ; j >= 0 ; j-- ) { // Count backwards to facilitate integration over Z
1659 for ( Int_t i = 0 ; i < rows ; i++ ) {
1660 Int_t index = 1 ; // Simpsons rule if N=odd. If N!=odd then add extra point by trapezoidal rule.
1661 deltaEz(i,j) = 0.0 ;
1662 for ( Int_t k = j ; k < columns ; k++ ) {
1663 deltaEz(i,j) += index*(gridSizeZ/3.0)*arrayE(i,k) ;
1664 if ( index != 4 ) index = 4; else index = 2 ;
1666 if ( index == 4 ) deltaEz(i,j) -= (gridSizeZ/3.0)*arrayE(i,columns-1) ;
1667 if ( index == 2 ) deltaEz(i,j) +=
1668 (gridSizeZ/3.0)*(0.5*arrayE(i,columns-2)-2.5*arrayE(i,columns-1)) ;
1669 if ( j == columns-2 ) deltaEz(i,j) =
1670 (gridSizeZ/3.0)*(1.5*arrayE(i,columns-2)+1.5*arrayE(i,columns-1)) ;
1671 if ( j == columns-1 ) deltaEz(i,j) = 0.0 ;
1675 // if ( m == 0 ) { TCanvas* c1 = new TCanvas("erOverEz","erOverEz",50,50,840,600) ; c1 -> cd() ;
1676 // eroverEz.Draw("surf") ; } // JT test
1678 // calculate z distortion from the integrated Delta Ez residuals
1679 // and include the aquivalence (Volt to cm) of the ROC shift !!
1681 for ( Int_t j = 0 ; j < columns ; j++ ) {
1682 for ( Int_t i = 0 ; i < rows ; i++ ) {
1684 // Scale the Ez distortions with the drift velocity pertubation -> delivers cm
1685 deltaEz(i,j) = deltaEz(i,j)*fgkdvdE;
1687 // ROC Potential in cm aquivalent
1688 Double_t dzROCShift = arrayV(i, columns -1)/ezField;
1689 if ( rocDisplacement ) deltaEz(i,j) = deltaEz(i,j) + dzROCShift; // add the ROC misaligment
1694 } // end loop over phi
1695 AliSysInfo::AddStamp("IntegrateEz", 140,0,0);
1698 for ( Int_t k = 0 ; k < phislices ; k++ )
1700 arrayofSumChargeDensities[k]->Delete() ;
1709 Int_t AliTPCCorrection::IsPowerOfTwo(Int_t i) const {
1711 // Helperfunction: Check if integer is a power of 2
1714 while( i > 0 ) { j += (i&1) ; i = (i>>1) ; }
1715 if ( j == 1 ) return(1) ; // True
1716 return(0) ; // False
1720 AliExternalTrackParam * AliTPCCorrection::FitDistortedTrack(AliExternalTrackParam & trackIn, Double_t refX, Int_t dir, TTreeSRedirector * const pcstream){
1722 // Fit the track parameters - without and with distortion
1723 // 1. Space points in the TPC are simulated along the trajectory
1724 // 2. Space points distorted
1725 // 3. Fits the non distorted and distroted track to the reference plane at refX
1726 // 4. For visualization and debugging purposes the space points and tracks can be stored in the tree - using the TTreeSRedirector functionality
1728 // trackIn - input track parameters
1729 // refX - reference X to fit the track
1730 // dir - direction - out=1 or in=-1
1731 // pcstream - debug streamer to check the results
1733 // see AliExternalTrackParam.h documentation:
1734 // track1.fP[0] - local y (rphi)
1736 // track1.fP[2] - sinus of local inclination angle
1737 // track1.fP[3] - tangent of deep angle
1738 // track1.fP[4] - 1/pt
1740 AliTPCROC * roc = AliTPCROC::Instance();
1741 const Int_t npoints0=roc->GetNRows(0)+roc->GetNRows(36);
1742 const Double_t kRTPC0 =roc->GetPadRowRadii(0,0);
1743 const Double_t kRTPC1 =roc->GetPadRowRadii(36,roc->GetNRows(36)-1);
1744 const Double_t kMaxSnp = 0.85;
1745 const Double_t kSigmaY=0.1;
1746 const Double_t kSigmaZ=0.1;
1747 const Double_t kMaxR=500;
1748 const Double_t kMaxZ=500;
1750 const Double_t kMaxZ0=220;
1751 const Double_t kZcut=3;
1752 const Double_t kMass = TDatabasePDG::Instance()->GetParticle("pi+")->Mass();
1756 AliExternalTrackParam track(trackIn); //
1758 AliTrackPointArray pointArray0(npoints0);
1759 AliTrackPointArray pointArray1(npoints0);
1761 if (!AliTrackerBase::PropagateTrackTo(&track,kRTPC0,kMass,5,kTRUE,kMaxSnp)) return 0;
1763 // simulate the track
1765 Float_t covPoint[6]={0,0,0, kSigmaY*kSigmaY,0,kSigmaZ*kSigmaZ}; //covariance at the local frame
1766 for (Double_t radius=kRTPC0; radius<kRTPC1; radius++){
1767 if (!AliTrackerBase::PropagateTrackTo(&track,radius,kMass,5,kTRUE,kMaxSnp)) return 0;
1769 xyz[0]+=gRandom->Gaus(0,0.000005);
1770 xyz[1]+=gRandom->Gaus(0,0.000005);
1771 xyz[2]+=gRandom->Gaus(0,0.000005);
1772 if (TMath::Abs(track.GetZ())>kMaxZ0) continue;
1773 if (TMath::Abs(track.GetX())<kRTPC0) continue;
1774 if (TMath::Abs(track.GetX())>kRTPC1) continue;
1775 AliTrackPoint pIn0; // space point
1777 Int_t sector= (xyz[2]>0)? 0:18;
1778 pointArray0.GetPoint(pIn0,npoints);
1779 pointArray1.GetPoint(pIn1,npoints);
1780 Double_t alpha = TMath::ATan2(xyz[1],xyz[0]);
1781 Float_t distPoint[3]={xyz[0],xyz[1],xyz[2]};
1782 DistortPoint(distPoint, sector);
1783 pIn0.SetXYZ(xyz[0], xyz[1],xyz[2]);
1784 pIn1.SetXYZ(distPoint[0], distPoint[1],distPoint[2]);
1786 track.Rotate(alpha);
1787 AliTrackPoint prot0 = pIn0.Rotate(alpha); // rotate to the local frame - non distoted point
1788 AliTrackPoint prot1 = pIn1.Rotate(alpha); // rotate to the local frame - distorted point
1789 prot0.SetXYZ(prot0.GetX(),prot0.GetY(), prot0.GetZ(),covPoint);
1790 prot1.SetXYZ(prot1.GetX(),prot1.GetY(), prot1.GetZ(),covPoint);
1791 pIn0=prot0.Rotate(-alpha); // rotate back to global frame
1792 pIn1=prot1.Rotate(-alpha); // rotate back to global frame
1793 pointArray0.AddPoint(npoints, &pIn0);
1794 pointArray1.AddPoint(npoints, &pIn1);
1796 if (npoints>=npoints0) break;
1798 if (npoints<npoints0/4.) return 0;
1802 AliExternalTrackParam *track0=0;
1803 AliExternalTrackParam *track1=0;
1804 AliTrackPoint point1,point2,point3;
1805 if (dir==1) { //make seed inner
1806 pointArray0.GetPoint(point1,1);
1807 pointArray0.GetPoint(point2,11);
1808 pointArray0.GetPoint(point3,21);
1810 if (dir==-1){ //make seed outer
1811 pointArray0.GetPoint(point1,npoints-21);
1812 pointArray0.GetPoint(point2,npoints-11);
1813 pointArray0.GetPoint(point3,npoints-1);
1815 if ((TMath::Abs(point1.GetX()-point3.GetX())+TMath::Abs(point1.GetY()-point3.GetY()))<10){
1816 printf("fit points not properly initialized\n");
1819 track0 = AliTrackerBase::MakeSeed(point1, point2, point3);
1820 track1 = AliTrackerBase::MakeSeed(point1, point2, point3);
1821 track0->ResetCovariance(10);
1822 track1->ResetCovariance(10);
1823 if (TMath::Abs(AliTrackerBase::GetBz())<0.01){
1824 ((Double_t*)track0->GetParameter())[4]= trackIn.GetParameter()[4];
1825 ((Double_t*)track1->GetParameter())[4]= trackIn.GetParameter()[4];
1827 for (Int_t jpoint=0; jpoint<npoints; jpoint++){
1828 Int_t ipoint= (dir>0) ? jpoint: npoints-1-jpoint;
1832 pointArray0.GetPoint(pIn0,ipoint);
1833 pointArray1.GetPoint(pIn1,ipoint);
1834 AliTrackPoint prot0 = pIn0.Rotate(track0->GetAlpha()); // rotate to the local frame - non distoted point
1835 AliTrackPoint prot1 = pIn1.Rotate(track1->GetAlpha()); // rotate to the local frame - distorted point
1836 if (TMath::Abs(prot0.GetX())<kRTPC0) continue;
1837 if (TMath::Abs(prot0.GetX())>kRTPC1) continue;
1839 if (!AliTrackerBase::PropagateTrackTo(track0,prot0.GetX(),kMass,5,kFALSE,kMaxSnp)) break;
1840 if (!AliTrackerBase::PropagateTrackTo(track1,prot0.GetX(),kMass,5,kFALSE,kMaxSnp)) break;
1841 if (TMath::Abs(track0->GetZ())>kMaxZ) break;
1842 if (TMath::Abs(track0->GetX())>kMaxR) break;
1843 if (TMath::Abs(track1->GetZ())>kMaxZ) break;
1844 if (TMath::Abs(track1->GetX())>kMaxR) break;
1845 if (dir>0 && track1->GetX()>refX) continue;
1846 if (dir<0 && track1->GetX()<refX) continue;
1847 if (TMath::Abs(track1->GetZ())<kZcut)continue;
1848 track.GetXYZ(xyz); // distorted track also propagated to the same reference radius
1850 Double_t pointPos[2]={0,0};
1851 Double_t pointCov[3]={0,0,0};
1852 pointPos[0]=prot0.GetY();//local y
1853 pointPos[1]=prot0.GetZ();//local z
1854 pointCov[0]=prot0.GetCov()[3];//simay^2
1855 pointCov[1]=prot0.GetCov()[4];//sigmayz
1856 pointCov[2]=prot0.GetCov()[5];//sigmaz^2
1857 if (!track0->Update(pointPos,pointCov)) break;
1859 Double_t deltaX=prot1.GetX()-prot0.GetX(); // delta X
1860 Double_t deltaYX=deltaX*TMath::Tan(TMath::ASin(track1->GetSnp())); // deltaY due delta X
1861 Double_t deltaZX=deltaX*track1->GetTgl(); // deltaZ due delta X
1863 pointPos[0]=prot1.GetY()-deltaYX;//local y is sign correct? should be minus
1864 pointPos[1]=prot1.GetZ()-deltaZX;//local z is sign correct? should be minus
1865 pointCov[0]=prot1.GetCov()[3];//simay^2
1866 pointCov[1]=prot1.GetCov()[4];//sigmayz
1867 pointCov[2]=prot1.GetCov()[5];//sigmaz^2
1868 if (!track1->Update(pointPos,pointCov)) break;
1872 if (npoints2<npoints/4.) return 0;
1873 AliTrackerBase::PropagateTrackTo(track0,refX,kMass,5.,kTRUE,kMaxSnp);
1874 AliTrackerBase::PropagateTrackTo(track0,refX,kMass,1.,kTRUE,kMaxSnp);
1875 track1->Rotate(track0->GetAlpha());
1876 AliTrackerBase::PropagateTrackTo(track1,track0->GetX(),kMass,5.,kFALSE,kMaxSnp);
1878 if (pcstream) (*pcstream)<<Form("fitDistort%s",GetName())<<
1879 "point0.="<<&pointArray0<< // points
1880 "point1.="<<&pointArray1<< // distorted points
1881 "trackIn.="<<&track<< // original track
1882 "track0.="<<track0<< // fitted track
1883 "track1.="<<track1<< // fitted distorted track
1885 new(&trackIn) AliExternalTrackParam(*track0);
1894 TTree* AliTPCCorrection::CreateDistortionTree(Double_t step){
1896 // create the distortion tree on a mesh with granularity given by step
1897 // return the tree with distortions at given position
1898 // Map is created on the mesh with given step size
1900 TTreeSRedirector *pcstream = new TTreeSRedirector(Form("correction%s.root",GetName()));
1902 for (Double_t x= -250; x<250; x+=step){
1903 for (Double_t y= -250; y<250; y+=step){
1904 Double_t r = TMath::Sqrt(x*x+y*y);
1906 if (r>250) continue;
1907 for (Double_t z= -250; z<250; z+=step){
1908 Int_t roc=(z>0)?0:18;
1912 Double_t phi = TMath::ATan2(y,x);
1913 DistortPoint(xyz,roc);
1914 Double_t r1 = TMath::Sqrt(xyz[0]*xyz[0]+xyz[1]*xyz[1]);
1915 Double_t phi1 = TMath::ATan2(xyz[1],xyz[0]);
1916 if ((phi1-phi)>TMath::Pi()) phi1-=TMath::Pi();
1917 if ((phi1-phi)<-TMath::Pi()) phi1+=TMath::Pi();
1918 Double_t dx = xyz[0]-x;
1919 Double_t dy = xyz[1]-y;
1920 Double_t dz = xyz[2]-z;
1922 Double_t drphi=(phi1-phi)*r;
1923 (*pcstream)<<"distortion"<<
1924 "x="<<x<< // original position
1929 "x1="<<xyz[0]<< // distorted position
1935 "dx="<<dx<< // delta position
1945 TFile f(Form("correction%s.root",GetName()));
1946 TTree * tree = (TTree*)f.Get("distortion");
1947 TTree * tree2= tree->CopyTree("1");
1948 tree2->SetName(Form("dist%s",GetName()));
1949 tree2->SetDirectory(0);
1957 void AliTPCCorrection::MakeTrackDistortionTree(TTree *tinput, Int_t dtype, Int_t ptype, const TObjArray * corrArray, Int_t step, Int_t offset, Bool_t debug ){
1960 // For each partial correction (specified in array) and given track topology (phi, theta, snp, refX)
1961 // calculates partial distortions
1962 // Partial distortion is stored in the resulting tree
1963 // Output is storred in the file distortion_<dettype>_<partype>.root
1964 // Partial distortion is stored with the name given by correction name
1967 // Parameters of function:
1968 // input - input tree
1969 // dtype - distortion type 0 - ITSTPC, 1 -TPCTRD, 2 - TPCvertex , 3 - TPC-TOF, 4 - TPCTPC track crossing
1970 // ppype - parameter type
1971 // corrArray - array with partial corrections
1972 // step - skipe entries - if 1 all entries processed - it is slow
1973 // debug 0 if debug on also space points dumped - it is slow
1975 const Double_t kMaxSnp = 0.85;
1976 const Double_t kcutSnp=0.25;
1977 const Double_t kcutTheta=1.;
1978 const Double_t kRadiusTPC=85;
1979 // AliTPCROC *tpcRoc =AliTPCROC::Instance();
1981 const Double_t kMass = TDatabasePDG::Instance()->GetParticle("pi+")->Mass();
1982 // const Double_t kB2C=-0.299792458e-3;
1983 const Int_t kMinEntries=20;
1984 Double_t phi,theta, snp, mean,rms, entries,sector,dsec;
1987 tinput->SetBranchAddress("run",&run);
1988 tinput->SetBranchAddress("theta",&theta);
1989 tinput->SetBranchAddress("phi", &phi);
1990 tinput->SetBranchAddress("snp",&snp);
1991 tinput->SetBranchAddress("mean",&mean);
1992 tinput->SetBranchAddress("rms",&rms);
1993 tinput->SetBranchAddress("entries",&entries);
1994 tinput->SetBranchAddress("sector",§or);
1995 tinput->SetBranchAddress("dsec",&dsec);
1996 tinput->SetBranchAddress("refX",&refX);
1997 TTreeSRedirector *pcstream = new TTreeSRedirector(Form("distortion%d_%d_%d.root",dtype,ptype,offset));
1999 Int_t nentries=tinput->GetEntries();
2000 Int_t ncorr=corrArray->GetEntries();
2001 Double_t corrections[100]={0}; //
2003 Double_t cov[15]={0,0,0,0,0,0,0,0,0,0,0,0,0,0};
2005 if (dtype==5 || dtype==6) dtype=4;
2006 if (dtype==0) { dir=-1;}
2007 if (dtype==1) { dir=1;}
2008 if (dtype==2) { dir=-1;}
2009 if (dtype==3) { dir=1;}
2010 if (dtype==4) { dir=-1;}
2012 for (Int_t ientry=offset; ientry<nentries; ientry+=step){
2013 tinput->GetEntry(ientry);
2014 if (TMath::Abs(snp)>kMaxSnp) continue;
2017 if (dtype==2) tPar[1]=theta*kRadiusTPC;
2020 tPar[4]=(gRandom->Rndm()-0.5)*0.02; // should be calculated - non equal to 0
2022 // tracks crossing CE
2023 tPar[1]=0; // track at the CE
2024 //if (TMath::Abs(theta) <0.05) continue; // deep cross
2027 if (TMath::Abs(snp) >kcutSnp) continue;
2028 if (TMath::Abs(theta) >kcutTheta) continue;
2029 printf("%f\t%f\t%f\t%f\t%f\t%f\n",entries, sector,theta,snp, mean,rms);
2030 Double_t bz=AliTrackerBase::GetBz();
2031 if (dtype !=4) { //exclude TPC - for TPC mainly non primary tracks
2032 if (dtype!=2 && TMath::Abs(bz)>0.1 ) tPar[4]=snp/(refX*bz*kB2C*2);
2034 if (dtype==2 && TMath::Abs(bz)>0.1 ) {
2035 tPar[4]=snp/(kRadiusTPC*bz*kB2C*2);//
2036 // snp at the TPC inner radius in case the vertex match used
2040 tPar[4]+=(gRandom->Rndm()-0.5)*0.02;
2041 AliExternalTrackParam track(refX,phi,tPar,cov);
2045 Double_t pt=1./tPar[4];
2046 Double_t dRrec=0; // dummy value - needed for points - e.g for laser
2047 //if (ptype==4 &&bz<0) mean*=-1; // interpret as curvature -- COMMENTED out - in lookup signed 1/pt used
2048 Double_t refXD=refX;
2049 (*pcstream)<<"fit"<<
2050 "run="<<run<< // run number
2051 "bz="<<bz<< // magnetic filed used
2052 "dtype="<<dtype<< // detector match type
2053 "ptype="<<ptype<< // parameter type
2054 "theta="<<theta<< // theta
2055 "phi="<<phi<< // phi
2056 "snp="<<snp<< // snp
2057 "mean="<<mean<< // mean dist value
2058 "rms="<<rms<< // rms
2061 "refX="<<refXD<< // referece X as double
2062 "gx="<<xyz[0]<< // global position at reference
2063 "gy="<<xyz[1]<< // global position at reference
2064 "gz="<<xyz[2]<< // global position at reference
2065 "dRrec="<<dRrec<< // delta Radius in reconstruction
2067 "id="<<id<< // track id
2068 "entries="<<entries;// number of entries in bin
2071 if (entries<kMinEntries) isOK=kFALSE;
2073 if (dtype!=4) for (Int_t icorr=0; icorr<ncorr; icorr++) {
2074 AliTPCCorrection *corr = (AliTPCCorrection*)corrArray->At(icorr);
2075 corrections[icorr]=0;
2076 if (entries>kMinEntries){
2077 AliExternalTrackParam trackIn(refX,phi,tPar,cov);
2078 AliExternalTrackParam *trackOut = 0;
2079 if (debug) trackOut=corr->FitDistortedTrack(trackIn, refX, dir,pcstream);
2080 if (!debug) trackOut=corr->FitDistortedTrack(trackIn, refX, dir,0);
2081 if (dtype==0) {dir= -1;}
2082 if (dtype==1) {dir= 1;}
2083 if (dtype==2) {dir= -1;}
2084 if (dtype==3) {dir= 1;}
2087 if (!AliTrackerBase::PropagateTrackTo(&trackIn,refX,kMass,5,kTRUE,kMaxSnp)) isOK=kFALSE;
2088 if (!trackOut->Rotate(trackIn.GetAlpha())) isOK=kFALSE;
2089 if (!AliTrackerBase::PropagateTrackTo(trackOut,trackIn.GetX(),kMass,5,kFALSE,kMaxSnp)) isOK=kFALSE;
2090 // trackOut->PropagateTo(trackIn.GetX(),AliTrackerBase::GetBz());
2092 corrections[icorr]= trackOut->GetParameter()[ptype]-trackIn.GetParameter()[ptype];
2095 corrections[icorr]=0;
2098 //if (ptype==4 &&bz<0) corrections[icorr]*=-1; // interpret as curvature - commented out
2100 (*pcstream)<<"fit"<<
2101 Form("%s=",corr->GetName())<<corrections[icorr]; // dump correction value
2104 if (dtype==4) for (Int_t icorr=0; icorr<ncorr; icorr++) {
2106 // special case of the TPC tracks crossing the CE
2108 AliTPCCorrection *corr = (AliTPCCorrection*)corrArray->At(icorr);
2109 corrections[icorr]=0;
2110 if (entries>kMinEntries){
2111 AliExternalTrackParam trackIn0(refX,phi,tPar,cov); //Outer - direction to vertex
2112 AliExternalTrackParam trackIn1(refX,phi,tPar,cov); //Inner - direction magnet
2113 AliExternalTrackParam *trackOut0 = 0;
2114 AliExternalTrackParam *trackOut1 = 0;
2116 if (debug) trackOut0=corr->FitDistortedTrack(trackIn0, refX, dir,pcstream);
2117 if (!debug) trackOut0=corr->FitDistortedTrack(trackIn0, refX, dir,0);
2118 if (debug) trackOut1=corr->FitDistortedTrack(trackIn1, refX, -dir,pcstream);
2119 if (!debug) trackOut1=corr->FitDistortedTrack(trackIn1, refX, -dir,0);
2121 if (trackOut0 && trackOut1){
2122 if (!AliTrackerBase::PropagateTrackTo(&trackIn0,refX,kMass,5,kTRUE,kMaxSnp)) isOK=kFALSE;
2123 if (!AliTrackerBase::PropagateTrackTo(&trackIn0,refX,kMass,1,kFALSE,kMaxSnp)) isOK=kFALSE;
2124 if (!trackOut0->Rotate(trackIn0.GetAlpha())) isOK=kFALSE;
2125 if (!AliTrackerBase::PropagateTrackTo(trackOut0,trackIn0.GetX(),kMass,5,kFALSE,kMaxSnp)) isOK=kFALSE;
2127 if (!AliTrackerBase::PropagateTrackTo(&trackIn1,refX,kMass,5,kTRUE,kMaxSnp)) isOK=kFALSE;
2128 if (!trackIn1.Rotate(trackIn0.GetAlpha())) isOK=kFALSE;
2129 if (!AliTrackerBase::PropagateTrackTo(&trackIn1,trackIn0.GetX(),kMass,1,kFALSE,kMaxSnp)) isOK=kFALSE;
2130 if (!trackOut1->Rotate(trackIn1.GetAlpha())) isOK=kFALSE;
2131 if (!AliTrackerBase::PropagateTrackTo(trackOut1,trackIn1.GetX(),kMass,5,kFALSE,kMaxSnp)) isOK=kFALSE;
2133 corrections[icorr] = (trackOut0->GetParameter()[ptype]-trackIn0.GetParameter()[ptype]);
2134 corrections[icorr]-= (trackOut1->GetParameter()[ptype]-trackIn1.GetParameter()[ptype]);
2136 if ((TMath::Abs(trackOut0->GetX()-trackOut1->GetX())>0.1)||
2137 (TMath::Abs(trackOut0->GetX()-trackIn1.GetX())>0.1)||
2138 (TMath::Abs(trackOut0->GetAlpha()-trackOut1->GetAlpha())>0.00001)||
2139 (TMath::Abs(trackOut0->GetAlpha()-trackIn1.GetAlpha())>0.00001)||
2140 (TMath::Abs(trackIn0.GetTgl()-trackIn1.GetTgl())>0.0001)||
2141 (TMath::Abs(trackIn0.GetSnp()-trackIn1.GetSnp())>0.0001)
2148 corrections[icorr]=0;
2152 //if (ptype==4 &&bz<0) corrections[icorr]*=-1; // interpret as curvature - commented out no in lookup
2154 (*pcstream)<<"fit"<<
2155 Form("%s=",corr->GetName())<<corrections[icorr]; // dump correction value
2158 (*pcstream)<<"fit"<<"isOK="<<isOK<<"\n";
2167 void AliTPCCorrection::MakeSectorDistortionTree(TTree *tinput, Int_t dtype, Int_t ptype, const TObjArray * corrArray, Int_t step, Int_t offset, Bool_t debug ){
2170 // For each partial correction (specified in array) and given track topology (phi, theta, snp, refX)
2171 // calculates partial distortions
2172 // Partial distortion is stored in the resulting tree
2173 // Output is storred in the file distortion_<dettype>_<partype>.root
2174 // Partial distortion is stored with the name given by correction name
2177 // Parameters of function:
2178 // input - input tree
2179 // dtype - distortion type 10 - IROC-OROC
2180 // ppype - parameter type
2181 // corrArray - array with partial corrections
2182 // step - skipe entries - if 1 all entries processed - it is slow
2183 // debug 0 if debug on also space points dumped - it is slow
2185 const Double_t kMaxSnp = 0.8;
2186 const Int_t kMinEntries=200;
2187 // AliTPCROC *tpcRoc =AliTPCROC::Instance();
2189 const Double_t kMass = TDatabasePDG::Instance()->GetParticle("pi+")->Mass();
2190 // const Double_t kB2C=-0.299792458e-3;
2191 Double_t phi,theta, snp, mean,rms, entries,sector,dsec,globalZ;
2196 tinput->SetBranchAddress("run",&run);
2197 tinput->SetBranchAddress("theta",&theta);
2198 tinput->SetBranchAddress("phi", &phi);
2199 tinput->SetBranchAddress("snp",&snp);
2200 tinput->SetBranchAddress("mean",&mean);
2201 tinput->SetBranchAddress("rms",&rms);
2202 tinput->SetBranchAddress("entries",&entries);
2203 tinput->SetBranchAddress("sector",§or);
2204 tinput->SetBranchAddress("dsec",&dsec);
2205 tinput->SetBranchAddress("refX",&refXD);
2206 tinput->SetBranchAddress("z",&globalZ);
2207 tinput->SetBranchAddress("isec0",&isec0);
2208 tinput->SetBranchAddress("isec1",&isec1);
2209 TTreeSRedirector *pcstream = new TTreeSRedirector(Form("distortionSector%d_%d_%d.root",dtype,ptype,offset));
2211 Int_t nentries=tinput->GetEntries();
2212 Int_t ncorr=corrArray->GetEntries();
2213 Double_t corrections[100]={0}; //
2215 Double_t cov[15]={0,0,0,0,0,0,0,0,0,0,0,0,0,0};
2218 for (Int_t ientry=offset; ientry<nentries; ientry+=step){
2219 tinput->GetEntry(ientry);
2222 if (TMath::Abs(TMath::Abs(isec0%18)-TMath::Abs(isec1%18))==0) id=1; // IROC-OROC - opposite side
2223 if (TMath::Abs(TMath::Abs(isec0%36)-TMath::Abs(isec1%36))==0) id=2; // IROC-OROC - same side
2224 if (dtype==10 && id==-1) continue;
2231 tPar[4]=(gRandom->Rndm()-0.1)*0.2; //
2232 Double_t pt=1./tPar[4];
2234 printf("%f\t%f\t%f\t%f\t%f\t%f\n",entries, sector,theta,snp, mean,rms);
2235 Double_t bz=AliTrackerBase::GetBz();
2236 AliExternalTrackParam track(refX,phi,tPar,cov);
2237 Double_t xyz[3],xyzIn[3],xyzOut[3];
2239 track.GetXYZAt(85,bz,xyzIn);
2240 track.GetXYZAt(245,bz,xyzOut);
2241 Double_t phiIn = TMath::ATan2(xyzIn[1],xyzIn[0]);
2242 Double_t phiOut = TMath::ATan2(xyzOut[1],xyzOut[0]);
2243 Double_t phiRef = TMath::ATan2(xyz[1],xyz[0]);
2244 Int_t sectorRef = TMath::Nint(9.*phiRef/TMath::Pi()-0.5);
2245 Int_t sectorIn = TMath::Nint(9.*phiIn/TMath::Pi()-0.5);
2246 Int_t sectorOut = TMath::Nint(9.*phiOut/TMath::Pi()-0.5);
2249 if (sectorIn!=sectorOut) isOK=kFALSE; // requironment - cluster in the same sector
2250 if (sectorIn!=sectorRef) isOK=kFALSE; // requironment - cluster in the same sector
2251 if (entries<kMinEntries/(1+TMath::Abs(globalZ/100.))) isOK=kFALSE; // requironment - minimal amount of tracks in bin
2253 if (TMath::Abs(theta)>1) isOK=kFALSE;
2255 Double_t dRrec=0; // dummy value - needed for points - e.g for laser
2257 (*pcstream)<<"fit"<<
2259 "bz="<<bz<< // magnetic filed used
2260 "dtype="<<dtype<< // detector match type
2261 "ptype="<<ptype<< // parameter type
2262 "theta="<<theta<< // theta
2263 "phi="<<phi<< // phi
2264 "snp="<<snp<< // snp
2265 "mean="<<mean<< // mean dist value
2266 "rms="<<rms<< // rms
2269 "refX="<<refXD<< // referece X
2270 "gx="<<xyz[0]<< // global position at reference
2271 "gy="<<xyz[1]<< // global position at reference
2272 "gz="<<xyz[2]<< // global position at reference
2273 "dRrec="<<dRrec<< // delta Radius in reconstruction
2275 "id="<<id<< // track id
2276 "entries="<<entries;// number of entries in bin
2278 AliExternalTrackParam *trackOut0 = 0;
2279 AliExternalTrackParam *trackOut1 = 0;
2280 AliExternalTrackParam *ptrackIn0 = 0;
2281 AliExternalTrackParam *ptrackIn1 = 0;
2283 for (Int_t icorr=0; icorr<ncorr; icorr++) {
2285 // special case of the TPC tracks crossing the CE
2287 AliTPCCorrection *corr = (AliTPCCorrection*)corrArray->At(icorr);
2288 corrections[icorr]=0;
2289 if (entries>kMinEntries &&isOK){
2290 AliExternalTrackParam trackIn0(refX,phi,tPar,cov);
2291 AliExternalTrackParam trackIn1(refX,phi,tPar,cov);
2292 ptrackIn1=&trackIn0;
2293 ptrackIn0=&trackIn1;
2295 if (debug) trackOut0=corr->FitDistortedTrack(trackIn0, refX, dir,pcstream);
2296 if (!debug) trackOut0=corr->FitDistortedTrack(trackIn0, refX, dir,0);
2297 if (debug) trackOut1=corr->FitDistortedTrack(trackIn1, refX, -dir,pcstream);
2298 if (!debug) trackOut1=corr->FitDistortedTrack(trackIn1, refX, -dir,0);
2300 if (trackOut0 && trackOut1){
2302 if (!AliTrackerBase::PropagateTrackTo(&trackIn0,refX,kMass,1,kTRUE,kMaxSnp)) isOK=kFALSE;
2303 if (!AliTrackerBase::PropagateTrackTo(&trackIn0,refX,kMass,1,kFALSE,kMaxSnp)) isOK=kFALSE;
2304 // rotate all tracks to the same frame
2305 if (!trackOut0->Rotate(trackIn0.GetAlpha())) isOK=kFALSE;
2306 if (!trackIn1.Rotate(trackIn0.GetAlpha())) isOK=kFALSE;
2307 if (!trackOut1->Rotate(trackIn0.GetAlpha())) isOK=kFALSE;
2309 if (!AliTrackerBase::PropagateTrackTo(trackOut0,refX,kMass,1,kFALSE,kMaxSnp)) isOK=kFALSE;
2310 if (!AliTrackerBase::PropagateTrackTo(&trackIn1,refX,kMass,1,kFALSE,kMaxSnp)) isOK=kFALSE;
2311 if (!AliTrackerBase::PropagateTrackTo(trackOut1,refX,kMass,1,kFALSE,kMaxSnp)) isOK=kFALSE;
2313 corrections[icorr] = (trackOut0->GetParameter()[ptype]-trackIn0.GetParameter()[ptype]);
2314 corrections[icorr]-= (trackOut1->GetParameter()[ptype]-trackIn1.GetParameter()[ptype]);
2315 (*pcstream)<<"fitDebug"<< // just to debug the correction
2317 "pIn0.="<<ptrackIn0<<
2318 "pIn1.="<<ptrackIn1<<
2319 "pOut0.="<<trackOut0<<
2320 "pOut1.="<<trackOut1<<
2326 corrections[icorr]=0;
2330 (*pcstream)<<"fit"<<
2331 Form("%s=",corr->GetName())<<corrections[icorr]; // dump correction value
2334 (*pcstream)<<"fit"<<"isOK="<<isOK<<"\n";
2341 void AliTPCCorrection::MakeLaserDistortionTreeOld(TTree* tree, TObjArray *corrArray, Int_t itype){
2343 // Make a laser fit tree for global minimization
2345 const Double_t cutErrY=0.1;
2346 const Double_t cutErrZ=0.1;
2347 const Double_t kEpsilon=0.00000001;
2348 const Double_t kMaxDist=1.; // max distance - space correction
2349 const Double_t kMaxRMS=0.05; // max distance -between point and local mean
2354 AliTPCLaserTrack *ltr=0;
2355 AliTPCLaserTrack::LoadTracks();
2356 tree->SetBranchAddress("dY.",&vecdY);
2357 tree->SetBranchAddress("dZ.",&vecdZ);
2358 tree->SetBranchAddress("eY.",&veceY);
2359 tree->SetBranchAddress("eZ.",&veceZ);
2360 tree->SetBranchAddress("LTr.",<r);
2361 Int_t entries= tree->GetEntries();
2362 TTreeSRedirector *pcstream= new TTreeSRedirector("distortionLaser_0.root");
2363 Double_t bz=AliTrackerBase::GetBz();
2366 for (Int_t ientry=0; ientry<entries; ientry++){
2367 tree->GetEntry(ientry);
2368 if (!ltr->GetVecGX()){
2369 ltr->UpdatePoints();
2371 TVectorD * delta= (itype==0)? vecdY:vecdZ;
2372 TVectorD * err= (itype==0)? veceY:veceZ;
2373 TLinearFitter fitter(2,"pol1");
2374 for (Int_t iter=0; iter<2; iter++){
2375 Double_t kfit0=0, kfit1=0;
2376 Int_t npoints=fitter.GetNpoints();
2379 kfit0=fitter.GetParameter(0);
2380 kfit1=fitter.GetParameter(1);
2382 for (Int_t irow=0; irow<159; irow++){
2385 Int_t nentries = 1000;
2386 if (veceY->GetMatrixArray()[irow]>cutErrY||veceZ->GetMatrixArray()[irow]>cutErrZ) nentries=0;
2387 if (veceY->GetMatrixArray()[irow]<kEpsilon||veceZ->GetMatrixArray()[irow]<kEpsilon) nentries=0;
2390 Int_t first3=TMath::Max(irow-3,0);
2391 Int_t last3 =TMath::Min(irow+3,159);
2393 if ((*ltr->GetVecSec())[irow]>=0 && err) {
2394 for (Int_t jrow=first3; jrow<=last3; jrow++){
2395 if ((*ltr->GetVecSec())[irow]!= (*ltr->GetVecSec())[jrow]) continue;
2396 if ((*err)[jrow]<kEpsilon) continue;
2397 array[counter]=(*delta)[jrow];
2404 rms3 = TMath::RMS(counter,array);
2405 mean3 = TMath::Mean(counter,array);
2409 Double_t phi =(*ltr->GetVecPhi())[irow];
2410 Double_t theta =ltr->GetTgl();
2411 Double_t mean=delta->GetMatrixArray()[irow];
2412 Double_t gx=0,gy=0,gz=0;
2413 Double_t snp = (*ltr->GetVecP2())[irow];
2415 // Double_t rms = err->GetMatrixArray()[irow];
2417 gx = (*ltr->GetVecGX())[irow];
2418 gy = (*ltr->GetVecGY())[irow];
2419 gz = (*ltr->GetVecGZ())[irow];
2421 // get delta R used in reconstruction
2422 AliTPCcalibDB* calib=AliTPCcalibDB::Instance();
2423 AliTPCCorrection * correction = calib->GetTPCComposedCorrection(AliTrackerBase::GetBz());
2424 // const AliTPCRecoParam * recoParam = calib->GetTransform()->GetCurrentRecoParam();
2425 //Double_t xyz0[3]={gx,gy,gz};
2426 Double_t oldR=TMath::Sqrt(gx*gx+gy*gy);
2427 Double_t fphi = TMath::ATan2(gy,gx);
2428 Double_t fsector = 9.*fphi/TMath::Pi();
2429 if (fsector<0) fsector+=18;
2430 Double_t dsec = fsector-Int_t(fsector)-0.5;
2432 Int_t id= ltr->GetId();
2436 Float_t xyz1[3]={gx,gy,gz};
2437 Int_t sector=(gz>0)?0:18;
2438 correction->CorrectPoint(xyz1, sector);
2439 refX=TMath::Sqrt(xyz1[0]*xyz1[0]+xyz1[1]*xyz1[1]);
2442 if (TMath::Abs(rms3)>kMaxRMS) isOK=kFALSE;
2443 if (TMath::Abs(mean-mean3)>kMaxRMS) isOK=kFALSE;
2444 if (counter<4) isOK=kFALSE;
2445 if (npoints<90) isOK=kFALSE;
2447 fitter.AddPoint(&refX,mean);
2449 Double_t deltaF=kfit0+kfit1*refX;
2451 (*pcstream)<<"fitFull"<< // dumpe also intermediate results
2452 "bz="<<bz<< // magnetic filed used
2453 "dtype="<<dtype<< // detector match type
2454 "ptype="<<itype<< // parameter type
2455 "theta="<<theta<< // theta
2456 "phi="<<phi<< // phi
2457 "snp="<<snp<< // snp
2458 "mean="<<mean3<< // mean dist value
2459 "rms="<<rms3<< // rms
2461 "npoints="<<npoints<< //number of points
2462 "mean3="<<mean3<< // mean dist value
2463 "rms3="<<rms3<< // rms
2464 "counter="<<counter<<
2465 "sector="<<fsector<<
2468 "refX="<<refX<< // reference radius
2469 "gx="<<gx<< // global position
2470 "gy="<<gy<< // global position
2471 "gz="<<gz<< // global position
2472 "dRrec="<<dRrec<< // delta Radius in reconstruction
2473 "id="<<id<< //bundle
2474 "entries="<<nentries<<// number of entries in bin
2477 if (iter==1) (*pcstream)<<"fit"<< // dump valus for fit
2478 "bz="<<bz<< // magnetic filed used
2479 "dtype="<<dtype<< // detector match type
2480 "ptype="<<itype<< // parameter type
2481 "theta="<<theta<< // theta
2482 "phi="<<phi<< // phi
2483 "snp="<<snp<< // snp
2484 "mean="<<mean3<< // mean dist value
2485 "rms="<<rms3<< // rms
2486 "sector="<<fsector<<
2489 "refX="<<refX<< // reference radius
2490 "gx="<<gx<< // global position
2491 "gy="<<gy<< // global position
2492 "gz="<<gz<< // global position
2493 "dRrec="<<dRrec<< // delta Radius in reconstruction
2495 "id="<<id<< //bundle
2496 "entries="<<nentries;// number of entries in bin
2499 Double_t ky = TMath::Tan(TMath::ASin(snp));
2500 Int_t ncorr = corrArray->GetEntries();
2501 Double_t r0 = TMath::Sqrt(gx*gx+gy*gy);
2502 Double_t phi0 = TMath::ATan2(gy,gx);
2503 Double_t distortions[1000]={0};
2504 Double_t distortionsR[1000]={0};
2506 for (Int_t icorr=0; icorr<ncorr; icorr++) {
2507 AliTPCCorrection *corr = (AliTPCCorrection*)corrArray->At(icorr);
2508 Float_t distPoint[3]={gx,gy,gz};
2509 Int_t sector= (gz>0)? 0:18;
2511 corr->DistortPoint(distPoint, sector);
2513 // Double_t value=distPoint[2]-gz;
2514 if (itype==0 && r0>1){
2515 Double_t r1 = TMath::Sqrt(distPoint[0]*distPoint[0]+distPoint[1]*distPoint[1]);
2516 Double_t phi1 = TMath::ATan2(distPoint[1],distPoint[0]);
2517 Double_t drphi= r0*(phi1-phi0);
2518 Double_t dr = r1-r0;
2519 distortions[icorr] = drphi-ky*dr;
2520 distortionsR[icorr] = dr;
2522 if (TMath::Abs(distortions[icorr])>kMaxDist) {isOKF=icorr+1; isOK=kFALSE; }
2523 if (TMath::Abs(distortionsR[icorr])>kMaxDist) {isOKF=icorr+1; isOK=kFALSE;}
2524 (*pcstream)<<"fit"<<
2525 Form("%s=",corr->GetName())<<distortions[icorr]; // dump correction value
2527 (*pcstream)<<"fit"<<"isOK="<<isOK<<"\n";
2537 void AliTPCCorrection::MakeDistortionMap(THnSparse * his0, TTreeSRedirector * const pcstream, const char* hname, Int_t run, Float_t refX, Int_t type, Int_t integ){
2539 // make a distortion map out ou fthe residual histogram
2540 // Results are written to the debug streamer - pcstream
2542 // his0 - input (4D) residual histogram
2543 // pcstream - file to write the tree
2545 // refX - track matching reference X
2546 // type - 0- y 1-z,2 -snp, 3-theta, 4=1/pt
2548 // OBJ: TAxis #Delta #Delta
2549 // OBJ: TAxis tanTheta tan(#Theta)
2550 // OBJ: TAxis phi #phi
2551 // OBJ: TAxis snp snp
2553 // marian.ivanov@cern.ch
2554 const Int_t kMinEntries=10;
2555 Double_t bz=AliTrackerBase::GetBz();
2556 Int_t idim[4]={0,1,2,3};
2560 Int_t nbins3=his0->GetAxis(3)->GetNbins();
2561 Int_t first3=his0->GetAxis(3)->GetFirst();
2562 Int_t last3 =his0->GetAxis(3)->GetLast();
2564 for (Int_t ibin3=first3; ibin3<last3; ibin3+=1){ // axis 3 - local angle
2565 his0->GetAxis(3)->SetRange(TMath::Max(ibin3-integ,1),TMath::Min(ibin3+integ,nbins3));
2566 Double_t x3= his0->GetAxis(3)->GetBinCenter(ibin3);
2567 THnSparse * his3= his0->Projection(3,idim); //projected histogram according selection 3
2569 Int_t nbins2 = his3->GetAxis(2)->GetNbins();
2570 Int_t first2 = his3->GetAxis(2)->GetFirst();
2571 Int_t last2 = his3->GetAxis(2)->GetLast();
2573 for (Int_t ibin2=first2; ibin2<last2; ibin2+=1){ // axis 2 - phi
2574 his3->GetAxis(2)->SetRange(TMath::Max(ibin2-integ,1),TMath::Min(ibin2+integ,nbins2));
2575 Double_t x2= his3->GetAxis(2)->GetBinCenter(ibin2);
2576 THnSparse * his2= his3->Projection(2,idim); //projected histogram according selection 2
2577 Int_t nbins1 = his2->GetAxis(1)->GetNbins();
2578 Int_t first1 = his2->GetAxis(1)->GetFirst();
2579 Int_t last1 = his2->GetAxis(1)->GetLast();
2580 for (Int_t ibin1=first1; ibin1<last1; ibin1++){ //axis 1 - theta
2582 Double_t x1= his2->GetAxis(1)->GetBinCenter(ibin1);
2583 his2->GetAxis(1)->SetRange(TMath::Max(ibin1-1,1),TMath::Min(ibin1+1,nbins1));
2584 if (TMath::Abs(x1)<0.1){
2585 if (x1<0) his2->GetAxis(1)->SetRange(TMath::Max(ibin1-1,1),TMath::Min(ibin1,nbins1));
2586 if (x1>0) his2->GetAxis(1)->SetRange(TMath::Max(ibin1,1),TMath::Min(ibin1+1,nbins1));
2588 if (TMath::Abs(x1)<0.06){
2589 his2->GetAxis(1)->SetRange(TMath::Max(ibin1,1),TMath::Min(ibin1,nbins1));
2591 TH1 * hisDelta = his2->Projection(0);
2593 Double_t entries = hisDelta->GetEntries();
2594 Double_t mean=0, rms=0;
2595 if (entries>kMinEntries){
2596 mean = hisDelta->GetMean();
2597 rms = hisDelta->GetRMS();
2599 Double_t sector = 9.*x2/TMath::Pi();
2600 if (sector<0) sector+=18;
2601 Double_t dsec = sector-Int_t(sector)-0.5;
2603 (*pcstream)<<hname<<
2608 "z="<<z<< // dummy z
2610 "entries="<<entries<<
2613 "refX="<<refX<< // track matching refernce plane
2619 //printf("%f\t%f\t%f\t%f\t%f\n",x3,x2,x1, entries,mean);
2630 void AliTPCCorrection::MakeDistortionMapCosmic(THnSparse * hisInput, TTreeSRedirector * const pcstream, const char* hname, Int_t run, Float_t refX, Int_t type){
2632 // make a distortion map out ou fthe residual histogram
2633 // Results are written to the debug streamer - pcstream
2635 // his0 - input (4D) residual histogram
2636 // pcstream - file to write the tree
2638 // refX - track matching reference X
2639 // type - 0- y 1-z,2 -snp, 3-theta, 4=1/pt
2640 // marian.ivanov@cern.ch
2643 // Collection name='TObjArray', class='TObjArray', size=16
2644 // 0. OBJ: TAxis #Delta #Delta
2645 // 1. OBJ: TAxis N_{cl} N_{cl}
2646 // 2. OBJ: TAxis dca_{r} (cm) dca_{r} (cm)
2647 // 3. OBJ: TAxis z (cm) z (cm)
2648 // 4. OBJ: TAxis sin(#phi) sin(#phi)
2649 // 5. OBJ: TAxis tan(#theta) tan(#theta)
2650 // 6. OBJ: TAxis 1/pt (1/GeV) 1/pt (1/GeV)
2651 // 7. OBJ: TAxis pt (GeV) pt (GeV)
2652 // 8. OBJ: TAxis alpha alpha
2653 const Int_t kMinEntries=10;
2655 // 1. make default selections
2658 Int_t idim0[4]={0 , 5, 8, 3}; // delta, theta, alpha, z
2659 hisInput->GetAxis(1)->SetRangeUser(110,190); //long tracks
2660 hisInput->GetAxis(2)->SetRangeUser(-10,35); //tracks close to beam pipe
2661 hisInput->GetAxis(4)->SetRangeUser(-0.3,0.3); //small snp at TPC entrance
2662 hisInput->GetAxis(7)->SetRangeUser(3,100); //"high pt tracks"
2663 hisDelta= hisInput->Projection(0);
2664 hisInput->GetAxis(0)->SetRangeUser(-6.*hisDelta->GetRMS(), +6.*hisDelta->GetRMS());
2666 THnSparse *his0= hisInput->Projection(4,idim0);
2668 // 2. Get mean in diferent bins
2670 Int_t nbins1=his0->GetAxis(1)->GetNbins();
2671 Int_t first1=his0->GetAxis(1)->GetFirst();
2672 Int_t last1 =his0->GetAxis(1)->GetLast();
2674 Double_t bz=AliTrackerBase::GetBz();
2675 Int_t idim[4]={0,1, 2, 3}; // delta, theta,alpha,z
2677 for (Int_t ibin1=first1; ibin1<=last1; ibin1++){ //axis 1 - theta
2679 Double_t x1= his0->GetAxis(1)->GetBinCenter(ibin1);
2680 his0->GetAxis(1)->SetRange(TMath::Max(ibin1-1,1),TMath::Min(ibin1+1,nbins1));
2682 THnSparse * his1 = his0->Projection(4,idim); // projected histogram according range1
2683 Int_t nbins3 = his1->GetAxis(3)->GetNbins();
2684 Int_t first3 = his1->GetAxis(3)->GetFirst();
2685 Int_t last3 = his1->GetAxis(3)->GetLast();
2687 for (Int_t ibin3=first3-1; ibin3<=last3; ibin3+=1){ // axis 3 - z at "vertex"
2688 his1->GetAxis(3)->SetRange(TMath::Max(ibin3-1,1),TMath::Min(ibin3+1,nbins3));
2689 Double_t x3= his1->GetAxis(3)->GetBinCenter(ibin3);
2691 his1->GetAxis(3)->SetRangeUser(-1,1);
2694 THnSparse * his3= his1->Projection(4,idim); //projected histogram according selection 3
2695 Int_t nbins2 = his3->GetAxis(2)->GetNbins();
2696 Int_t first2 = his3->GetAxis(2)->GetFirst();
2697 Int_t last2 = his3->GetAxis(2)->GetLast();
2699 for (Int_t ibin2=first2; ibin2<=last2; ibin2+=1){
2700 his3->GetAxis(2)->SetRange(TMath::Max(ibin2-1,1),TMath::Min(ibin2+1,nbins2));
2701 Double_t x2= his3->GetAxis(2)->GetBinCenter(ibin2);
2702 hisDelta = his3->Projection(0);
2704 Double_t entries = hisDelta->GetEntries();
2705 Double_t mean=0, rms=0;
2706 if (entries>kMinEntries){
2707 mean = hisDelta->GetMean();
2708 rms = hisDelta->GetRMS();
2710 Double_t sector = 9.*x2/TMath::Pi();
2711 if (sector<0) sector+=18;
2712 Double_t dsec = sector-Int_t(sector)-0.5;
2713 Double_t snp=0; // dummy snp - equal 0
2714 (*pcstream)<<hname<<
2716 "bz="<<bz<< // magnetic field
2717 "theta="<<x1<< // theta
2718 "phi="<<x2<< // phi (alpha)
2719 "z="<<x3<< // z at "vertex"
2720 "snp="<<snp<< // dummy snp
2721 "entries="<<entries<< // entries in bin
2722 "mean="<<mean<< // mean
2724 "refX="<<refX<< // track matching refernce plane
2725 "type="<<type<< // parameter type
2726 "sector="<<sector<< // sector
2727 "dsec="<<dsec<< // dummy delta sector
2730 printf("%f\t%f\t%f\t%f\t%f\n",x1,x3,x2, entries,mean);
2741 void AliTPCCorrection::MakeDistortionMapSector(THnSparse * hisInput, TTreeSRedirector * const pcstream, const char* hname, Int_t run, Int_t type){
2743 // make a distortion map out of the residual histogram
2744 // Results are written to the debug streamer - pcstream
2746 // his0 - input (4D) residual histogram
2747 // pcstream - file to write the tree
2749 // type - 0- y 1-z,2 -snp, 3-theta
2750 // marian.ivanov@cern.ch
2752 //Collection name='TObjArray', class='TObjArray', size=16
2753 //0 OBJ: TAxis delta delta
2754 //1 OBJ: TAxis phi phi
2755 //2 OBJ: TAxis localX localX
2756 //3 OBJ: TAxis kY kY
2757 //4 OBJ: TAxis kZ kZ
2758 //5 OBJ: TAxis is1 is1
2759 //6 OBJ: TAxis is0 is0
2761 //8. OBJ: TAxis IsPrimary IsPrimary
2763 const Int_t kMinEntries=10;
2764 THnSparse * hisSector0=0;
2765 TH1 * htemp=0; // histogram to calculate mean value of parameter
2766 Double_t bz=AliTrackerBase::GetBz();
2769 // Loop over pair of sector:
2780 for (Int_t isec0=0; isec0<72; isec0++){
2781 Int_t index0[9]={0, 4, 3, 7, 1, 2, 5, 6,8}; //regroup indeces
2783 //hisInput->GetAxis(8)->SetRangeUser(-0.1,0.4); // select secondaries only ? - get out later ?
2784 hisInput->GetAxis(6)->SetRangeUser(isec0-0.1,isec0+0.1);
2785 hisSector0=hisInput->Projection(7,index0);
2788 for (Int_t isec1=isec0+1; isec1<72; isec1++){
2789 //if (isec1!=isec0+36) continue;
2790 if ( TMath::Abs((isec0%18)-(isec1%18))>1.5 && TMath::Abs((isec0%18)-(isec1%18))<16.5) continue;
2791 printf("Sectors %d\t%d\n",isec1,isec0);
2792 hisSector0->GetAxis(6)->SetRangeUser(isec1-0.1,isec1+0.1);
2793 TH1 * hisX=hisSector0->Projection(5);
2794 Double_t refX= hisX->GetMean();
2796 TH1 *hisDelta=hisSector0->Projection(0);
2797 Double_t dmean = hisDelta->GetMean();
2798 Double_t drms = hisDelta->GetRMS();
2799 hisSector0->GetAxis(0)->SetRangeUser(dmean-5.*drms, dmean+5.*drms);
2802 // 1. make default selections
2804 Int_t idim0[5]={0 , 1, 2, 3, 4}; // {delta, theta, snp, z, phi }
2805 THnSparse *hisSector1= hisSector0->Projection(5,idim0);
2807 // 2. Get mean in diferent bins
2809 Int_t idim[5]={0, 1, 2, 3, 4}; // {delta, theta-1,snp-2 ,z-3, phi-4}
2811 // Int_t nbinsPhi=hisSector1->GetAxis(4)->GetNbins();
2812 Int_t firstPhi=hisSector1->GetAxis(4)->GetFirst();
2813 Int_t lastPhi =hisSector1->GetAxis(4)->GetLast();
2815 for (Int_t ibinPhi=firstPhi; ibinPhi<=lastPhi; ibinPhi+=1){ //axis 4 - phi
2819 Double_t xPhi= hisSector1->GetAxis(4)->GetBinCenter(ibinPhi);
2820 Double_t psec = (9*xPhi/TMath::Pi());
2821 if (psec<0) psec+=18;
2822 Bool_t isOK0=kFALSE;
2823 Bool_t isOK1=kFALSE;
2824 if (TMath::Abs(psec-isec0%18-0.5)<1. || TMath::Abs(psec-isec0%18-17.5)<1.) isOK0=kTRUE;
2825 if (TMath::Abs(psec-isec1%18-0.5)<1. || TMath::Abs(psec-isec1%18-17.5)<1.) isOK1=kTRUE;
2826 if (!isOK0) continue;
2827 if (!isOK1) continue;
2829 hisSector1->GetAxis(4)->SetRange(TMath::Max(ibinPhi-2,firstPhi),TMath::Min(ibinPhi+2,lastPhi));
2830 if (isec1!=isec0+36) {
2831 hisSector1->GetAxis(4)->SetRange(TMath::Max(ibinPhi-3,firstPhi),TMath::Min(ibinPhi+3,lastPhi));
2834 htemp = hisSector1->Projection(4);
2835 xPhi=htemp->GetMean();
2837 THnSparse * hisPhi = hisSector1->Projection(4,idim);
2838 //Int_t nbinsZ = hisPhi->GetAxis(3)->GetNbins();
2839 Int_t firstZ = hisPhi->GetAxis(3)->GetFirst();
2840 Int_t lastZ = hisPhi->GetAxis(3)->GetLast();
2842 for (Int_t ibinZ=firstZ; ibinZ<=lastZ; ibinZ+=1){ // axis 3 - z
2846 hisPhi->GetAxis(3)->SetRange(TMath::Max(ibinZ,firstZ),TMath::Min(ibinZ,lastZ));
2847 if (isec1!=isec0+36) {
2848 hisPhi->GetAxis(3)->SetRange(TMath::Max(ibinZ-1,firstZ),TMath::Min(ibinZ-1,lastZ));
2850 htemp = hisPhi->Projection(3);
2851 Double_t xZ= htemp->GetMean();
2853 THnSparse * hisZ= hisPhi->Projection(3,idim);
2854 //projected histogram according selection 3 -z
2857 //Int_t nbinsSnp = hisZ->GetAxis(2)->GetNbins();
2858 Int_t firstSnp = hisZ->GetAxis(2)->GetFirst();
2859 Int_t lastSnp = hisZ->GetAxis(2)->GetLast();
2860 for (Int_t ibinSnp=firstSnp; ibinSnp<=lastSnp; ibinSnp+=2){ // axis 2 - snp
2864 hisZ->GetAxis(2)->SetRange(TMath::Max(ibinSnp-1,firstSnp),TMath::Min(ibinSnp+1,lastSnp));
2865 if (isec1!=isec0+36) {
2866 hisZ->GetAxis(2)->SetRange(TMath::Max(ibinSnp-2,firstSnp),TMath::Min(ibinSnp+2,lastSnp));
2868 htemp = hisZ->Projection(2);
2869 Double_t xSnp= htemp->GetMean();
2871 THnSparse * hisSnp= hisZ->Projection(2,idim);
2872 //projected histogram according selection 2 - snp
2874 //Int_t nbinsTheta = hisSnp->GetAxis(1)->GetNbins();
2875 Int_t firstTheta = hisSnp->GetAxis(1)->GetFirst();
2876 Int_t lastTheta = hisSnp->GetAxis(1)->GetLast();
2878 for (Int_t ibinTheta=firstTheta; ibinTheta<=lastTheta; ibinTheta+=2){ // axis1 theta
2881 hisSnp->GetAxis(1)->SetRange(TMath::Max(ibinTheta-2,firstTheta),TMath::Min(ibinTheta+2,lastTheta));
2882 if (isec1!=isec0+36) {
2883 hisSnp->GetAxis(1)->SetRange(TMath::Max(ibinTheta-3,firstTheta),TMath::Min(ibinTheta+3,lastTheta));
2885 htemp = hisSnp->Projection(1);
2886 Double_t xTheta=htemp->GetMean();
2888 hisDelta = hisSnp->Projection(0);
2890 Double_t entries = hisDelta->GetEntries();
2891 Double_t mean=0, rms=0;
2892 if (entries>kMinEntries){
2893 mean = hisDelta->GetMean();
2894 rms = hisDelta->GetRMS();
2896 Double_t sector = 9.*xPhi/TMath::Pi();
2897 if (sector<0) sector+=18;
2898 Double_t dsec = sector-Int_t(sector)-0.5;
2899 Int_t dtype=1; // TPC alignment type
2900 (*pcstream)<<hname<<
2902 "bz="<<bz<< // magnetic field
2903 "ptype="<<type<< // parameter type
2904 "dtype="<<dtype<< // parameter type
2905 "isec0="<<isec0<< // sector 0
2906 "isec1="<<isec1<< // sector 1
2907 "sector="<<sector<< // sector as float
2908 "dsec="<<dsec<< // delta sector
2910 "theta="<<xTheta<< // theta
2911 "phi="<<xPhi<< // phi (alpha)
2913 "snp="<<xSnp<< // snp
2915 "entries="<<entries<< // entries in bin
2916 "mean="<<mean<< // mean
2917 "rms="<<rms<< // rms
2918 "refX="<<refX<< // track matching reference plane
2921 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);
2942 void AliTPCCorrection::StoreInOCDB(Int_t startRun, Int_t endRun, const char *comment){
2944 // Store object in the OCDB
2945 // By default the object is stored in the current directory
2946 // default comment consit of user name and the date
2948 TString ocdbStorage="";
2949 ocdbStorage+="local://"+gSystem->GetFromPipe("pwd")+"/OCDB";
2950 AliCDBMetaData *metaData= new AliCDBMetaData();
2951 metaData->SetObjectClassName("AliTPCCorrection");
2952 metaData->SetResponsible("Marian Ivanov");
2953 metaData->SetBeamPeriod(1);
2954 metaData->SetAliRootVersion("05-25-01"); //root version
2955 TString userName=gSystem->GetFromPipe("echo $USER");
2956 TString date=gSystem->GetFromPipe("date");
2958 if (!comment) metaData->SetComment(Form("Space point distortion calibration\n User: %s\n Data%s",userName.Data(),date.Data()));
2959 if (comment) metaData->SetComment(comment);
2961 id1=new AliCDBId("TPC/Calib/Correction", startRun, endRun);
2962 AliCDBStorage* gStorage = AliCDBManager::Instance()->GetStorage(ocdbStorage);
2963 gStorage->Put(this, (*id1), metaData);
2967 void AliTPCCorrection::FastSimDistortedVertex(Double_t orgVertex[3], Int_t nTracks, AliESDVertex &aV, AliESDVertex &avOrg, AliESDVertex &cV, AliESDVertex &cvOrg, TTreeSRedirector * const pcstream, Double_t etaCuts){
2969 // Fast method to simulate the influence of the given distortion on the vertex reconstruction
2972 AliMagF* magF= (AliMagF*)TGeoGlobalMagField::Instance()->GetField();
2973 if (!magF) AliError("Magneticd field - not initialized");
2974 Double_t bz = magF->SolenoidField(); //field in kGauss
2975 printf("bz: %f\n",bz);
2976 AliVertexerTracks *vertexer = new AliVertexerTracks(bz); // bz in kGauss
2978 TObjArray aTrk; // Original Track array of Aside
2979 TObjArray daTrk; // Distorted Track array of A side
2980 UShort_t *aId = new UShort_t[nTracks]; // A side Track ID
2983 UShort_t *cId = new UShort_t [nTracks];
2985 Double_t mass = TDatabasePDG::Instance()->GetParticle("pi+")->Mass();
2986 TF1 fpt("fpt",Form("x*(1+(sqrt(x*x+%f^2)-%f)/([0]*[1]))^(-[0])",mass,mass),0.4,10);
2987 fpt.SetParameters(7.24,0.120);
2989 for(Int_t nt=0; nt<nTracks; nt++){
2990 Double_t phi = gRandom->Uniform(0.0, 2*TMath::Pi());
2991 Double_t eta = gRandom->Uniform(-etaCuts, etaCuts);
2992 Double_t pt = fpt.GetRandom(); // momentum for f1
2993 // printf("phi %lf eta %lf pt %lf\n",phi,eta,pt);
2995 if(gRandom->Rndm() < 0.5){
3001 Double_t theta = 2*TMath::ATan(TMath::Exp(-eta))-TMath::Pi()/2.;
3003 pxyz[0]=pt*TMath::Cos(phi);
3004 pxyz[1]=pt*TMath::Sin(phi);
3005 pxyz[2]=pt*TMath::Tan(theta);
3006 Double_t cv[21]={0};
3007 AliExternalTrackParam *t= new AliExternalTrackParam(orgVertex, pxyz, cv, sign);
3011 AliExternalTrackParam *td = FitDistortedTrack(*t, refX, dir, NULL);
3013 if (pcstream) (*pcstream)<<"track"<<
3019 if(( eta>0.07 )&&( eta<etaCuts )) { // - log(tan(0.5*theta)), theta = 0.5*pi - ATan(5.0/80.0)
3023 Int_t nn=aTrk.GetEntriesFast();
3026 }else if(( eta<-0.07 )&&( eta>-etaCuts )){
3030 Int_t nn=cTrk.GetEntriesFast();
3035 }// end of track loop
3037 vertexer->SetTPCMode();
3038 vertexer->SetConstraintOff();
3040 aV = *((AliESDVertex*)vertexer->FindPrimaryVertex(&daTrk,aId));
3041 avOrg = *((AliESDVertex*)vertexer->FindPrimaryVertex(&aTrk,aId));
3042 cV = *((AliESDVertex*)vertexer->FindPrimaryVertex(&dcTrk,cId));
3043 cvOrg = *((AliESDVertex*)vertexer->FindPrimaryVertex(&cTrk,cId));
3044 if (pcstream) (*pcstream)<<"vertex"<<
3045 "x="<<orgVertex[0]<<
3046 "y="<<orgVertex[1]<<
3047 "z="<<orgVertex[2]<<
3048 "av.="<<&aV<< // distorted vertex A side
3049 "cv.="<<&cV<< // distroted vertex C side
3050 "avO.="<<&avOrg<< // original vertex A side
3057 void AliTPCCorrection::AddVisualCorrection(AliTPCCorrection* corr, Int_t position){
3059 // make correction available for visualization using
3060 // TFormula, TFX and TTree::Draw
3061 // important in order to check corrections and also compute dervied variables
3062 // e.g correction partial derivatives
3064 // NOTE - class is not owner of correction
3066 if (!fgVisualCorrection) fgVisualCorrection=new TObjArray(10000);
3067 if (position>=fgVisualCorrection->GetEntriesFast())
3068 fgVisualCorrection->Expand((position+10)*2);
3069 fgVisualCorrection->AddAt(corr, position);
3074 Double_t AliTPCCorrection::GetCorrSector(Double_t sector, Double_t r, Double_t kZ, Int_t axisType, Int_t corrType){
3076 // calculate the correction at given position - check the geffCorr
3078 // corrType return values
3084 if (!fgVisualCorrection) return 0;
3085 AliTPCCorrection *corr = (AliTPCCorrection*)fgVisualCorrection->At(corrType);
3086 if (!corr) return 0;
3088 Double_t phi=sector*TMath::Pi()/9.;
3089 Double_t gx = r*TMath::Cos(phi);
3090 Double_t gy = r*TMath::Sin(phi);
3092 Int_t nsector=(gz>=0) ? 0:18;
3096 Float_t distPoint[3]={gx,gy,gz};
3097 corr->DistortPoint(distPoint, nsector);
3098 Double_t r0=TMath::Sqrt(gx*gx+gy*gy);
3099 Double_t r1=TMath::Sqrt(distPoint[0]*distPoint[0]+distPoint[1]*distPoint[1]);
3100 Double_t phi0=TMath::ATan2(gy,gx);
3101 Double_t phi1=TMath::ATan2(distPoint[1],distPoint[0]);
3102 if (axisType==0) return r1-r0;
3103 if (axisType==1) return (phi1-phi0)*r0;
3104 if (axisType==2) return distPoint[2]-gz;
3105 if (axisType==3) return (TMath::Cos(phi)*(distPoint[0]-gx)+ TMath::Cos(phi)*(distPoint[1]-gy));
3109 Double_t AliTPCCorrection::GetCorrXYZ(Double_t gx, Double_t gy, Double_t gz, Int_t axisType, Int_t corrType){
3111 // return correction at given x,y,z
3113 if (!fgVisualCorrection) return 0;
3114 AliTPCCorrection *corr = (AliTPCCorrection*)fgVisualCorrection->At(corrType);
3115 if (!corr) return 0;
3116 Double_t phi0= TMath::ATan2(gy,gx);
3117 Int_t nsector=(gz>=0) ? 0:18;
3118 Float_t distPoint[3]={gx,gy,gz};
3119 corr->CorrectPoint(distPoint, nsector);
3120 Double_t r0=TMath::Sqrt(gx*gx+gy*gy);
3121 Double_t r1=TMath::Sqrt(distPoint[0]*distPoint[0]+distPoint[1]*distPoint[1]);
3122 Double_t phi1=TMath::ATan2(distPoint[1],distPoint[0]);
3123 if (axisType==0) return r1-r0;
3124 if (axisType==1) return (phi1-phi0)*r0;
3125 if (axisType==2) return distPoint[2]-gz;
3129 Double_t AliTPCCorrection::GetCorrXYZDz(Double_t gx, Double_t gy, Double_t gz, Int_t axisType, Int_t corrType,Double_t delta){
3131 // return correction at given x,y,z
3133 if (!fgVisualCorrection) return 0;
3134 AliTPCCorrection *corr = (AliTPCCorrection*)fgVisualCorrection->At(corrType);
3135 if (!corr) return 0;
3136 Double_t phi0= TMath::ATan2(gy,gx);
3137 Int_t nsector=(gz>=0) ? 0:18;
3138 Float_t distPoint[3]={gx,gy,gz};
3139 Float_t dxyz[3]={gx,gy,gz};
3141 corr->GetCorrectionDz(distPoint, nsector,dxyz,delta);
3142 distPoint[0]+=dxyz[0];
3143 distPoint[1]+=dxyz[1];
3144 distPoint[2]+=dxyz[2];
3145 Double_t r0=TMath::Sqrt(gx*gx+gy*gy);
3146 Double_t r1=TMath::Sqrt(distPoint[0]*distPoint[0]+distPoint[1]*distPoint[1]);
3147 Double_t phi1=TMath::ATan2(distPoint[1],distPoint[0]);
3148 if (axisType==0) return r1-r0;
3149 if (axisType==1) return (phi1-phi0)*r0;
3150 if (axisType==2) return distPoint[2]-gz;
3154 Double_t AliTPCCorrection::GetCorrXYZIntegrateZ(Double_t gx, Double_t gy, Double_t gz, Int_t axisType, Int_t corrType,Double_t delta){
3156 // return correction at given x,y,z
3158 if (!fgVisualCorrection) return 0;
3159 AliTPCCorrection *corr = (AliTPCCorrection*)fgVisualCorrection->At(corrType);
3160 if (!corr) return 0;
3161 Double_t phi0= TMath::ATan2(gy,gx);
3162 Int_t nsector=(gz>=0) ? 0:18;
3163 Float_t distPoint[3]={gx,gy,gz};
3164 Float_t dxyz[3]={gx,gy,gz};
3166 corr->GetCorrectionIntegralDz(distPoint, nsector,dxyz,delta);
3167 distPoint[0]+=dxyz[0];
3168 distPoint[1]+=dxyz[1];
3169 distPoint[2]+=dxyz[2];
3170 Double_t r0=TMath::Sqrt(gx*gx+gy*gy);
3171 Double_t r1=TMath::Sqrt(distPoint[0]*distPoint[0]+distPoint[1]*distPoint[1]);
3172 Double_t phi1=TMath::ATan2(distPoint[1],distPoint[0]);
3173 if (axisType==0) return r1-r0;
3174 if (axisType==1) return (phi1-phi0)*r0;
3175 if (axisType==2) return distPoint[2]-gz;
3180 Double_t AliTPCCorrection::GetDistXYZ(Double_t gx, Double_t gy, Double_t gz, Int_t axisType, Int_t corrType){
3182 // return correction at given x,y,z
3184 if (!fgVisualCorrection) return 0;
3185 AliTPCCorrection *corr = (AliTPCCorrection*)fgVisualCorrection->At(corrType);
3186 if (!corr) return 0;
3187 Double_t phi0= TMath::ATan2(gy,gx);
3188 Int_t nsector=(gz>=0) ? 0:18;
3189 Float_t distPoint[3]={gx,gy,gz};
3190 corr->DistortPoint(distPoint, nsector);
3191 Double_t r0=TMath::Sqrt(gx*gx+gy*gy);
3192 Double_t r1=TMath::Sqrt(distPoint[0]*distPoint[0]+distPoint[1]*distPoint[1]);
3193 Double_t phi1=TMath::ATan2(distPoint[1],distPoint[0]);
3194 if (axisType==0) return r1-r0;
3195 if (axisType==1) return (phi1-phi0)*r0;
3196 if (axisType==2) return distPoint[2]-gz;
3200 Double_t AliTPCCorrection::GetDistXYZDz(Double_t gx, Double_t gy, Double_t gz, Int_t axisType, Int_t corrType,Double_t delta){
3202 // return correction at given x,y,z
3204 if (!fgVisualCorrection) return 0;
3205 AliTPCCorrection *corr = (AliTPCCorrection*)fgVisualCorrection->At(corrType);
3206 if (!corr) return 0;
3207 Double_t phi0= TMath::ATan2(gy,gx);
3208 Int_t nsector=(gz>=0) ? 0:18;
3209 Float_t distPoint[3]={gx,gy,gz};
3210 Float_t dxyz[3]={gx,gy,gz};
3212 corr->GetDistortionDz(distPoint, nsector,dxyz,delta);
3213 distPoint[0]+=dxyz[0];
3214 distPoint[1]+=dxyz[1];
3215 distPoint[2]+=dxyz[2];
3216 Double_t r0=TMath::Sqrt(gx*gx+gy*gy);
3217 Double_t r1=TMath::Sqrt(distPoint[0]*distPoint[0]+distPoint[1]*distPoint[1]);
3218 Double_t phi1=TMath::ATan2(distPoint[1],distPoint[0]);
3219 if (axisType==0) return r1-r0;
3220 if (axisType==1) return (phi1-phi0)*r0;
3221 if (axisType==2) return distPoint[2]-gz;
3225 Double_t AliTPCCorrection::GetDistXYZIntegrateZ(Double_t gx, Double_t gy, Double_t gz, Int_t axisType, Int_t corrType,Double_t delta){
3227 // return correction at given x,y,z
3229 if (!fgVisualCorrection) return 0;
3230 AliTPCCorrection *corr = (AliTPCCorrection*)fgVisualCorrection->At(corrType);
3231 if (!corr) return 0;
3232 Double_t phi0= TMath::ATan2(gy,gx);
3233 Int_t nsector=(gz>=0) ? 0:18;
3234 Float_t distPoint[3]={gx,gy,gz};
3235 Float_t dxyz[3]={gx,gy,gz};
3237 corr->GetDistortionIntegralDz(distPoint, nsector,dxyz,delta);
3238 distPoint[0]+=dxyz[0];
3239 distPoint[1]+=dxyz[1];
3240 distPoint[2]+=dxyz[2];
3241 Double_t r0=TMath::Sqrt(gx*gx+gy*gy);
3242 Double_t r1=TMath::Sqrt(distPoint[0]*distPoint[0]+distPoint[1]*distPoint[1]);
3243 Double_t phi1=TMath::ATan2(distPoint[1],distPoint[0]);
3244 if (axisType==0) return r1-r0;
3245 if (axisType==1) return (phi1-phi0)*r0;
3246 if (axisType==2) return distPoint[2]-gz;
3252 void AliTPCCorrection::MakeLaserDistortionTree(TTree* tree, TObjArray */*corrArray*/, Int_t /*itype*/){
3254 // Make a laser fit tree for global minimization
3256 AliTPCcalibDB* calib=AliTPCcalibDB::Instance();
3257 AliTPCCorrection * correction = calib->GetTPCComposedCorrection();
3258 if (!correction) correction = calib->GetTPCComposedCorrection(AliTrackerBase::GetBz());
3259 correction->AddVisualCorrection(correction,0); //register correction
3261 // AliTPCTransform *transform = AliTPCcalibDB::Instance()->GetTransform() ;
3262 //AliTPCParam *param = AliTPCcalibDB::Instance()->GetParameters();
3264 const Double_t cutErrY=0.05;
3265 const Double_t kSigmaCut=4;
3266 // const Double_t cutErrZ=0.03;
3267 const Double_t kEpsilon=0.00000001;
3268 // const Double_t kMaxDist=1.; // max distance - space correction
3273 AliTPCLaserTrack *ltr=0;
3274 AliTPCLaserTrack::LoadTracks();
3275 tree->SetBranchAddress("dY.",&vecdY);
3276 tree->SetBranchAddress("dZ.",&vecdZ);
3277 tree->SetBranchAddress("eY.",&veceY);
3278 tree->SetBranchAddress("eZ.",&veceZ);
3279 tree->SetBranchAddress("LTr.",<r);
3280 Int_t entries= tree->GetEntries();
3281 TTreeSRedirector *pcstream= new TTreeSRedirector("distortionLaser_0.root");
3282 Double_t bz=AliTrackerBase::GetBz();
3284 // Double_t globalXYZ[3];
3285 //Double_t globalXYZCorr[3];
3286 for (Int_t ientry=0; ientry<entries; ientry++){
3287 tree->GetEntry(ientry);
3288 if (!ltr->GetVecGX()){
3289 ltr->UpdatePoints();
3294 printf("Entry\t%d\n",ientry);
3295 for (Int_t irow0=0; irow0<158; irow0+=1){
3297 TLinearFitter fitter10(4,"hyp3");
3298 TLinearFitter fitter5(2,"hyp1");
3299 Int_t sector= (Int_t)(*ltr->GetVecSec())[irow0];
3300 if (sector<0) continue;
3301 //if (TMath::Abs(vecdY->GetMatrixArray()[irow0])<kEpsilon) continue;
3303 Double_t refX= (*ltr->GetVecLX())[irow0];
3304 Int_t firstRow1 = TMath::Max(irow0-10,0);
3305 Int_t lastRow1 = TMath::Min(irow0+10,158);
3306 Double_t padWidth=(irow0<64)?0.4:0.6;
3307 // make long range fit
3308 for (Int_t irow1=firstRow1; irow1<=lastRow1; irow1++){
3309 if (TMath::Abs((*ltr->GetVecSec())[irow1]-sector)>kEpsilon) continue;
3310 if (veceY->GetMatrixArray()[irow1]>cutErrY) continue;
3311 if (TMath::Abs(vecdY->GetMatrixArray()[irow1])<kEpsilon) continue;
3312 Double_t idealX= (*ltr->GetVecLX())[irow1];
3313 Double_t idealY= (*ltr->GetVecLY())[irow1];
3314 // Double_t idealZ= (*ltr->GetVecLZ())[irow1];
3315 Double_t gx= (*ltr->GetVecGX())[irow1];
3316 Double_t gy= (*ltr->GetVecGY())[irow1];
3317 Double_t gz= (*ltr->GetVecGZ())[irow1];
3318 Double_t measY=(*vecdY)[irow1]+idealY;
3319 Double_t deltaR = GetCorrXYZ(gx, gy, gz, 0,0);
3320 // deltaR = R distorted -R ideal
3321 Double_t xxx[4]={idealX+deltaR-refX,TMath::Cos(idealY/padWidth), TMath::Sin(idealY/padWidth)};
3322 fitter10.AddPoint(xxx,measY,1);
3325 Double_t rms10=0;//TMath::Sqrt(fitter10.GetChisquare()/(fitter10.GetNpoints()-4));
3326 Double_t mean10 =0;// fitter10.GetParameter(0);
3327 Double_t slope10 =0;// fitter10.GetParameter(0);
3328 Double_t cosPart10 = 0;// fitter10.GetParameter(2);
3329 Double_t sinPart10 =0;// fitter10.GetParameter(3);
3331 if (fitter10.GetNpoints()>10){
3333 rms10=TMath::Sqrt(fitter10.GetChisquare()/(fitter10.GetNpoints()-4));
3334 mean10 = fitter10.GetParameter(0);
3335 slope10 = fitter10.GetParameter(1);
3336 cosPart10 = fitter10.GetParameter(2);
3337 sinPart10 = fitter10.GetParameter(3);
3339 // make short range fit
3341 for (Int_t irow1=firstRow1+5; irow1<=lastRow1-5; irow1++){
3342 if (TMath::Abs((*ltr->GetVecSec())[irow1]-sector)>kEpsilon) continue;
3343 if (veceY->GetMatrixArray()[irow1]>cutErrY) continue;
3344 if (TMath::Abs(vecdY->GetMatrixArray()[irow1])<kEpsilon) continue;
3345 Double_t idealX= (*ltr->GetVecLX())[irow1];
3346 Double_t idealY= (*ltr->GetVecLY())[irow1];
3347 // Double_t idealZ= (*ltr->GetVecLZ())[irow1];
3348 Double_t gx= (*ltr->GetVecGX())[irow1];
3349 Double_t gy= (*ltr->GetVecGY())[irow1];
3350 Double_t gz= (*ltr->GetVecGZ())[irow1];
3351 Double_t measY=(*vecdY)[irow1]+idealY;
3352 Double_t deltaR = GetCorrXYZ(gx, gy, gz, 0,0);
3353 // deltaR = R distorted -R ideal
3354 Double_t expY= mean10+slope10*(idealX+deltaR-refX);
3355 if (TMath::Abs(measY-expY)>kSigmaCut*rms10) continue;
3357 Double_t corr=cosPart10*TMath::Cos(idealY/padWidth)+sinPart10*TMath::Sin(idealY/padWidth);
3358 Double_t xxx[4]={idealX+deltaR-refX,TMath::Cos(idealY/padWidth), TMath::Sin(idealY/padWidth)};
3359 fitter5.AddPoint(xxx,measY-corr,1);
3364 if (fitter5.GetNpoints()<8) isOK=kFALSE;
3366 Double_t rms5=0;//TMath::Sqrt(fitter5.GetChisquare()/(fitter5.GetNpoints()-4));
3367 Double_t offset5 =0;// fitter5.GetParameter(0);
3368 Double_t slope5 =0;// fitter5.GetParameter(0);
3371 rms5=TMath::Sqrt(fitter5.GetChisquare()/(fitter5.GetNpoints()-4));
3372 offset5 = fitter5.GetParameter(0);
3373 slope5 = fitter5.GetParameter(0);
3378 Double_t phi =(*ltr->GetVecPhi())[irow0];
3379 Double_t theta =ltr->GetTgl();
3380 Double_t mean=(vecdY)->GetMatrixArray()[irow0];
3381 Double_t gx=0,gy=0,gz=0;
3382 Double_t snp = (*ltr->GetVecP2())[irow0];
3383 Int_t bundle= ltr->GetBundle();
3384 Int_t id= ltr->GetId();
3385 // Double_t rms = err->GetMatrixArray()[irow];
3387 gx = (*ltr->GetVecGX())[irow0];
3388 gy = (*ltr->GetVecGY())[irow0];
3389 gz = (*ltr->GetVecGZ())[irow0];
3390 Double_t dRrec = GetCorrXYZ(gx, gy, gz, 0,0);
3391 fitter10.GetParameters(fit10);
3392 fitter5.GetParameters(fit5);
3393 Double_t idealY= (*ltr->GetVecLY())[irow0];
3394 Double_t measY=(*vecdY)[irow0]+idealY;
3395 Double_t corr=cosPart10*TMath::Cos(idealY/padWidth)+sinPart10*TMath::Sin(idealY/padWidth);
3396 if (TMath::Max(rms5,rms10)>0.06) isOK=kFALSE;
3398 (*pcstream)<<"fitFull"<< // dumpe also intermediate results
3399 "bz="<<bz<< // magnetic filed used
3400 "dtype="<<dtype<< // detector match type
3401 "ptype="<<ptype<< // parameter type
3402 "theta="<<theta<< // theta
3403 "phi="<<phi<< // phi
3404 "snp="<<snp<< // snp
3407 // // "dsec="<<dsec<<
3408 "refX="<<refX<< // reference radius
3409 "gx="<<gx<< // global position
3410 "gy="<<gy<< // global position
3411 "gz="<<gz<< // global position
3412 "dRrec="<<dRrec<< // delta Radius in reconstruction
3413 "id="<<id<< //bundle