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