StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
StarMagField.h
1 /***********************************************************************
2  *
3  * $Id: StarMagField.h,v 1.17 2017/04/26 21:11:25 perev Exp $
4  *
5  * Author: Jim Thomas 11/1/2000
6  *
7  ***********************************************************************
8  *
9  * Description: the STAR Magnetic Field
10  *
11  ***********************************************************************
12  *
13  * $Log: StarMagField.h,v $
14  * Revision 1.17 2017/04/26 21:11:25 perev
15  * Add setConstBz
16  *
17  * Revision 1.16 2014/06/26 21:50:17 fisyak
18  * New Tpc Alignment, v632
19  *
20  * Revision 1.15 2013/03/26 13:38:18 fisyak
21  * restore back modififcations as not related to drop in no. of reconstructed tracks
22  *
23  * Revision 1.13 2013/01/17 15:11:33 fisyak
24  * More clear handling ROOT and non ROOT versions
25  *
26  * Revision 1.11 2013/01/16 00:05:27 fisyak
27  * Add TObject
28  *
29  * Revision 1.10 2013/01/15 23:45:02 fisyak
30  * Account ROOT version with TVirtualMagField
31  *
32  * Revision 1.9 2013/01/15 17:35:23 fisyak
33  * Create clean versions of ROOT and non ROOT StarMagField
34  *
35  * Revision 1.8 2009/12/11 14:19:07 fisyak
36  * switch from define to enum
37  *
38  * Revision 1.7 2009/12/08 15:33:57 fisyak
39  * Hold replacement defines via enum till StMagUtilities will be updated
40  *
41  * Revision 1.6 2009/12/07 23:38:15 fisyak
42  * Move size definition from #define to enumerations
43  *
44  * Revision 1.5 2009/01/13 03:19:44 perev
45  * Mag field nou controlled from starsim. BugFix
46  *
47  * Revision 1.4 2007/09/21 21:07:08 fisyak
48  * Remove ClassDef and ClassImp
49  *
50  * Revision 1.3 2007/09/13 00:00:27 fisyak
51  * add mag.field in steel, from Lijuan Ruan
52  *
53  * Revision 1.2 2005/07/28 19:46:01 fisyak
54  * Add:
55  * - frindge magnetic field from P.Nevski extrapolation,
56  * - lock mechanism for Magnetic Field parameters
57  * Remove:
58  * - Electric Field (Still part of StDbUtilitie/StMagUtilities)
59  *
60  * Comparision between mfldgeo and StarMagField calculation for Z and R components of mag.field
61  * is at http://www.star.bnl.gov/~fisyak/star/MagField/
62  *
63  * Revision 1.1.1.1 2005/07/07 14:13:55 fisyak
64  * The version of STAR mag. field extracted from StDbUtilities/StMagUtilities to be used in Simulation and Reconstruction instead of agufld
65  *
66  * Revision 1.2 2005/07/07 14:07:55 fisyak
67  * Final version before moving to official repository
68  *
69  * Revision 1.1 2004/03/12 13:26:24 fisyak
70  * Singleton for STAR magnetic field
71  *
72  ***********************************************************************/
74 // //
75 // StarMagField Class //
76 // //
78 #ifndef StarMagField_H
79 #define StarMagField_H
80 
81 #include <math.h>
82 #include <stdio.h>
83 #include <stdlib.h>
84 #include <Stiostream.h>
85 #include <Rtypes.h>
86 #if defined (__ROOT__)
87 #include "TH2.h"
88 #include "TGeoMatrix.h"
89 #if ROOT_VERSION_CODE >= 335360 /* ROOT_VERSION(5,30,0) */
90 #include "TVirtualMagField.h"
91 class StarMagField : public TVirtualMagField
92 #else
93 #include "TObject.h"
94 class StarMagField : public TObject
95 #endif
96 #else
98 #endif
99 {
100  public:
101  enum EBField { kUndefined = 0, kConstant = 1, kMapped = 2, kChain = 3 } ;
102  enum ESmFSizes {nZ = 57, nR = 28, nPhi = 37, nZSteel = 16, nRSteel = 115, nPhiSteel = 25};
103  static void Search ( Int_t N, const Float_t Xarray[], Float_t x, Int_t &low ) ;
104  virtual Float_t Interpolate ( const Float_t Xarray[], const Float_t Yarray[],
105  const Int_t ORDER, const Float_t x ) ;
106  virtual void Interpolate2DBfield ( const Float_t r, const Float_t z,
107  Float_t &Br_value, Float_t &Bz_value ) ;
108  virtual void Interpolate2ExtDBfield ( const Float_t r, const Float_t z,
109  Float_t &Br_value, Float_t &Bz_value ) ;
110  virtual void Interpolate3DBfield ( const Float_t r, const Float_t z, const Float_t phi,
111  Float_t &Br_value, Float_t &Bz_value, Float_t &Bphi_value ) ;
112 
113  private:
114 
115  virtual void ReadField ( ) ;
116  static StarMagField *fgInstance;
117 #if defined (__ROOT__)
118  TGeoRotation fStarMagFieldRotation;
119  TH2F *fBzdZCorrection; // correction due to endcap calomiter
120  TH2F *fBrdZCorrection; // correction due to endcap calomiter
121 #endif
122  public:
123  //added by Lijuan
124 
125  virtual void Interpolate3DBSteelfield ( const Float_t r, const Float_t z, const Float_t phi,
126  Float_t &Br_value, Float_t &Bz_value, Float_t &Bphi_value );
127  //end added by Lijuan
128 
129 
130 
131  EBField fMap; // (D) = kMapped; Global flag to indicate static arrays are full
132  Float_t fFactor; // (D) = 1.0 ; Multiplicative factor (allows scaling and sign reversal)
133  Float_t fRescale; // (D) = 1.0 ; Multiplicative factor (allows re-scaling wrt which map read)
134  Float_t fBDipole; // (D) = -42.67 ; field value (kG)
135  Float_t fRmaxDip; // (D) = 15.34 ; Inner field volume radius
136  Float_t fZminDip; // (D) = 980.0 ; StArt of the DX mAgnet in Z
137  Float_t fZmaxDip; // (D) = 1350.0 ; End of the DX mAgnet in Z
138  Bool_t fLock; // (D) = kFALSE ; Set kTRUE if lock above values
139 
140  Float_t Bz[nZ][nR], Br[nZ][nR] ;
141  Float_t Radius[nR], ZList[nZ] ;
142  Float_t Bz3D[nPhi][nZ][nR], Br3D[nPhi][nZ][nR], Bphi3D[nPhi][nZ][nR] ;
143  Float_t R3D[nR], Z3D[nZ], Phi3D[nPhi] ;
144  Float_t R3DSteel[nRSteel], Z3DSteel[nZSteel], Phi3DSteel[nPhiSteel] ;
145  Float_t Bz3DSteel[nPhiSteel][nZSteel][nRSteel];
146  Float_t Bx3DSteel[nPhiSteel][nZSteel][nRSteel], By3DSteel[nPhiSteel][nZSteel][nRSteel] ;
147  static bool mConstBz;
148 
149  //added by Lijuan
150  Float_t Br3DSteel[nPhiSteel][nZSteel][nRSteel], Bphi3DSteel[nPhiSteel][nZSteel][nRSteel] ;
151  //end added by Lijuan
152 
153  public:
154 
155  StarMagField ( EBField map = kMapped, Float_t Factor = 1,
156  Bool_t Lock = kFALSE, Float_t Rescale = 1,
157  Float_t Bdipole = -42.67, Float_t Rmaxdip = 15.34,
158  Float_t Zmindip = 980.0, Float_t Zmaxdip = 1350.0) ;
159  virtual ~StarMagField () {
160  fgInstance = 0;
161 #if defined (__ROOT__)
162  SafeDelete(fBzdZCorrection);
163  SafeDelete(fBrdZCorrection);
164 #endif
165  }
166  static StarMagField *Instance();
167 
168  static void setConstBz( bool state ){ mConstBz = state; }
169 
170  virtual void BField ( const Float_t x[], Float_t B[] ) ;
171  virtual void BField ( const Double_t x[], Double_t B[] ) ;
172  virtual void Field ( const Float_t x[], Float_t B[] ) {BField(x,B);}
173  virtual void Field ( const Double_t x[], Double_t B[] ) {BField(x,B);}
174  virtual void BrBzField( const Float_t r, const Float_t z, Float_t &Br_value, Float_t &Bz_value ) ;
175  virtual void B3DField ( const Float_t x[], Float_t B[] ) ;
176  virtual void B3DField ( const Double_t x[], Double_t B[] ) ;
177  virtual void BrBz3DField ( const Float_t r, const Float_t z, const Float_t phi,
178  Float_t &Br_value, Float_t &Bz_value, Float_t &Bphi_value ) ;
179  virtual void SetFactor (Float_t factor = 1);
180  virtual void SetRescale(Float_t factor = 1);
181  virtual void SetBDipole(Float_t m = -42.67);
182  virtual void SetRmaxDip(Float_t m = 15.3);
183  virtual void SetZminDip(Float_t m = 980.0);
184  virtual void SetZmaxDip(Float_t m = 1350.0);
185  virtual void SetLock();
186  virtual EBField GetMap() {return fMap;}
187  virtual Float_t GetFactor() {return fFactor;}
188  virtual Float_t GetRescale() {return fRescale;}
189  virtual Bool_t IsLocked() {return fLock;}
190  virtual void Print(Option_t* opt="") const;
191 #ifdef __ROOT__
192  void SetStarMagFieldRotation(TGeoRotation &rot);
193  void SetStarMagFieldRotation(Double_t *rot);
194  const TGeoRotation &StarMagFieldRotation() {return *&fStarMagFieldRotation;}
195  ClassDef(StarMagField,1) // Base class for all STAR MagField
196 #endif
197 };
198 
199 #endif
virtual void Interpolate3DBSteelfield(const Float_t r, const Float_t z, const Float_t phi, Float_t &Br_value, Float_t &Bz_value, Float_t &Bphi_value)
Interpolate the B field map - 3D interpolation.
virtual void Interpolate2DBfield(const Float_t r, const Float_t z, Float_t &Br_value, Float_t &Bz_value)
Interpolate the B field map - 2D interpolation.
virtual void B3DField(const Float_t x[], Float_t B[])
Bfield in Cartesian coordinates - 3D field.
static void Search(Int_t N, const Float_t Xarray[], Float_t x, Int_t &low)
Search an ordered table by starting at the most recently used point.
virtual Float_t Interpolate(const Float_t Xarray[], const Float_t Yarray[], const Int_t ORDER, const Float_t x)
Interpolate a 3x2 table (quadratic) or a 2x2 table (linear)
virtual void Interpolate3DBfield(const Float_t r, const Float_t z, const Float_t phi, Float_t &Br_value, Float_t &Bz_value, Float_t &Bphi_value)
Interpolate the B field map - 3D interpolation.
virtual void BrBzField(const Float_t r, const Float_t z, Float_t &Br_value, Float_t &Bz_value)
B field in Radial coordinates - 2D field (ie Phi symmetric)