FairRoot/PandaRoot
PndFieldMap.cxx
Go to the documentation of this file.
1 // -------------------------------------------------------------------------
2 // ----- PndFieldMap source file -----
3 // ----- Created 12/01/04 by M. Al/Turany (FairField.cxx) -----
4 // -------------------------------------------------------------------------
5 
6 
7 #include <iomanip>
8 #include <iostream>
9 #include <fstream>
10 #include "stdlib.h"
11 
12 #include "TArrayF.h"
13 #include "TFile.h"
14 #include "TMath.h"
15 
16 #include "PndFieldMap.h"
17 #include "PndFieldMapData.h"
18 #include "PndFieldPar.h"
19 
20 using namespace std;
21 
22 // ------------- Default constructor ----------------------------------
24  : FairField(),
25  fFileName(""),
26  fScale(1.0),
27  funit(10.0),
28  fPosX(0), fPosY(0), fPosZ(0),
29  fXmin(0), fXmax(0), fXstep(0),
30  fYmin(0), fYmax(0), fYstep(0),
31  fZmin(0), fZmax(0), fZstep(0),
32  fNx(0),fNy(0),fNz(0),
33  fBx(NULL), fBy(NULL), fBz(NULL)
34 {
35  SetName("");
36  fType = 1;
37 }
38 // ------------------------------------------------------------------------
39 
40 
41 
42 // ------------- Standard constructor ---------------------------------
43 PndFieldMap::PndFieldMap(const char* mapName, const char* fileType)
44  : FairField(mapName),
45  fFileName(TString("")),
46  fScale(1.0),
47  funit(10.0),
48  fPosX(0), fPosY(0), fPosZ(0),
49  fXmin(0), fXmax(0), fXstep(0),
50  fYmin(0), fYmax(0), fYstep(0),
51  fZmin(0), fZmax(0), fZstep(0),
52  fNx(0),fNy(0),fNz(0),
53  fBx(NULL), fBy(NULL), fBz(NULL)
54 {
55  SetName(mapName);
56  TString dir = getenv("VMCWORKDIR");
57  fFileName = dir + "/input/fieldmaps/" + mapName;
58  if ( fileType[0] == 'R' ) fFileName += ".root";
59  else fFileName += ".dat";
60  fType = 1;
61 }
62 // ------------------------------------------------------------------------
63 
64 
65 
66 // ------------ Constructor from PndFieldPar --------------------------
68  : FairField(),
69  fFileName(TString("")),
70  fScale(1.0),
71  funit(10.0),
72  fPosX(0), fPosY(0), fPosZ(0),
73  fXmin(0), fXmax(0), fXstep(0),
74  fYmin(0), fYmax(0), fYstep(0),
75  fZmin(0), fZmax(0), fZstep(0),
76  fNx(0),fNy(0),fNz(0),
77  fBx(NULL), fBy(NULL), fBz(NULL)
78 {
79  fType = 1;
80  if ( ! fieldPar ) {
81  LOG(ERROR) << "PndConstField::PndConstField: empty parameter container!" ;
82  SetName("");
83  fType = -1;
84  }
85  else {
86  TString Name=GetName();
87  fieldPar->MapName(Name);
88  fPosX = fieldPar->GetPositionX();
89  fPosY = fieldPar->GetPositionY();
90  fPosZ = fieldPar->GetPositionZ();
91  fScale = fieldPar->GetScale();
92  TString dir = getenv("VMCWORKDIR");
93  fFileName = dir + "/input/fieldmaps/" + Name + ".root";
94  fType = fieldPar->GetType();
95  }
96 }
97 // ------------------------------------------------------------------------
98 
99 // ------------- ----Copy constructor ---------------------------------
101  : FairField(),
102  fFileName(L.fFileName),
103  fScale(L.fScale),
104  funit(L.fScale),
105  fPosX(L.fPosX), fPosY(L.fPosY), fPosZ(L.fPosZ),
106  fXmin(L.fXmin), fXmax(L.fXmax), fXstep(L.fXstep),
107  fYmin(L.fYmin), fYmax(L.fYmax), fYstep(L.fYstep),
108  fZmin(L.fZmin), fZmax(L.fZmax), fZstep(L.fZstep),
109  fNx(L.fNx),fNy(L.fNy),fNz(L.fNz),
110  fBx(L.fBx), fBy(L.fBy), fBz(L.fBz)
111 {
112  fType = L.fType;
113  for (Int_t ii=0; ii<2; ii++)
114  {
115  fHc[ii] = L.fHc[ii];
116  for (Int_t jj=0; jj<2; jj++)
117  {
118  fHb[ii][jj] = L.fHb[ii][jj];
119  for (Int_t kk=0; kk<2; kk++)
120  {
121  fHa[ii][jj][kk] = L.fHa[ii][jj][kk];
122  }
123  }
124  }
125 
126 }
127 // ------------------------------------------------------------------------
128 
129 
130 
131 // ------------ Destructor --------------------------------------------
133 
134  //printf("PndFieldMap::~PndFieldMap() \n");
135  if ( fBx ) delete fBx;
136  if ( fBy ) delete fBy;
137  if ( fBz ) delete fBz;
138 }
139 // ------------------------------------------------------------------------
140 
141 
142 
143 // ----------- Intialisation ------------------------------------------
145  if (fFileName.EndsWith(".root")) ReadRootFile(fFileName, fName);
146  else if (fFileName.EndsWith(".dat")) ReadAsciiFile(fFileName);
147  else {
148  LOG(ERROR) << "-E- PndFieldMap::Init: No proper file name defined! "<<fFileName.Data() ;
149  Fatal("Init", "No proper file name");
150  }
151 }
152 // ------------------------------------------------------------------------
153 
154 
155 
156 // ----------- Get x component of the field ---------------------------
158 
159  Int_t ix = 0;
160  Int_t iy = 0;
161  Int_t iz = 0;
162  Double_t dx = 0.;
163  Double_t dy = 0.;
164  Double_t dz = 0.;
165 
166  if ( IsInside(x, y, z, ix, iy, iz, dx, dy, dz) ) {
167 
168  // Get Bx field values at grid cell corners
169  fHa[0][0][0] = fBx->At(ix *fNy*fNz + iy *fNz + iz);
170  fHa[1][0][0] = fBx->At((ix+1)*fNy*fNz + iy *fNz + iz);
171  fHa[0][1][0] = fBx->At(ix *fNy*fNz + (iy+1)*fNz + iz);
172  fHa[1][1][0] = fBx->At((ix+1)*fNy*fNz + (iy+1)*fNz + iz);
173  fHa[0][0][1] = fBx->At(ix *fNy*fNz + iy *fNz + (iz+1));
174  fHa[1][0][1] = fBx->At((ix+1)*fNy*fNz + iy *fNz + (iz+1));
175  fHa[0][1][1] = fBx->At(ix *fNy*fNz + (iy+1)*fNz + (iz+1));
176  fHa[1][1][1] = fBx->At((ix+1)*fNy*fNz + (iy+1)*fNz + (iz+1));
177 
178  // Return interpolated field value
179  return Interpolate(dx, dy, dz);
180 
181  }
182 
183  return 0.;
184 }
185 // ------------------------------------------------------------------------
186 
187 
188 
189 // ----------- Get y component of the field ---------------------------
191 
192  Int_t ix = 0;
193  Int_t iy = 0;
194  Int_t iz = 0;
195  Double_t dx = 0.;
196  Double_t dy = 0.;
197  Double_t dz = 0.;
198 
199  if ( IsInside(x, y, z, ix, iy, iz, dx, dy, dz) ) {
200 
201  // Get By field values at grid cell corners
202  fHa[0][0][0] = fBy->At(ix *fNy*fNz + iy *fNz + iz);
203  fHa[1][0][0] = fBy->At((ix+1)*fNy*fNz + iy *fNz + iz);
204  fHa[0][1][0] = fBy->At(ix *fNy*fNz + (iy+1)*fNz + iz);
205  fHa[1][1][0] = fBy->At((ix+1)*fNy*fNz + (iy+1)*fNz + iz);
206  fHa[0][0][1] = fBy->At(ix *fNy*fNz + iy *fNz + (iz+1));
207  fHa[1][0][1] = fBy->At((ix+1)*fNy*fNz + iy *fNz + (iz+1));
208  fHa[0][1][1] = fBy->At(ix *fNy*fNz + (iy+1)*fNz + (iz+1));
209  fHa[1][1][1] = fBy->At((ix+1)*fNy*fNz + (iy+1)*fNz + (iz+1));
210 
211  // Return interpolated field value
212  return Interpolate(dx, dy, dz);
213 
214  }
215 
216  return 0.;
217 }
218 // ------------------------------------------------------------------------
219 
220 
221 
222 // ----------- Get z component of the field ---------------------------
224 
225  Int_t ix = 0;
226  Int_t iy = 0;
227  Int_t iz = 0;
228  Double_t dx = 0.;
229  Double_t dy = 0.;
230  Double_t dz = 0.;
231 
232  if ( IsInside(x, y, z, ix, iy, iz, dx, dy, dz) ) {
233 
234  // Get Bz field values at grid cell corners
235  fHa[0][0][0] = fBz->At(ix *fNy*fNz + iy *fNz + iz);
236  fHa[1][0][0] = fBz->At((ix+1)*fNy*fNz + iy *fNz + iz);
237  fHa[0][1][0] = fBz->At(ix *fNy*fNz + (iy+1)*fNz + iz);
238  fHa[1][1][0] = fBz->At((ix+1)*fNy*fNz + (iy+1)*fNz + iz);
239  fHa[0][0][1] = fBz->At(ix *fNy*fNz + iy *fNz + (iz+1));
240  fHa[1][0][1] = fBz->At((ix+1)*fNy*fNz + iy *fNz + (iz+1));
241  fHa[0][1][1] = fBz->At(ix *fNy*fNz + (iy+1)*fNz + (iz+1));
242  fHa[1][1][1] = fBz->At((ix+1)*fNy*fNz + (iy+1)*fNz + (iz+1));
243 
244  // Return interpolated field value
245  return Interpolate(dx, dy, dz);
246 
247  }
248 
249  return 0.;
250 }
251 // ------------------------------------------------------------------------
252 
253 
254 
255 // ----------- Check whether a point is inside the map ----------------
257  Int_t& ix, Int_t& iy, Int_t& iz,
258  Double_t& dx, Double_t& dy, Double_t& dz) {
259 
260  // --- Transform into local coordinate system
261  Double_t xl = x - fPosX;
262  Double_t yl = y - fPosY;
263  Double_t zl = z - fPosZ;
264 
265  // --- Check for being outside the map range
266  if ( ! ( xl >= fXmin && xl < fXmax && yl >= fYmin && yl < fYmax &&
267  zl >= fZmin && zl < fZmax ) ) {
268  ix = iy = iz = 0;
269  dx = dy = dz = 0.;
270  return kFALSE;
271  }
272 
273  // --- Determine grid cell
274  ix = Int_t( xl / fXstep );
275  iy = Int_t( yl / fYstep );
276  iz = Int_t( zl / fZstep );
277 
278 
279  // Relative distance from grid point (in units of cell size)
280  dx = xl / fXstep - Double_t(ix);
281  dy = yl / fYstep - Double_t(iy);
282  dz = zl / fZstep - Double_t(iz);
283 
284  return kTRUE;
285 
286 }
287 // ------------------------------------------------------------------------
288 
289 
290 
291 // ---------- Write the map to an ASCII file --------------------------
292 void PndFieldMap::WriteAsciiFile(const char* fileName) {
293 
294  // Open file
295  LOG(INFO) <<"PndFieldMap: Writing field map to ASCII file " <<fileName ;
296  ofstream mapFile(fileName);
297  if ( ! mapFile.is_open() ) {
298  LOG(ERROR) << "PndFieldMap:ReadAsciiFile: Could not open file! " ;
299  return;
300  }
301 
302  // Write field map grid parameters
303  mapFile.precision(4);
304  mapFile << showpoint;
305  if ( fType == 1 ) mapFile << "nosym" << endl;
306  if ( fType == 2 ) mapFile << "Solenoid" << endl;
307  if ( fType == 3 ) mapFile << "Dipole" << endl;
308  if ( fType == 4 ) mapFile << "Trans" << endl;
309  if ( funit == 10.0 ) mapFile << "T" << endl;
310  if ( funit == 0.001 ) mapFile << "G" << endl;
311  if ( funit == 1.0 ) mapFile << "kG" << endl;
312 
313  mapFile << fXmin << " " << fXmax << " " << fNx << endl;
314  mapFile << fYmin << " " << fYmax << " " << fNy << endl;
315  mapFile << fZmin << " " << fZmax << " " << fNz << endl;
316 
317  // Write field values
318  Double_t factor = funit * fScale; // Takes out scaling
319  cout << right;
320  Int_t nTot = fNx * fNy * fNz;
321  cout << "-I- PndFieldMap: " << fNx*fNy*fNz << " entries to write... "
322  << setw(3) << 0 << " % ";
323  Int_t index=0;
324  div_t modul;
325  Int_t iDiv = TMath::Nint(nTot/100.);
326  for(Int_t ix=0; ix<fNx; ix++) {
327  for(Int_t iy=0; iy<fNy; iy++) {
328  for(Int_t iz=0; iz<fNz; iz++) {
329  index =ix*fNy*fNz + iy*fNz + iz;
330  if (iDiv!=0)
331  {
332  modul = div(index,iDiv);
333  if ( modul.rem == 0 ) {
334  Double_t perc = TMath::Nint(100.*index/nTot);
335  cout << "\b\b\b\b\b\b" << setw(3) << perc << " % " << flush;
336  }
337  }
338  mapFile << fBx->At(index)/factor << " " << fBy->At(index)/factor
339  << " " << fBz->At(index)/factor << endl;
340  } // z-Loop
341  } // y-Loop
342  } // x-Loop
343  cout << " " << index+1 << " written" << endl;
344  mapFile.close();
345 
346 }
347 // ------------------------------------------------------------------------
348 
349 
350 
351 // ------- Write field map to a ROOT file -----------------------------
352 void PndFieldMap::WriteRootFile(const char* fileName,
353  const char* mapName) {
354 
355  PndFieldMapData* data = new PndFieldMapData(mapName, *this);
356  TFile* oldFile = gFile;
357  TFile* file = new TFile(fileName, "RECREATE");
358  data->Write();
359  file->Close();
360  if(oldFile) oldFile->cd();
361 
362 }
363 // ------------------------------------------------------------------------
364 
365 
366 
367 // ----- Set the position of the field centre in global coordinates -----
369  fPosX = x;
370  fPosY = y;
371  fPosZ = z;
372 }
373 // ------------------------------------------------------------------------
374 
375 
376 
377 // --------- Screen output --------------------------------------------
379  TString type = "Map";
380  if ( fType == 2 ) type = "Soleniod Map ";
381  if ( fType == 3 ) type = "Dipole Map ";
382  if ( fType == 4 ) type = "Trans Map ";
383  cout << "======================================================" << endl;
384  cout.precision(4);
385  cout << showpoint;
386  cout << "---- " << fTitle << " : " << fName << endl;
387  cout << "----" << endl;
388  cout << "---- Field type : " << type << endl;
389  cout << "----" << endl;
390  cout << "---- Field map grid : " << endl;
391  cout << "---- x = " << setw(4) << fXmin << " to " << setw(4) << fXmax
392  << " cm, " << fNx << " grid points, dx = " << fXstep << " cm" << endl;
393  cout << "---- y = " << setw(4) << fYmin << " to " << setw(4) << fYmax
394  << " cm, " << fNy << " grid points, dy = " << fYstep << " cm" << endl;
395  cout << "---- z = " << setw(4) << fZmin << " to " << setw(4) << fZmax
396  << " cm, " << fNz << " grid points, dz = " << fZstep << " cm" << endl;
397  cout << endl;
398  cout << "---- Field centre position: ( " << setw(6) << fPosX << ", "
399  << setw(6) << fPosY << ", " << setw(6) << fPosZ << ") cm" << endl;
400  cout << "---- Field scaling factor: " << fScale << endl;
401  Double_t bx = GetBx(0.,0.,0.);
402  Double_t by = GetBy(0.,0.,0.);
403  Double_t bz = GetBz(0.,0.,0.);
404  cout << "----" << endl;
405  cout << "---- Field at origin is ( " << setw(6) << bx << ", " << setw(6)
406  << by << ", " << setw(6) << bz << ") kG" << endl;
407  cout << "======================================================" << endl;
408 }
409 // ------------------------------------------------------------------------
410 
411 
412 
413 // --------- Reset parameters and data (private) ----------------------
415  fPosX = fPosY = fPosZ = 0.;
416  fXmin = fYmin = fZmin = 0.;
417  fXmax = fYmax = fZmax = 0.;
418  fXstep = fYstep = fZstep = 0.;
419  fNx = fNy = fNz = 0;
420  fScale = 1.;
421  funit = 10.0;
422  if ( fBx ) { delete fBx; fBx = NULL; }
423  if ( fBy ) { delete fBy; fBy = NULL; }
424  if ( fBz ) { delete fBz; fBz = NULL; }
425 }
426 // ------------------------------------------------------------------------
427 
428 
429 
430 // ----- Read field map from ASCII file (private) ---------------------
431 void PndFieldMap::ReadAsciiFile(const char* fileName) {
432 
433  Double_t bx=0., by=0., bz=0.;
434  // Open file
435  cout << "-I- PndFieldMap: Reading field map from ASCII file "
436  << fileName << endl;
437  ifstream mapFile(fileName);
438  if ( ! mapFile.is_open() ) {
439  cerr << "-E- PndFieldMap:ReadAsciiFile: Could not open file! " << endl;
440  Fatal("ReadAsciiFile","Could not open file");
441  }
442 
443  // Read map type
444  TString type;
445  mapFile >> type;
446 
447  Int_t iType = 0;
448  if ( type == "nosym" ) iType = 1;
449  if ( type == "Solenoid") iType = 2;
450  if ( type == "Dipole" ) iType = 3;
451  if ( type == "Trans" ) iType = 4;
452  if ( fType != iType ) {
453  cout << "-E- PndFieldMap::ReadAsciiFile: Incompatible map types!"
454  << endl;
455  cout << " Field map is of type " << fType
456  << " but map on file is of type " << iType << endl;
457  Fatal("ReadAsciiFile","Incompatible map types");
458  }
459  // Read Units
460  TString unit;
461  mapFile >> unit;
462  if ( unit == "G" ) funit = 0.001;
463  else if ( unit == "T" ) funit = 10.0;
464  else if ( unit == "kG" ) funit=1.0;
465  else {
466  cout << "-E- FieldMap::ReadAsciiFile: No units!"
467  << endl;
468  Fatal("ReadAsciiFile","No units defined");
469  }
470 
471 
472  // Read grid parameters
473 
474  mapFile >>fXmin >> fXmax >> fNx;
475  mapFile >>fYmin >> fYmax >> fNy;
476  mapFile >>fZmin >> fZmax >> fNz;
477  fXstep = ( fXmax - fXmin ) / Double_t( fNx - 1 );
478  fYstep = ( fYmax - fYmin ) / Double_t( fNy - 1 );
479  fZstep = ( fZmax - fZmin ) / Double_t( fNz - 1 );
480 
481  // Create field arrays
482  fBx = new TArrayF(fNx * fNy * fNz);
483  fBy = new TArrayF(fNx * fNy * fNz);
484  fBz = new TArrayF(fNx * fNy * fNz);
485 
486  // Read the field values
487  Double_t factor = fScale * funit; // Factor 1/1000 for G -> kG
488  cout << right;
489  Int_t nTot = fNx * fNy * fNz;
490  cout << "-I- PndFieldMap: " << nTot << " entries to read... "
491  << setw(3) << 0 << " % ";
492  Int_t index = 0;
493  div_t modul;
494  Int_t iDiv = TMath::Nint(nTot/100.);
495  for (Int_t ix=0; ix<fNx; ix++) {
496  for (Int_t iy = 0; iy<fNy; iy++) {
497  for (Int_t iz = 0; iz<fNz; iz++) {
498  if (! mapFile.good()) cerr << "-E- PndFieldMap::ReadAsciiFile: "
499  << "I/O Error at " << ix << " "
500  << iy << " " << iz << endl;
501  index = ix*fNy*fNz + iy*fNz + iz;
502  if (iDiv!=0)
503  {
504  modul = div(index,iDiv);
505  if ( modul.rem == 0 ) {
506  Double_t perc = TMath::Nint(100.*index/nTot);
507  cout << "\b\b\b\b\b\b" << setw(3) << perc << " % " << flush;
508  }
509  }
510  //mapFile >> xx>>yy>>zz>> bx >> by >> bz ;
511  mapFile >> bx >> by >> bz ;
512  //cout << " x= " <<xx <<" y= " << yy<<" z= " << zz<<" bx= " << bx <<" by= " <<by <<" bz= " << bz<< endl;
513  fBx->AddAt(factor*bx, index);
514  fBy->AddAt(factor*by, index);
515  fBz->AddAt(factor*bz, index);
516  if ( mapFile.eof() ) {
517  cerr << endl << "-E- PndFieldMap::ReadAsciiFile: EOF"
518  << " reached at " << ix << " " << iy << " " << iz << endl;
519  mapFile.close();
520  break;
521  }
522  } // z-Loop
523  } // y-Loop0)
524  } // x-Loop
525 
526  cout << " " << index+1 << " read" << endl;
527 
528  mapFile.close();
529 
530 }
531 // ------------------------------------------------------------------------
532 
533 
534 
535 // ------------- Read field map from ROOT file (private) ---------------
536 void PndFieldMap::ReadRootFile(const char* fileName,
537  const char* mapName) {
538 
539  // Store gFile pointer
540  TFile* oldFile = gFile;
541 
542  // Open root file
543  LOG(INFO) <<"PndFieldMap: Reading field map from ROOT file "<<fileName ;
544  TFile* file = new TFile(fileName, "READ");
545  if (file->IsZombie()) {
546  LOG(ERROR) << "-E- PndFieldMap::ReadRootfile: Cannot read from file! " ;
547  Fatal("ReadRootFile","Cannot read from file");
548  }
549 
550  // Get the field data object
551  PndFieldMapData* data = NULL;
552  file->GetObject(mapName, data);
553  if ( ! data ) {
554  LOG(ERROR)<<"PndFieldMap::ReadRootFile: data object %s not found in file! "<< fileName ;
555  exit(-1);
556  }
557 
558  // Get the field parameters
559  SetField(data);
560 
561  // Close the root file and delete the data object
562  file->Close();
563  delete data;
564  delete file;
565  if ( oldFile ) oldFile->cd();
566 
567 }
568 // ------------------------------------------------------------------------
569 
570 
571 
572 // ------------ Set field parameters and data (private) ----------------
574 
575  // Check compatibility
576  if ( data->GetType() != fType ) {
577  LOG(ERROR)<<"PndFieldMap::SetField: Incompatible map types Field map is of type "<<fType<<" \n but map on file is of type "<<data->GetType() ;
578  Fatal("SetField","Incompatible map types");
579  }
580 
581 
582  fXmin = data->GetXmin();
583  fYmin = data->GetYmin();
584  fZmin = data->GetZmin();
585  fXmax = data->GetXmax();
586  fYmax = data->GetYmax();
587  fZmax = data->GetZmax();
588  fNx = data->GetNx();
589  fNy = data->GetNy();
590  fNz = data->GetNz();
591  fXstep = ( fXmax - fXmin ) / Double_t( fNx - 1 );
592  fYstep = ( fYmax - fYmin ) / Double_t( fNy - 1 );
593  fZstep = ( fZmax - fZmin ) / Double_t( fNz - 1 );
594  if ( fBx ) delete fBx;
595  if ( fBy ) delete fBy;
596  if ( fBz ) delete fBz;
597  fBx = new TArrayF(*(data->GetBx()));
598  fBy = new TArrayF(*(data->GetBy()));
599  fBz = new TArrayF(*(data->GetBz()));
600 
601  // Scale and convert from G(or T) to kG
602  Double_t factor = fScale * funit;
603  Int_t index = 0;
604  for (Int_t ix=0; ix<fNx; ix++) {
605  for (Int_t iy=0; iy<fNy; iy++) {
606  for (Int_t iz=0; iz<fNz; iz++) {
607  index = ix*fNy*fNz + iy*fNz + iz;
608  if ( fBx ) (*fBx)[index] = (*fBx)[index] * factor;
609  if ( fBy ) (*fBy)[index] = (*fBy)[index] * factor;
610  if ( fBz ) (*fBz)[index] = (*fBz)[index] * factor;
611  }
612  }
613  }
614 
615 }
616 // ------------------------------------------------------------------------
617 
618 
619 
620 // ------------ Interpolation in a grid cell (private) -----------------
622 
623  // Interpolate in x coordinate
624  fHb[0][0] = fHa[0][0][0] + ( fHa[1][0][0]-fHa[0][0][0] ) * dx;
625  fHb[1][0] = fHa[0][1][0] + ( fHa[1][1][0]-fHa[0][1][0] ) * dx;
626  fHb[0][1] = fHa[0][0][1] + ( fHa[1][0][1]-fHa[0][0][1] ) * dx;
627  fHb[1][1] = fHa[0][1][1] + ( fHa[1][1][1]-fHa[0][1][1] ) * dx;
628 
629  // Interpolate in y coordinate
630  fHc[0] = fHb[0][0] + ( fHb[1][0] - fHb[0][0] ) * dy;
631  fHc[1] = fHb[0][1] + ( fHb[1][1] - fHb[0][1] ) * dy;
632 
633  // Interpolate in z coordinate
634  return fHc[0] + ( fHc[1] - fHc[0] ) * dz;
635 
636 }
637 // ------------------------------------------------------------------------
638 
639 
640 
Double_t fHa[2][2][2]
Definition: PndFieldMap.h:189
Double_t funit
Definition: PndFieldMap.h:164
double dy
virtual void Init()
Int_t GetType() const
Definition: PndFieldPar.h:37
Double_t fYmin
Definition: PndFieldMap.h:173
cout<< "-----------------------------------------------> Quarter VOLUME<<endl;name="QuarterShape";QuarterShape=newTGeoArb8(name,dz,vertQuar);name="Quarter4Vol";TStringmedium="air";QuarterVol=newTGeoVolumeAssembly(name);name="SubunitShape";SubunitShape=newTGeoArb8(name,dz,vertSub);TStringmedium="air";name="SubunitVol";name1="SubunitVol1";name2="SubunitVol2";name3="SubunitVol3";name4="SubunitVol4";name5="SubunitVol5";name6="SubunitVol6";name7="SubunitVol7";name8="SubunitVol8";name9="SubunitVol9";SubunitVol=newTGeoVolumeAssembly(name);SubunitVol1=newTGeoVolumeAssembly(name1);SubunitVol2=newTGeoVolumeAssembly(name2);SubunitVol3=newTGeoVolumeAssembly(name3);SubunitVol4=newTGeoVolumeAssembly(name4);SubunitVol5=newTGeoVolumeAssembly(name5);SubunitVol6=newTGeoVolumeAssembly(name6);SubunitVol7=newTGeoVolumeAssembly(name7);SubunitVol8=newTGeoVolumeAssembly(name8);SubunitVol9=newTGeoVolumeAssembly(name9);name="BoxShape";BoxShape=newTGeoArb8(name,dz,vertBox);TStringmedium="air";name="BoxVol";BoxVol=newTGeoVolumeAssembly(name);name1="BoxVol1";name2="BoxVol2";name3="BoxVol3";name4="BoxVol4";BoxVol1=newTGeoVolumeAssembly(name1);BoxVol2=newTGeoVolumeAssembly(name2);BoxVol3=newTGeoVolumeAssembly(name3);BoxVol4=newTGeoVolumeAssembly(name4);for(Int_tb=0;b<kNumOfBoxes;b++){cout<<""<<endl;cout<<"---------------->BOXnumber:"<<b<<endl;if(b==0){trBox=newTGeoTranslation(tr+kSpaceInSub+kAlveoleThickness+0.5*kSpaceInBox,tr+kSpaceInSub+kAlveoleThickness+0.5*kSpaceInBox,0.);rotBox=newTGeoRotation();rotBox.RotateX(0.465518);rotBox.RotateY(-0.465518);}if(b==1){trBox=newTGeoTranslation(-tr-kSpaceInSub-kAlveoleThickness-0.5*kSpaceInBox,-tr-kSpaceInSub-kAlveoleThickness-0.5*kSpaceInBox,0.);rotBox=newTGeoRotation();rotBox.RotateX(-0.465518);rotBox.RotateY(0.465518);}if(b==2){trBox=newTGeoTranslation(tr+kSpaceInSub+kAlveoleThickness+0.5*kSpaceInBox,-tr-kSpaceInSub-kAlveoleThickness-0.5*kSpaceInBox,0.);rotBox=newTGeoRotation();rotBox.RotateX(-0.465518);rotBox.RotateY(-0.465518);}if(b==3){trBox=newTGeoTranslation(-tr-kSpaceInSub-kAlveoleThickness-0.5*kSpaceInBox,tr+kSpaceInSub+kAlveoleThickness+0.5*kSpaceInBox,0.);rotBox=newTGeoRotation();rotBox.RotateX(0.465518);rotBox.RotateY(0.465518);}TGeoCombiTrans*trrotBox=newTGeoCombiTrans(trBox,rotBox);name="BoxVol";name+=b;trrotBox->SetName(name);trrotBox->RegisterYourself();SubunitVol->AddNode(BoxVol,b,trrotBox);if(b==1){name+=b;trrotBox->SetName(name);trrotBox->RegisterYourself();SubunitVol1->AddNode(BoxVol,b,trrotBox);}if(b==2){name+=b;trrotBox->SetName(name);trrotBox->RegisterYourself();SubunitVol2->AddNode(BoxVol1,b,trrotBox);}if(b==0){name+=b;trrotBox-> SetName(name)
Double_t GetYmin() const
void WriteAsciiFile(const char *fileName)
Double_t GetZmin() const
TFile * file
Double_t GetXmax() const
exit(0)
void MapName(TString &name)
Definition: PndFieldPar.h:47
Double_t GetPositionZ() const
Definition: PndFieldPar.h:50
Double_t fHb[2][2]
Field at corners of a grid cell.
Definition: PndFieldMap.h:190
Double_t fZstep
Definition: PndFieldMap.h:174
void WriteRootFile(const char *fileName, const char *mapName)
Int_t GetNy() const
Double_t fZmin
Definition: PndFieldMap.h:174
Double_t GetZmax() const
virtual Bool_t IsInside(Double_t x, Double_t y, Double_t z, Int_t &ix, Int_t &iy, Int_t &iz, Double_t &dx, Double_t &dy, Double_t &dz)
Double_t fYmax
Definition: PndFieldMap.h:173
TArrayF * fBz
Definition: PndFieldMap.h:184
Double_t fXstep
Definition: PndFieldMap.h:172
Double_t fXmax
Definition: PndFieldMap.h:172
TString fFileName
Definition: PndFieldMap.h:157
Double_t fHc[2]
Interpolated field (2-dim)
Definition: PndFieldMap.h:191
Double_t fPosZ
Definition: PndFieldMap.h:168
TArrayF * GetBz() const
Double_t fScale
Definition: PndFieldMap.h:161
Double_t GetPositionY() const
Definition: PndFieldPar.h:49
void SetPosition(Double_t x, Double_t y, Double_t z)
Double_t Interpolate(Double_t dx, Double_t dy, Double_t dz)
Double_t
TArrayF * GetBz() const
Definition: PndFieldMap.h:116
Double_t GetXmin() const
PndMultiFieldPar * fieldPar
Definition: sim_ftof.C:102
TArrayF * GetBy() const
Definition: PndFieldMap.h:115
Double_t z
Double_t fPosY
Definition: PndFieldMap.h:168
TArrayF * fBy
Definition: PndFieldMap.h:183
virtual ~PndFieldMap()
Int_t GetNx() const
double dx
void SetField(const PndFieldMapData *data)
TArrayF * fBx
Definition: PndFieldMap.h:182
Double_t GetYmax() const
Double_t fPosX
Definition: PndFieldMap.h:168
Double_t x
Double_t fYstep
Definition: PndFieldMap.h:173
Int_t GetNz() const
ClassImp(PndAnaContFact)
TArrayF * GetBy() const
TArrayF * GetBx() const
Double_t fZmax
Definition: PndFieldMap.h:174
Double_t GetPositionX() const
Definition: PndFieldPar.h:48
TArrayF * GetBx() const
Definition: PndFieldMap.h:114
Double_t y
Double_t GetScale() const
Definition: PndFieldPar.h:51
Int_t GetType() const
virtual void Print()
void ReadAsciiFile(const char *fileName)
int Nint(float x)
Definition: PndCAMath.h:117
void ReadRootFile(const char *fileName, const char *mapName)
Double_t fXmin
Definition: PndFieldMap.h:172