]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCExBBShape.cxx
M AliTPCcalibCalib.cxx - start refit form TPC out instead of track...
[u/mrichter/AliRoot.git] / TPC / AliTPCExBBShape.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 ////////////////////////////////////////////////////////////////////////////
17 //                                                                        //
18 // AliExBBShape class                                                     //
19 // The class calculates the space point distortions due to the B field    //
20 // shape imperfections using a second order technique based on integrals  //
21 // over Bz (e.g. int By/Bz) obtained via the AliMagF class                //
22 // The class allows "effective Omega Tau" corrections.                    //
23 //                                                                        //
24 // date: 27/04/2010                                                       //
25 // Authors: Magnus Mager, Jim Thomas, Stefan Rossegger                    //
26 //                                                                        //
27 // Example usage:                                                         //
28 //  AliMagF mag("mag","mag");                                             //
29 //  AliTPCExBBShape exb;                                                  //
30 //  exb.SetBField(&mag);             // use Bfield from AliMagF           //
31 //  exb.SetOmegaTauT1T2(0.32,1.,1.); // values ideally from OCDB          //
32 //  // plot dRPhi distortions ...                                         //
33 //  exb.CreateHistoDRPhiinZR(0,100,100)->Draw("surf2");                   //
34 ////////////////////////////////////////////////////////////////////////////
35
36 #include <AliMagF.h>
37
38 #include "AliTPCExBBShape.h"
39
40 AliTPCExBBShape::AliTPCExBBShape()
41   : AliTPCCorrection("exb_bshape","ExB B-shape"),
42     fC1(0.),fC2(0.),
43     fBField(0)
44 {
45   //
46   // default constructor
47   //
48 }
49
50 AliTPCExBBShape::~AliTPCExBBShape() {
51   //
52   // virtual destructor
53   //
54 }
55
56 void AliTPCExBBShape::Init() {
57   //
58   // Initialization funtion (not used at the moment)
59   //
60   
61   // Set default parameters
62   // FIXME: Ask the database for these entries
63   
64
65   AliMagF * mag = new AliMagF("mag","mag"); // from database (GRP?)
66   SetBField(mag);
67
68   Double_t vdrift = 2.6; // [cm/us]   // From dataBase: to be updated: per second (ideally)
69   Double_t bzField = -0.5; // [Tesla] // From dataBase: to be updated: per run
70
71   Double_t ezField = 400; // [V/cm]   // to be updated: never (hopefully)
72   Double_t wt = -10.0 * (bzField*10) * vdrift / ezField ; 
73
74   // Correction Terms for effective omegaTau; obtained by a laser calibration run
75   Double_t t1 = 0.9;   // ideally from database
76   Double_t t2 = 1.5;   // ideally from database
77
78   SetOmegaTauT1T2(wt,t1,t2);
79
80
81 }
82
83 void AliTPCExBBShape::Update(const TTimeStamp &/*timeStamp*/) {
84   //
85   // Update function 
86   //
87
88   Double_t vdrift = 2.6; // [cm/us]   // From dataBase: to be updated: per second (ideally)
89   Double_t bzField = -0.5; // [Tesla] // From dataBase: to be updated: per run
90
91   Double_t ezField = 400; // [V/cm]   // to be updated: never (hopefully)
92   Double_t wt = -10.0 * (bzField*10) * vdrift / ezField ; 
93
94   // Correction Terms for effective omegaTau; obtained by a laser calibration run
95   Double_t t1 = 0.9;   // ideally from database
96   Double_t t2 = 1.5;   // ideally from database
97
98   SetOmegaTauT1T2(wt,t1,t2);
99
100
101 }
102
103
104
105 void AliTPCExBBShape::GetCorrection(const Float_t x[],const Short_t roc,Float_t dx[]) {
106   //
107   // Calculates the space point corrections of the B field inperfections (B field shape) 
108   //
109
110   if (!fBField) {
111     for (Int_t j=0;j<3;++j) dx[j]=0.;
112     return;
113   }
114
115   const Double_t xStart[3]={ x[0], x[1], x[2] };
116   const Double_t xEnd[3]={ x[0],  x[1],  roc%36<18?fgkTPC_Z0:-fgkTPC_Z0 };
117
118   Double_t intBStart[3];
119   Double_t intBEnd[3];
120
121   fBField->GetTPCRatInt(xStart,intBStart);
122   fBField->GetTPCRatInt(xEnd,  intBEnd  );
123
124   const Float_t intBxOverBz=(intBEnd[0]-intBStart[0]);
125   const Float_t intByOverBz=(intBEnd[1]-intBStart[1]);
126   
127   dx[0]=fC2*intBxOverBz-fC1*intByOverBz;
128   dx[1]=fC1*intBxOverBz+fC2*intByOverBz;
129   dx[2]=0.;
130
131
132 }
133
134 void AliTPCExBBShape::GetBxAndByOverBz(const Float_t x[],const Short_t roc,Float_t BxByOverBz[]) {
135   //
136   // This function is purely for calibration purposes
137   // Returns the via AliMagF obtaind B field integrals  
138   // 
139
140   if (!fBField) {
141     for (Int_t j=0;j<3;++j) BxByOverBz[j]=0.;
142     return;
143   }
144
145   const Double_t xStart[3]={ x[0], x[1], x[2] };
146   const Double_t xEnd[3]={ x[0],  x[1],  roc%36<18?fgkTPC_Z0:-fgkTPC_Z0 };
147
148   Double_t intBStart[3];
149   Double_t intBEnd[3];
150
151   fBField->GetTPCRatInt(xStart,intBStart);
152   fBField->GetTPCRatInt(xEnd,  intBEnd  );
153
154   const Float_t intBxOverBz=(intBEnd[0]-intBStart[0]);
155   const Float_t intByOverBz=(intBEnd[1]-intBStart[1]);
156   
157   BxByOverBz[0]=intBxOverBz;
158   BxByOverBz[1]=intByOverBz;
159
160 }
161
162 void AliTPCExBBShape::Print(Option_t* option) const {
163   //
164   // Print function to check the settings (e.g. voltage offsets)
165   // option=="a" prints details of the B field settings and the 
166   // C0 and C1 coefficents (for calibration purposes)
167   //
168   TString opt = option; opt.ToLower();
169   printf("%s\n - B field settings:\n",GetTitle());
170   fBField->Print(option);
171   //  printf(" - B field: X-Twist: %1.5lf rad, Y-Twist: %1.5lf rad \n",fBField->Print(option));
172   if (opt.Contains("a")) { // Print all details
173     printf(" - C1: %1.4f, C2: %1.4f \n",fC1,fC2);
174   }    
175  
176 }