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