]> git.uio.no Git - u/mrichter/AliRoot.git/blame - STEER/AliMagF.cxx
Introduce factor for magnetic field
[u/mrichter/AliRoot.git] / STEER / AliMagF.cxx
CommitLineData
4c039060 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$Log$
dd045843 18Revision 1.3 1999/09/29 09:24:29 fca
19Introduction of the Copyright and cvs Log
20
4c039060 21*/
22
fe4da5cc 23
24#include "AliMagF.h"
25#include "TSystem.h"
26#include <stdlib.h>
27#include <stdio.h>
28
29//ZDC part -------------------------------------------------------------------
30
31 static const Float_t G1=20.03;
32 static const Float_t FDIP=-37.34;
33 static const Float_t FDIMU=6.;
34 static const Float_t FCORN=11.72;
35//
36// ZBEG Beginning of the inner triplet
37// D1BEG Beginning of separator dipole 1
38// D2BEG Beginning of separator dipole 2
39// CORBEG Corrector dipole beginning (because of dimuon arm)
40//
41 static const Float_t CORBEG=1920,COREND=CORBEG+190, CORRA2=4.5*4.5;
42//
43 static const Float_t ZBEG=2300;
44 static const Float_t Z1BEG=ZBEG+ 0,Z1END=Z1BEG+630,Z1RA2=3.5*3.5;
45 static const Float_t Z2BEG=ZBEG+ 880,Z2END=Z2BEG+550,Z2RA2=3.5*3.5;
46 static const Float_t Z3BEG=ZBEG+1530,Z3END=Z3BEG+550,Z3RA2=3.5*3.5;
47 static const Float_t Z4BEG=ZBEG+2430,Z4END=Z4BEG+630,Z4RA2=3.5*3.5;
48 static const Float_t D1BEG=5843.5 ,D1END=D1BEG+945,D1RA2=4.5*4.5;
49 static const Float_t D2BEG=12113.2 ,D2END=D2BEG+945,D2RA2=4.5*.5;
50
51//ZDC part -------------------------------------------------------------------
52
53ClassImp(AliMagF)
54
55//________________________________________
56AliMagF::AliMagF(const char *name, const char *title, const Int_t integ, const Int_t map,
57 const Float_t factor, const Float_t fmax)
58 : TNamed(name,title)
59{
60 fMap = map;
61 fType = Undef;
62 fInteg = integ;
63 fFactor = factor;
64 fMax = fmax;
65}
66
67//________________________________________
68void AliMagF::Field(Float_t*, Float_t *b)
69{
70 printf("Undefined MagF Field called, returning 0\n");
71 b[0]=b[1]=b[2]=0;
72}
73
74ClassImp(AliMagFC)
75
76//________________________________________
77AliMagFC::AliMagFC(const char *name, const char *title, const Int_t integ, const Int_t map,
78 const Float_t factor, const Float_t fmax)
79 : AliMagF(name,title,integ,map,factor,fmax)
80{
81 printf("Constant Field %s created: map= %d, factor= %f\n",fName.Data(),map,factor);
82 fType = Const;
83}
84
85//________________________________________
86void AliMagFC::Field(Float_t *x, Float_t *b)
87{
88 b[0]=b[1]=b[2]=0;
89 if(fMap==1) {
90 if(TMath::Abs(x[2])<700 && x[0]*x[0]+(x[1]+30)*(x[1]+30) < 560*560) {
91 b[2]=2;
92 } else {
93 if ( 725 <= x[2] && x[2] <= 1225 ) {
94 Float_t dz = TMath::Abs(975-x[2])*0.01;
95 b[0]=(1-0.1*dz*dz)*7;
96 }
97 else {
98//This is the ZDC part
99 Float_t rad2=x[0]*x[0]+x[1]*x[1];
100 if(rad2<D2RA2) {
101 if(x[2]>D2BEG) {
102
103// Separator Dipole D2
104 if(x[2]<D2END) b[1]=FDIP;
105 } else if(x[2]>D1BEG) {
106
107// Separator Dipole D1
108 if(x[2]<D1END) b[1]=-FDIP;
109 }
110 if(rad2<CORRA2) {
111
112// First quadrupole of inner triplet de-focussing in x-direction
113// Inner triplet
114 if(x[2]>Z4BEG) {
115 if(x[2]<Z4END) {
116
117// 2430 <-> 3060
118 b[0]=-G1*x[1];
119 b[1]=-G1*x[0];
120 }
121 } else if(x[2]>Z3BEG) {
122 if(x[2]<Z3END) {
123
124// 1530 <-> 2080
125 b[0]=G1*x[1];
126 b[1]=G1*x[0];
127 }
128 } else if(x[2]>Z2BEG) {
129 if(x[2]<Z2END) {
130
131// 890 <-> 1430
132 b[0]=G1*x[1];
133 b[1]=G1*x[0];
134 }
135 } else if(x[2]>Z1BEG) {
136 if(x[2]<Z1END) {
137
138// 0 <-> 630
139 b[0]=-G1*x[1];
140 b[1]=-G1*x[0];
141 }
142 } else if(x[2]>CORBEG) {
143 if(x[2]<COREND) {
144// Corrector dipole (because of dimuon arm)
145 b[0]=FCORN;
146 }
147 }
148 }
149 }
150 }
151 }
dd045843 152 if(fFactor!=1) {
153 b[0]*=fFactor;
154 b[1]*=fFactor;
155 b[2]*=fFactor;
156 }
fe4da5cc 157 } else {
158 printf("Invalid field map for constant field %d\n",fMap);
159 exit(1);
160 }
161}
162
163ClassImp(AliMagFCM)
164
165//________________________________________
166AliMagFCM::AliMagFCM(const char *name, const char *title, const Int_t integ, const Int_t map,
167 const Float_t factor, const Float_t fmax)
168 : AliMagF(name,title,integ,map,factor,fmax)
169{
170 fType = ConMesh;
171 printf("Constant Mesh Field %s created: map= %d, factor= %f, file= %s\n",fName.Data(),map,factor,fTitle.Data());
172}
173
174//________________________________________
175void AliMagFCM::Field(Float_t *x, Float_t *b)
176{
177 Double_t ratx, raty, ratz, hix, hiy, hiz, ratx1, raty1, ratz1,
178 bhyhz, bhylz, blyhz, blylz, bhz, blz, xl[3];
179 const Double_t one=1;
180 Int_t ix, iy, iz;
181
182 // --- find the position in the grid ---
183
184 b[0]=b[1]=b[2]=0;
185 if(-700<x[2] && x[2]<fZbeg && x[0]*x[0]+(x[1]+30)*(x[1]+30) < 560*560) {
186 b[2]=2;
187 } else {
11f6d36c 188 Bool_t infield=(fZbeg<=x[2] && x[2]<fZbeg+fZdel*(fZn-1)
fe4da5cc 189 && ( fXbeg <= TMath::Abs(x[0]) && TMath::Abs(x[0]) < fXbeg+fXdel*(fXn-1) )
11f6d36c 190 && ( fYbeg <= TMath::Abs(x[1]) && TMath::Abs(x[1]) < fYbeg+fYdel*(fYn-1) ));
191 if(infield) {
fe4da5cc 192 xl[0]=TMath::Abs(x[0])-fXbeg;
193 xl[1]=TMath::Abs(x[1])-fYbeg;
194 xl[2]=x[2]-fZbeg;
195
196 // --- start with x
197
198 hix=xl[0]*fXdeli;
199 ratx=hix-int(hix);
200 ix=int(hix);
201
202 hiy=xl[1]*fYdeli;
203 raty=hiy-int(hiy);
204 iy=int(hiy);
205
206 hiz=xl[2]*fZdeli;
207 ratz=hiz-int(hiz);
208 iz=int(hiz);
209
210 if(fMap==2) {
211 // ... simple interpolation
212 ratx1=one-ratx;
213 raty1=one-raty;
214 ratz1=one-ratz;
215 bhyhz = Bx(ix ,iy+1,iz+1)*ratx1+Bx(ix+1,iy+1,iz+1)*ratx;
216 bhylz = Bx(ix ,iy+1,iz )*ratx1+Bx(ix+1,iy+1,iz )*ratx;
217 blyhz = Bx(ix ,iy ,iz+1)*ratx1+Bx(ix+1,iy ,iz+1)*ratx;
218 blylz = Bx(ix ,iy ,iz )*ratx1+Bx(ix+1,iy ,iz )*ratx;
219 bhz = blyhz *raty1+bhyhz *raty;
220 blz = blylz *raty1+bhylz *raty;
221 b[0] = blz *ratz1+bhz *ratz;
222 //
223 bhyhz = By(ix ,iy+1,iz+1)*ratx1+By(ix+1,iy+1,iz+1)*ratx;
224 bhylz = By(ix ,iy+1,iz )*ratx1+By(ix+1,iy+1,iz )*ratx;
225 blyhz = By(ix ,iy ,iz+1)*ratx1+By(ix+1,iy ,iz+1)*ratx;
226 blylz = By(ix ,iy ,iz )*ratx1+By(ix+1,iy ,iz )*ratx;
227 bhz = blyhz *raty1+bhyhz *raty;
228 blz = blylz *raty1+bhylz *raty;
229 b[1] = blz *ratz1+bhz *ratz;
230 //
231 bhyhz = Bz(ix ,iy+1,iz+1)*ratx1+Bz(ix+1,iy+1,iz+1)*ratx;
232 bhylz = Bz(ix ,iy+1,iz )*ratx1+Bz(ix+1,iy+1,iz )*ratx;
233 blyhz = Bz(ix ,iy ,iz+1)*ratx1+Bz(ix+1,iy ,iz+1)*ratx;
234 blylz = Bz(ix ,iy ,iz )*ratx1+Bz(ix+1,iy ,iz )*ratx;
235 bhz = blyhz *raty1+bhyhz *raty;
236 blz = blylz *raty1+bhylz *raty;
237 b[2] = blz *ratz1+bhz *ratz;
238 //printf("ratx,raty,ratz,b[0],b[1],b[2] %f %f %f %f %f %f\n",
239 //ratx,raty,ratz,b[0],b[1],b[2]);
240 //
241 // ... use the dipole symmetry
242 if (x[0]*x[1] < 0) b[1]=-b[1];
243 if (x[0]<0) b[2]=-b[2];
244 } else {
245 printf("Invalid field map for constant mesh %d\n",fMap);
246 }
247 } else {
248//This is the ZDC part
249 Float_t rad2=x[0]*x[0]+x[1]*x[1];
250 if(rad2<D2RA2) {
251 if(x[2]>D2BEG) {
252
253// Separator Dipole D2
254 if(x[2]<D2END) b[1]=FDIP;
255 } else if(x[2]>D1BEG) {
256
257// Separator Dipole D1
258 if(x[2]<D1END) b[1]=-FDIP;
259 }
260 if(rad2<CORRA2) {
261
262// First quadrupole of inner triplet de-focussing in x-direction
263// Inner triplet
264 if(x[2]>Z4BEG) {
265 if(x[2]<Z4END) {
266
267// 2430 <-> 3060
268 b[0]=-G1*x[1];
269 b[1]=-G1*x[0];
270 }
271 } else if(x[2]>Z3BEG) {
272 if(x[2]<Z3END) {
273
274// 1530 <-> 2080
275 b[0]=G1*x[1];
276 b[1]=G1*x[0];
277 }
278 } else if(x[2]>Z2BEG) {
279 if(x[2]<Z2END) {
280
281// 890 <-> 1430
282 b[0]=G1*x[1];
283 b[1]=G1*x[0];
284 }
285 } else if(x[2]>Z1BEG) {
286 if(x[2]<Z1END) {
287
288// 0 <-> 630
289 b[0]=-G1*x[1];
290 b[1]=-G1*x[0];
291 }
292 } else if(x[2]>CORBEG) {
293 if(x[2]<COREND) {
294// Corrector dipole (because of dimuon arm)
295 b[0]=FCORN;
296 }
297 }
298 }
299 }
300 }
301 }
dd045843 302 if(fFactor!=1) {
303 b[0]*=fFactor;
304 b[1]*=fFactor;
305 b[2]*=fFactor;
306 }
fe4da5cc 307}
308
309//________________________________________
310void AliMagFCM::ReadField()
311{
312 FILE *magfile;
313 Int_t ix, iy, iz, ipx, ipy, ipz;
314 Float_t bx, by, bz;
315 char *fname;
316 printf("Reading Magnetic Field %s from file %s\n",fName.Data(),fTitle.Data());
317 fname = gSystem->ExpandPathName(fTitle.Data());
318 magfile=fopen(fname,"r");
319 delete [] fname;
320 if (magfile) {
321 fscanf(magfile,"%d %d %d %f %f %f %f %f %f",
322 &fXn, &fYn, &fZn, &fXdel, &fYdel, &fZdel, &fXbeg, &fYbeg, &fZbeg);
323 printf("fXn %d, fYn %d, fZn %d, fXdel %f, fYdel %f, fZdel %f, fXbeg %f, fYbeg %f, fZbeg %f\n",
324 fXn, fYn, fZn, fXdel, fYdel, fZdel, fXbeg, fYbeg, fZbeg);
325 fXdeli=1./fXdel;
326 fYdeli=1./fYdel;
327 fZdeli=1./fZdel;
328 fB = new TVector(3*fXn*fYn*fZn);
329 for (iz=0; iz<fZn; iz++) {
330 ipz=iz*3*(fXn*fYn);
331 for (iy=0; iy<fYn; iy++) {
332 ipy=ipz+iy*3*fXn;
333 for (ix=0; ix<fXn; ix++) {
334 ipx=ipy+ix*3;
335 fscanf(magfile,"%f %f %f",&bz,&by,&bx);
336 (*fB)(ipx+2)=bz;
337 (*fB)(ipx+1)=by;
338 (*fB)(ipx )=bx;
339 }
340 }
341 }
342 } else {
343 printf("File %s not found !\n",fTitle.Data());
344 exit(1);
345 }
346}
347
348
349