GEMC  2.3
Geant4 Monte-Carlo Framework
phi-segmented.cc
Go to the documentation of this file.
1 // gemc headers
2 #include "asciiField.h"
3 #include "utils.h"
4 
5 
6 // phi-segmented 3D field in cylindrical coordinates. Field itself is in cartesian coordinates.
7 // Dependent on 3 cartesian coordinates (transverse, azimuthal, longitudinal) expressed in cylindrical coordinate
8 // The values can be loaded from the map in any order as long as their speed is ordered.
9 // The values are indexed as B1_3D[transverse][longi]
10 // The field is two dimensional, ordered in the class as B1=BT, B2=BL
11 
12 // from fieldFactory: load field map
14 {
15  setlocale(LC_NUMERIC, "en_US");
16 
17  int np_1 = map->getCoordinateWithSpeed(0).np;
18  int np_2 = map->getCoordinateWithSpeed(1).np;
19  int np_3 = map->getCoordinateWithSpeed(2).np;
20  double min1 = map->getCoordinateWithSpeed(0).min;
21  double min2 = map->getCoordinateWithSpeed(1).min;
22  double min3 = map->getCoordinateWithSpeed(2).min;
23  double cell1 = (map->getCoordinateWithSpeed(0).max - min1)/(np_1 - 1);
24  double cell2 = (map->getCoordinateWithSpeed(1).max - min2)/(np_2 - 1);
25  double cell3 = (map->getCoordinateWithSpeed(2).max - min3)/(np_3 - 1);
26 
27  // Allocate memory. [AZI][TRANSVERSE][LONGI]
28  // as initialized in the map
29  map->B1_3D = new double**[map->np[0]];
30  map->B2_3D = new double**[map->np[0]];
31  map->B3_3D = new double**[map->np[0]];
32  for (unsigned i = 0; i < map->np[0]; ++i)
33  {
34  map->B1_3D[i] = new double*[map->np[1]];
35  map->B2_3D[i] = new double*[map->np[1]];
36  map->B3_3D[i] = new double*[map->np[1]];
37  for (unsigned j = 0; j < map->np[1]; ++j)
38  {
39  map->B1_3D[i][j] = new double[map->np[2]];
40  map->B2_3D[i][j] = new double[map->np[2]];
41  map->B3_3D[i][j] = new double[map->np[2]];
42 
43  }
44  }
45 
46  double unit1 = get_number("1*" + map->getCoordinateWithSpeed(0).unit);
47  double unit2 = get_number("1*" + map->getCoordinateWithSpeed(1).unit);
48  double unit3 = get_number("1*" + map->getCoordinateWithSpeed(2).unit);
49 
50  double scale = map->scaleFactor*get_number("1*" + map->unit);
51 
52  double d1, d2, d3, b1, b2, b3;
53 
54  // progress bar
55  int barWidth = 50;
56 
57  // using fscanf instead of c++ to read file, this is a lot faster
58  char ctmp[100];
59  string tmp;
60  FILE *fp = fopen (map->identifier.c_str(), "r");
61 
62  // ignoring header
63  while(tmp != "</mfield>")
64  {
65  fscanf(fp, "%s", ctmp);
66  tmp = string(ctmp);
67  }
68 
69  // now reading map
70  // values as read from map
71  for(int i1 = 0; i1<np_1 ; i1++)
72  {
73  for(int i2 = 0; i2<np_2 ; i2++)
74  {
75  for(int i3 = 0; i3<np_3 ; i3++)
76  {
77  fscanf(fp, "%lg %lg %lg %lg %lg %lg", &d1, &d2, &d3, &b1, &b2, &b3);
78 
79  d1 *= unit1;
80  d2 *= unit2;
81  d3 *= unit3;
82  b1 *= scale;
83  b2 *= scale;
84  b3 *= scale;
85 
86  // checking map consistency for first coordinate
87  if( (min1 + i1*cell1 - d1)/d1 > 0.001)
88  {
89  cout << " !! Error: coordinate index wrong. Map point should be " << min1 + i1*cell1
90  << " but it's " << d1 << " instead." << endl;
91  }
92  // checking map consistency for second coordinate
93  if( (min2 + i2*cell2 - d2)/d2 > 0.001)
94  {
95  cout << " !! Error: coordinate index wrong. Map point should be " << min2 + i2*cell2
96  << " but it's " << d2 << " instead." << endl;
97  }
98 
99  // checking map consistency for third coordinate
100  if( (min3 + i3*cell3 - d3)/d3 > 0.001)
101  {
102  cout << " !! Error: coordinate index wrong. Map point should be " << min2 + i2*cell2
103  << " but it's " << d2 << " instead." << endl;
104  }
105 
106  if(verbosity>4)
107  cout << " Loading Map: coordinates (" << d1 << ", " << d2 << ", " << d3 << ") values: (" << b1 << ", " << b2 << ", " << b3 << ")." << endl;
108 
109  // calculating index
110  unsigned t1 = (unsigned) floor( ( d1 - min1 + cell1/2 ) / ( cell1 ) ) ;
111  unsigned t2 = (unsigned) floor( ( d2 - min2 + cell2/2 ) / ( cell2 ) ) ;
112  unsigned t3 = (unsigned) floor( ( d3 - min3 + cell3/2 ) / ( cell3 ) ) ;
113 
114  // The values are indexed as B1_3D[AZI][TRANSVERSE][LONGI]
115  if( map->getCoordinateWithSpeed(0).name == "azimuthal"
116  && map->getCoordinateWithSpeed(1).name == "transverse"
117  && map->getCoordinateWithSpeed(2).name == "longitudinal" )
118  {
119  map->B1_3D[t1][t2][t3] = b1;
120  map->B2_3D[t1][t2][t3] = b2;
121  map->B3_3D[t1][t2][t3] = b3;
122  }
123  if( map->getCoordinateWithSpeed(0).name == "azimuthal"
124  && map->getCoordinateWithSpeed(1).name == "longitudinal"
125  && map->getCoordinateWithSpeed(2).name == "transverse" )
126  {
127  map->B1_3D[t1][t3][t2] = b1;
128  map->B2_3D[t1][t3][t2] = b2;
129  map->B3_3D[t1][t3][t2] = b3;
130  }
131  if( map->getCoordinateWithSpeed(0).name == "longitudinal"
132  && map->getCoordinateWithSpeed(1).name == "azimuthal"
133  && map->getCoordinateWithSpeed(2).name == "transverse" )
134  {
135  map->B1_3D[t3][t1][t2] = b1;
136  map->B2_3D[t3][t1][t2] = b2;
137  map->B3_3D[t3][t1][t2] = b3;
138  }
139  if( map->getCoordinateWithSpeed(0).name == "longitudinal"
140  && map->getCoordinateWithSpeed(1).name == "transverse"
141  && map->getCoordinateWithSpeed(2).name == "azimuthal" )
142  {
143  map->B1_3D[t3][t2][t1] = b1;
144  map->B2_3D[t3][t2][t1] = b2;
145  map->B3_3D[t3][t2][t1] = b3;
146  }
147  if( map->getCoordinateWithSpeed(0).name == "transverse"
148  && map->getCoordinateWithSpeed(1).name == "longitudinal"
149  && map->getCoordinateWithSpeed(2).name == "azimuthal" )
150  {
151  map->B1_3D[t2][t3][t1] = b1;
152  map->B2_3D[t2][t3][t1] = b2;
153  map->B3_3D[t2][t3][t1] = b3;
154  }
155  if( map->getCoordinateWithSpeed(0).name == "transverse"
156  && map->getCoordinateWithSpeed(1).name == "azimuthal"
157  && map->getCoordinateWithSpeed(2).name == "longitudinal" )
158  {
159  map->B1_3D[t2][t1][t3] = b1;
160  map->B2_3D[t2][t1][t3] = b2;
161  map->B3_3D[t2][t1][t3] = b3;
162  }
163 
164  }
165  }
166 
167  cout << " [";
168  double progress = (double)i1/(double)np_1;
169  int pos = progress*barWidth;
170  for (int i = 0; i < barWidth; ++i)
171  {
172  if (i < pos) cout << "=";
173  else if (i == pos) cout << ">";
174  else cout << " ";
175  }
176  cout << "] " << int(progress * 100.0) << " %\r";
177  cout.flush();
178  }
179  fclose(fp);
180 
181  cout << endl;
182 }
183 
184 // from mappedField: GetFieldValue
185 void gMappedField::GetFieldValue_phiSegmented( const double x[3], double *Bfield, int FIRST_ONLY) const
186 {
187 
188  // this routine will rotate the coordinate to segment local
189  // coordinates to find the map indexes.
190  // Then it will rotate the field back to the original position
191 
192  double mfield[3] = {0,0,0};
193 
194  double aC = 0; // azimuthal
195  double tC = 0; // transverse
196  double lC = 0; // longitudinal
197 
198  double aLC;
199  unsigned aI, tI, lI;
200 
201 
202  aC = atan2( x[1], x[0] )*rad;
203  if( aC < 0 ) aC += 360*deg;
204 
205  tC = sqrt(x[0]*x[0] + x[1]*x[1]);
206  lC = x[2];
207 
208  // Rotating the point to within the map limit
209  int segment = 1;
210  aLC = aC;
211  while (aLC/deg > 30)
212  {
213  aLC -= 60*deg;
214  segment++;
215  }
216 
217  // need the delta phi from the local coordinate, so we can rotate the field back to the original point
218  double dphi = aC - aLC;
219 
220  // for the map index we need the absolute value
221  double aaLC = fabs(aLC);
222 
223  // map indexes, bottom of the cell
224  aI = floor( ( aaLC - startMap[0] ) / cellSize[0] );
225  tI = floor( ( tC - startMap[1] ) / cellSize[1] );
226  lI = floor( ( lC - startMap[2] ) / cellSize[2] );
227 
228  // checking if the point is closer to the top of the cell
229  if( fabs( startMap[0] + aI*cellSize[0] - aaLC) > fabs( startMap[0] + (aI + 1)*cellSize[0] - aaLC) ) aI++;
230  if( fabs( startMap[1] + tI*cellSize[1] - tC) > fabs( startMap[1] + (tI + 1)*cellSize[1] - tC) ) tI++;
231  if( fabs( startMap[2] + lI*cellSize[2] - lC) > fabs( startMap[2] + (lI + 1)*cellSize[2] - lC) ) lI++;
232 
233  // outside map, returning no field
234  if(aI >= np[0] || tI >= np[1] || lI >= np[2]) return;
235 
236 
237  // positive on the right side of the segment
238  int sign = (aLC >= 0 ? 1 : -1);
239 
240 
241  // no interpolation
242  if(interpolation == "none")
243  {
244  // Field at local point
245  mfield[0] = B1_3D[aI][tI][lI];
246  mfield[1] = B2_3D[aI][tI][lI];
247  mfield[2] = B3_3D[aI][tI][lI];
248 
249 
250  // Rotating the field back to original point
251  Bfield[0] = sign*mfield[0] * cos(dphi/rad) - mfield[1] * sin(dphi/rad);
252  Bfield[1] = sign*mfield[0] * sin(dphi/rad) + mfield[1] * cos(dphi/rad);
253  Bfield[2] = sign*mfield[2];
254 
255  if(verbosity>3 && FIRST_ONLY != 99)
256  {
257 
258  cout << " > Track position in magnetic field: "
259  << "(" << (x[0] + mapOrigin[0])/cm << ", "
260  << (x[1] + mapOrigin[1])/cm << ", "
261  << (x[2] + mapOrigin[2])/cm << ") cm, " << endl;
262 
263  cout << " > Track position in magnetic field (phi, r, z): "
264  << "(" << aC/deg << " deg, "
265  << tC/cm << " cm, "
266  << lC/cm << " cm) " << endl;
267 
268  cout << " > Local coordinates (phi, segment, delta phi) "
269  << "(" << aaLC/deg << " deg, "
270  << segment << ", "
271  << dphi/deg << " deg) " << endl;
272 
273  cout << " > Map Indexes (phi, r, z) "
274  << "(" << aI << ", "
275  << tI << ", "
276  << lI << ") " << endl;
277 
278  cout << " > Field Values (local) (Bx, By, Bz) "
279  << "(" << mfield[0]/tesla << ", "
280  << mfield[1]/tesla << ", "
281  << mfield[2]/tesla << ") tesla " << endl;
282 
283  cout << " > Field Values (absolute) (Bx, By, Bz) "
284  << "(" << Bfield[0]/tesla << ", "
285  << Bfield[1]/tesla << ", "
286  << Bfield[2]/tesla << ") tesla " << endl;
287 
288  cout << endl;
289  }
290  }
291 }
292 
293 
294 
295 
296 
297 
unsigned int np
Definition: mappedField.h:40
string unit
field unit in the map
Definition: mappedField.h:97
double verbosity
string name
Definition: mappedField.h:39
unsigned int * np
Definition: mappedField.h:116
double *** B2_3D
Definition: mappedField.h:103
double scaleFactor
copy of the gfield scaleFactor
Definition: mappedField.h:96
double min
Definition: mappedField.h:41
double get_number(string v, int warn_no_unit)
Returns number with dimension from string, i.e. 100*cm.
double *** B3_3D
Definition: mappedField.h:104
double max
Definition: mappedField.h:42
void loadFieldMap_phiSegmented(gMappedField *, double)
gcoord getCoordinateWithSpeed(int speed)
return coordinate based on speed
Definition: mappedField.cc:38
string unit
Definition: mappedField.h:43
void GetFieldValue_phiSegmented(const double x[3], double *Bfield, int FIRST_ONLY) const
double *** B1_3D
Definition: mappedField.h:102
string identifier
Pointer to map in factory (for example, hostname / filename with path / date)
Definition: mappedField.h:92