]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TPC/TPCbase/AliTPCExBTwist.cxx
TPC module
[u/mrichter/AliRoot.git] / TPC / TPCbase / AliTPCExBTwist.cxx
CommitLineData
0116859c 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////////////////////////////////////////////////////////////////////////////
0116859c 17// AliTPCExBTwist class //
0116859c 18////////////////////////////////////////////////////////////////////////////
b4caed64 19
20
534fd34a 21#include "AliMagF.h"
22#include "TGeoGlobalMagField.h"
23#include "AliTPCcalibDB.h"
24#include "AliTPCParam.h"
25#include "AliLog.h"
0116859c 26
27#include "AliTPCExBTwist.h"
28
29AliTPCExBTwist::AliTPCExBTwist()
30 : AliTPCCorrection("exb_twist","ExB twist"),
31 fC1(0.),fC2(0.),
32 fXTwist(0.),fYTwist(0.)
33{
34 //
35 // default constructor
36 //
37}
38
39AliTPCExBTwist::~AliTPCExBTwist() {
40 //
41 // default destructor
42 //
43}
44
69d03c4d 45Bool_t AliTPCExBTwist::AddCorrectionCompact(AliTPCCorrection* corr, Double_t weight){
46 //
47 // Add correction and make them compact
48 // Assumptions:
49 // - origin of distortion/correction are additive
50 // - only correction ot the same type supported ()
51 if (corr==NULL) {
52 AliError("Zerro pointer - correction");
53 return kFALSE;
54 }
55 AliTPCExBTwist * corrC = dynamic_cast< AliTPCExBTwist*>(corr);
56 if (corrC == NULL) return kFALSE;
57 fXTwist+=weight*corrC->fXTwist; // Twist of E to B field in X-Z [rad]
58 fYTwist+=weight*corrC->fYTwist; // Twist of E to B field in Y-Z [rad]
59 return kTRUE;
60}
61
534fd34a 62
63
e527a1b9 64void AliTPCExBTwist::Init() {
65 //
534fd34a 66 // Initialization funtion
e527a1b9 67 //
68
534fd34a 69 AliMagF* magF= (AliMagF*)TGeoGlobalMagField::Instance()->GetField();
70 if (!magF) AliError("Magneticd field - not initialized");
71 Double_t bzField = magF->SolenoidField()/10.; //field in T
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)
e527a1b9 75 Double_t ezField = 400; // [V/cm] // to be updated: never (hopefully)
76 Double_t wt = -10.0 * (bzField*10) * vdrift / ezField ;
e527a1b9 77 // Correction Terms for effective omegaTau; obtained by a laser calibration run
534fd34a 78 SetOmegaTauT1T2(wt,fT1,fT2);
e527a1b9 79
e527a1b9 80
e527a1b9 81}
82
83void AliTPCExBTwist::Update(const TTimeStamp &/*timeStamp*/) {
84 //
85 // Update function
86 //
534fd34a 87 AliMagF* magF= (AliMagF*)TGeoGlobalMagField::Instance()->GetField();
88 if (!magF) AliError("Magneticd field - not initialized");
89 Double_t bzField = magF->SolenoidField()/10.; //field in T
90 AliTPCParam *param= AliTPCcalibDB::Instance()->GetParameters();
91 if (!param) AliError("Parameters - not initialized");
92 Double_t vdrift = param->GetDriftV()/1000000.; // [cm/us] // From dataBase: to be updated: per second (ideally)
e527a1b9 93 Double_t ezField = 400; // [V/cm] // to be updated: never (hopefully)
94 Double_t wt = -10.0 * (bzField*10) * vdrift / ezField ;
e527a1b9 95 // Correction Terms for effective omegaTau; obtained by a laser calibration run
534fd34a 96 SetOmegaTauT1T2(wt,fT1,fT2);
97
e527a1b9 98
e527a1b9 99}
100
101
102
0116859c 103void AliTPCExBTwist::GetCorrection(const Float_t x[],const Short_t roc,Float_t dx[]) {
104 //
105 // Calculates the correction of a mismatch between the E and B field axis
106 //
107
108 const Float_t zstart=x[2];
b1f0a2a5 109 const Float_t zend =(roc%36<18?fgkTPCZ0:-fgkTPCZ0);
0116859c 110 const Float_t zdrift=zstart-zend;
111
112 dx[0]=(fC2*fXTwist-fC1*fYTwist)*zdrift;
113 dx[1]=(fC1*fXTwist+fC2*fYTwist)*zdrift;
114 dx[2]=0.;
115}
116
b1f0a2a5 117void AliTPCExBTwist::Print(const Option_t* option) const {
0116859c 118 //
119 // Print function to check the settings (e.g. the twist in the X direction)
120 // option=="a" prints the C0 and C1 coefficents for calibration purposes
121 //
122
123 TString opt = option; opt.ToLower();
124 printf("%s\n",GetTitle());
125
126 printf(" - Twist settings: X-Twist: %1.5f rad, Y-Twist: %1.5f rad \n",fXTwist,fYTwist);
127 if (opt.Contains("a")) { // Print all details
534fd34a 128 printf(" - T1: %1.4f, T2: %1.4f \n",fT1,fT2);
0116859c 129 printf(" - C1: %1.4f, C2: %1.4f \n",fC1,fC2);
130 }
131
132
133}