]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCExBBShape.cxx
coveritiy
[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 // AliTPCExBBShape class                                                  //
18 ////////////////////////////////////////////////////////////////////////////
19
20 #include <AliMagF.h>
21 #include "TGeoGlobalMagField.h"
22 #include "AliTPCcalibDB.h"
23 #include "AliTPCParam.h"
24 #include "AliLog.h"
25
26 #include "AliTPCExBBShape.h"
27
28 AliTPCExBBShape::AliTPCExBBShape()
29   : AliTPCCorrection("exb_bshape","ExB B-shape"),
30     fC1(0.),fC2(0.),
31     fBField(0)
32 {
33   //
34   // default constructor
35   //
36 }
37
38 AliTPCExBBShape::~AliTPCExBBShape() {
39   //
40   // virtual destructor
41   //
42 }
43
44 void AliTPCExBBShape::Init() {
45   //
46   // Initialization funtion
47   //
48   
49   AliMagF* magF= (AliMagF*)TGeoGlobalMagField::Instance()->GetField();
50   if (!magF) AliError("Magneticd field - not initialized");
51   Double_t bzField = magF->SolenoidField()/10.; //field in T
52   SetBField(magF);
53   AliTPCParam *param= AliTPCcalibDB::Instance()->GetParameters();
54   if (!param) AliError("Parameters - not initialized");
55   Double_t vdrift = param->GetDriftV()/1000000.; // [cm/us]   // From dataBase: to be updated: per second (ideally)
56   Double_t ezField = 400; // [V/cm]   // to be updated: never (hopefully)
57   Double_t wt = -10.0 * (bzField*10) * vdrift / ezField ; 
58   // Correction Terms for effective omegaTau; obtained by a laser calibration run
59   SetOmegaTauT1T2(wt,fT1,fT2);
60
61
62 }
63
64 void AliTPCExBBShape::Update(const TTimeStamp &/*timeStamp*/) {
65   //
66   // Update function 
67   //
68   AliMagF* magF= (AliMagF*)TGeoGlobalMagField::Instance()->GetField();
69   if (!magF) AliError("Magneticd field - not initialized");
70   Double_t bzField = magF->SolenoidField()/10.; //field in T
71   SetBField(magF);
72   AliTPCParam *param= AliTPCcalibDB::Instance()->GetParameters();
73   if (!param) AliError("Parameters - not initialized");
74   Double_t vdrift = param->GetDriftV()/1000000.; // [cm/us]   // From dataBase: to be updated: per second (ideally)
75   Double_t ezField = 400; // [V/cm]   // to be updated: never (hopefully)
76   Double_t wt = -10.0 * (bzField*10) * vdrift / ezField ; 
77   // Correction Terms for effective omegaTau; obtained by a laser calibration run
78   SetOmegaTauT1T2(wt,fT1,fT2);
79
80
81 }
82
83
84
85 void AliTPCExBBShape::GetCorrection(const Float_t x[],const Short_t roc,Float_t dx[]) {
86   //
87   // Calculates the space point corrections of the B field inperfections (B field shape) 
88   //
89
90   if (!fBField) {
91     for (Int_t j=0;j<3;++j) dx[j]=0.;
92     return;
93   }
94
95   const Double_t xStart[3]={ x[0], x[1], x[2] };
96   const Double_t xEnd[3]={ x[0],  x[1],  roc%36<18?fgkTPCZ0:-fgkTPCZ0 };
97
98   Double_t intBStart[3];
99   Double_t intBEnd[3];
100
101   fBField->GetTPCRatInt(xStart,intBStart);
102   fBField->GetTPCRatInt(xEnd,  intBEnd  );
103
104   const Float_t intBxOverBz=(intBEnd[0]-intBStart[0]);
105   const Float_t intByOverBz=(intBEnd[1]-intBStart[1]);
106   
107   dx[0]=fC2*intBxOverBz-fC1*intByOverBz;
108   dx[1]=fC1*intBxOverBz+fC2*intByOverBz;
109   dx[2]=0.;
110
111
112 }
113
114 void AliTPCExBBShape::GetBxAndByOverBz(const Float_t x[],const Short_t roc,Float_t BxByOverBz[]) {
115   //
116   // This function is purely for calibration purposes
117   // Returns the via AliMagF obtaind B field integrals  
118   // 
119
120   if (!fBField) {
121     for (Int_t j=0;j<3;++j) BxByOverBz[j]=0.;
122     return;
123   }
124
125   const Double_t xStart[3]={ x[0], x[1], x[2] };
126   const Double_t xEnd[3]={ x[0],  x[1],  roc%36<18?fgkTPCZ0:-fgkTPCZ0 };
127
128   Double_t intBStart[3];
129   Double_t intBEnd[3];
130
131   fBField->GetTPCRatInt(xStart,intBStart);
132   fBField->GetTPCRatInt(xEnd,  intBEnd  );
133
134   const Float_t intBxOverBz=(intBEnd[0]-intBStart[0]);
135   const Float_t intByOverBz=(intBEnd[1]-intBStart[1]);
136   
137   BxByOverBz[0]=intBxOverBz;
138   BxByOverBz[1]=intByOverBz;
139
140 }
141
142 void AliTPCExBBShape::Print(Option_t* option) const {
143   //
144   // Print function to check the settings (e.g. voltage offsets)
145   // option=="a" prints details of the B field settings and the 
146   // C0 and C1 coefficents (for calibration purposes)
147   //
148   TString opt = option; opt.ToLower();
149   printf("%s\t%s\n - B field settings:\n",GetTitle(),GetName());
150   fBField->Print(option);
151   //  printf(" - B field: X-Twist: %1.5lf rad, Y-Twist: %1.5lf rad \n",fBField->Print(option));
152   if (opt.Contains("a")) { // Print all details  
153     printf(" - T1: %1.4f, T2: %1.4f \n",fT1,fT2);
154     printf(" - C1: %1.4f, C2: %1.4f \n",fC1,fC2);
155   }    
156 }
157
158 Double_t AliTPCExBBShape::GetBFieldXYZ(Double_t gx, Double_t gy, Double_t gz, Int_t axisType){
159   //
160   // return B field at given x,y,z
161   // 
162   AliMagF* field = (AliMagF*)TGeoGlobalMagField::Instance()->GetField();
163   if (!field) return 0;
164   Double_t xyz[3]={gx,gy,gz};
165   Double_t bxyz[3]={0};
166   field->Field(xyz,bxyz);
167   //
168   Double_t r=TMath::Sqrt(gx*gx+gy*gy);
169   //  Double_t b=TMath::Sqrt(bxyz[0]*bxyz[0]+bxyz[1]*bxyz[1]);
170   if (axisType==0) {
171     return (xyz[0]*bxyz[1]-xyz[1]*bxyz[0])/(bxyz[2]*r);
172   }
173   if (axisType==1){
174     return (xyz[0]*bxyz[0]+xyz[1]*bxyz[1])/(bxyz[2]*r);
175   }
176   if (axisType==2) return bxyz[2];
177   if (axisType==3) return bxyz[0];
178   if (axisType==4) return bxyz[1];
179   if (axisType==5) return bxyz[2];
180   return bxyz[2];
181 }