![]() |
AliRoot Core
d69033e (d69033e)
|
Calculation of the electric field. More...
#include <AliTPCEfield.h>
Public Member Functions | |
AliTPCEfield () | |
AliTPCEfield (const char *name, Int_t maxFreq, Bool_t is2D, Bool_t useLinear=kTRUE) | |
virtual | ~AliTPCEfield () |
void | SetRange (Double_t x0, Double_t x1, Double_t y0, Double_t y1, Double_t z00, Double_t z1=0) |
void | AddBoundaryLine (Double_t x0, Double_t y0, Double_t z0, Double_t v0, Double_t x1, Double_t y1, Double_t z1, Double_t v1, Int_t id=0, Int_t npoints=100) |
TTree * | GetTree (const char *tname="Boundary") |
void | MakeFitFunctions (Int_t maxFreq) |
void | FitField () |
void | DumpField (Double_t gridSize=5, Double_t step=0.5) |
Double_t | EvalField (Int_t ifun, Double_t x, Double_t y, Double_t z, Int_t type=0) |
Double_t | Eval (Double_t x, Double_t y, Double_t z, Int_t type=0) |
Double_t | Field (Int_t ftype, Double_t ifx, Double_t ify, Double_t ifz, Double_t x, Double_t y, Double_t z) |
Double_t | FieldDn (Int_t ftype, Double_t ifx, Double_t ify, Double_t ifz, Int_t dn, Double_t x, Double_t y, Double_t z) |
TMatrixD * | MakeCorrelation (TMatrixD &matrix) |
Double_t | SinHNorm (Double_t x, Double_t norm) |
Double_t | CosHNorm (Double_t x, Double_t norm) |
Static Public Member Functions | |
static Double_t | EvalS (Double_t x, Double_t y, Double_t z, Int_t type=0) |
Public Attributes | |
Double_t | fMin [3] |
range of coordinates from Min to Max More... | |
Double_t | fMax [3] |
Double_t | fScale |
scaling factor More... | |
Int_t | fMaxFreq |
maximal frequency of expansion More... | |
Bool_t | fIs2D |
flag for 2D field More... | |
Bool_t | fUseLinear |
flag to use also linear term of the field More... | |
TTreeSRedirector * | fWorkspace |
! workspace More... | |
TMatrixD * | fFitFunctions |
fit function description More... | |
TVectorD * | fFitParam |
fit parameters - coeficients More... | |
TMatrixD * | fFitCovar |
fit covariance More... | |
TLinearFitter * | fFitter |
linear fitter - temporary solution - integrals to be calculated More... | |
Static Public Attributes | |
static AliTPCEfield * | fgInstance =0 |
instance of fied - for visualization More... | |
Calculation of the electric field.
Solution of Laplace equation in cartesian system, with boundary condition.
See details: http://web.mit.edu/6.013_book/www/chapter5/5.10.html
Definition at line 16 of file AliTPCEfield.h.
AliTPCEfield::AliTPCEfield | ( | ) |
Definition at line 36 of file AliTPCEfield.cxx.
AliTPCEfield::AliTPCEfield | ( | const char * | name, |
Int_t | maxFreq, | ||
Bool_t | is2D, | ||
Bool_t | useLinear = kTRUE |
||
) |
Definition at line 51 of file AliTPCEfield.cxx.
|
virtual |
Destructor
Definition at line 124 of file AliTPCEfield.cxx.
void AliTPCEfield::AddBoundaryLine | ( | Double_t | x0, |
Double_t | y0, | ||
Double_t | z0, | ||
Double_t | v0, | ||
Double_t | x1, | ||
Double_t | y1, | ||
Double_t | z1, | ||
Double_t | v1, | ||
Int_t | id = 0 , |
||
Int_t | npoints = 100 |
||
) |
Add a e field boundary line From point (x0,y0) to point (x1,y1) Linear decrease of potential is assumed Boundary can be identified using boundary ID The line is written into tree Boundary
Definition at line 142 of file AliTPCEfield.cxx.
Referenced by MakeTPC2DExample().
|
inline |
Definition at line 42 of file AliTPCEfield.h.
void AliTPCEfield::DumpField | ( | Double_t | gridSize = 5 , |
Double_t | step = 0.5 |
||
) |
Definition at line 380 of file AliTPCEfield.cxx.
Double_t AliTPCEfield::Eval | ( | Double_t | x, |
Double_t | y, | ||
Double_t | z, | ||
Int_t | type = 0 |
||
) |
Evaluate function ifun at position gx amd gy type == 0 - field == 1 - Ex == 2 - Ey == 3 - Ez
Definition at line 278 of file AliTPCEfield.cxx.
Referenced by DumpField(), and EvalS().
Double_t AliTPCEfield::EvalField | ( | Int_t | ifun, |
Double_t | x, | ||
Double_t | y, | ||
Double_t | z, | ||
Int_t | type = 0 |
||
) |
Evaluate function ifun at position gx amd gy type == 0 - field == 1 - Ex == 2 - Ey == 3 - Ez
Definition at line 260 of file AliTPCEfield.cxx.
Referenced by Eval(), and FitField().
|
static |
static evaluation - possible to use it in the TF1
Definition at line 298 of file AliTPCEfield.cxx.
Double_t AliTPCEfield::Field | ( | Int_t | ftype, |
Double_t | ifx, | ||
Double_t | ify, | ||
Double_t | ifz, | ||
Double_t | x, | ||
Double_t | y, | ||
Double_t | z | ||
) |
Field component in f frequency
Definition at line 175 of file AliTPCEfield.cxx.
Referenced by EvalField().
Double_t AliTPCEfield::FieldDn | ( | Int_t | ftype, |
Double_t | ifx, | ||
Double_t | ify, | ||
Double_t | ifz, | ||
Int_t | dn, | ||
Double_t | x, | ||
Double_t | y, | ||
Double_t | z | ||
) |
Field component in f frequency
Definition at line 213 of file AliTPCEfield.cxx.
Referenced by EvalField().
void AliTPCEfield::FitField | ( | ) |
Fit the e field Minimize chi2 residuals at the boundary points ?Tempoary sollution - integrals can be calculated analytically -
Definition at line 304 of file AliTPCEfield.cxx.
TTree * AliTPCEfield::GetTree | ( | const char * | tname = "Boundary" | ) |
Definition at line 169 of file AliTPCEfield.cxx.
Referenced by FitField().
TMatrixD * AliTPCEfield::MakeCorrelation | ( | TMatrixD & | matrix | ) |
Definition at line 365 of file AliTPCEfield.cxx.
void AliTPCEfield::MakeFitFunctions | ( | Int_t | maxFreq | ) |
fit functions = \( f(x,y,z) = fx(x)*fy(y)*fz(z) \) function can be of following types: 0 - constant 1 - linear 2 - hypx 3 - hypy 4 - hypz
Definition at line 69 of file AliTPCEfield.cxx.
Referenced by AliTPCEfield().
void AliTPCEfield::SetRange | ( | Double_t | x0, |
Double_t | x1, | ||
Double_t | y0, | ||
Double_t | y1, | ||
Double_t | z00, | ||
Double_t | z1 = 0 |
||
) |
Set the ranges - coordinates are rescaled in order to use proper cos,sin expansion in scaled space
Definition at line 130 of file AliTPCEfield.cxx.
Referenced by MakeTPC2DExample().
|
inline |
Definition at line 41 of file AliTPCEfield.h.
TMatrixD* AliTPCEfield::fFitCovar |
fit covariance
Definition at line 54 of file AliTPCEfield.h.
Referenced by FitField(), and MakeFitFunctions().
TMatrixD* AliTPCEfield::fFitFunctions |
fit function description
Definition at line 52 of file AliTPCEfield.h.
Referenced by Eval(), EvalField(), FitField(), and MakeFitFunctions().
TVectorD* AliTPCEfield::fFitParam |
fit parameters - coeficients
Definition at line 53 of file AliTPCEfield.h.
Referenced by FitField(), and MakeFitFunctions().
TLinearFitter* AliTPCEfield::fFitter |
linear fitter - temporary solution - integrals to be calculated
Definition at line 55 of file AliTPCEfield.h.
Referenced by FitField().
|
static |
instance of fied - for visualization
Definition at line 56 of file AliTPCEfield.h.
Referenced by AliTPCEfield(), and EvalS().
Bool_t AliTPCEfield::fIs2D |
flag for 2D field
Definition at line 48 of file AliTPCEfield.h.
Referenced by MakeFitFunctions(), and SetRange().
Double_t AliTPCEfield::fMax[3] |
Definition at line 45 of file AliTPCEfield.h.
Referenced by AliTPCEfield(), DumpField(), Eval(), FitField(), and SetRange().
Int_t AliTPCEfield::fMaxFreq |
maximal frequency of expansion
Definition at line 47 of file AliTPCEfield.h.
Referenced by DumpField().
Double_t AliTPCEfield::fMin[3] |
range of coordinates from Min to Max
Definition at line 44 of file AliTPCEfield.h.
Referenced by AliTPCEfield(), DumpField(), Eval(), FitField(), and SetRange().
Double_t AliTPCEfield::fScale |
scaling factor
Definition at line 46 of file AliTPCEfield.h.
Referenced by DumpField(), Eval(), FitField(), and SetRange().
Bool_t AliTPCEfield::fUseLinear |
flag to use also linear term of the field
Definition at line 49 of file AliTPCEfield.h.
Referenced by MakeFitFunctions().
TTreeSRedirector* AliTPCEfield::fWorkspace |
! workspace
Definition at line 51 of file AliTPCEfield.h.
Referenced by AliTPCEfield(), GetTree(), and ~AliTPCEfield().