EGSnrc C++ class library  Report PIRS-898 (2021)
Iwan Kawrakow, Ernesto Mainegra-Hing, Frederic Tessier, Reid Townson and Blake Walters
egs_vhp_geometry.h
Go to the documentation of this file.
1 /*
2 ###############################################################################
3 #
4 # EGSnrc egs++ voxelized human phantom geometry headers
5 # Copyright (C) 2015 National Research Council Canada
6 #
7 # This file is part of EGSnrc.
8 #
9 # EGSnrc is free software: you can redistribute it and/or modify it under
10 # the terms of the GNU Affero General Public License as published by the
11 # Free Software Foundation, either version 3 of the License, or (at your
12 # option) any later version.
13 #
14 # EGSnrc is distributed in the hope that it will be useful, but WITHOUT ANY
15 # WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
16 # FOR A PARTICULAR PURPOSE. See the GNU Affero General Public License for
17 # more details.
18 #
19 # You should have received a copy of the GNU Affero General Public License
20 # along with EGSnrc. If not, see <http://www.gnu.org/licenses/>.
21 #
22 ###############################################################################
23 #
24 # Author: Iwan Kawrakow, 2008
25 #
26 # Contributors: Frederic Tessier
27 # Reid Townson
28 # Hubert Ho
29 # Marc Chamberland
30 #
31 ###############################################################################
32 #
33 # A "Voxelized Human Phantom" (VHP) geometry.
34 #
35 # Although a VHP geometry is just a XYZ geometry, due to the huge number of
36 # voxels one has to be more careful with memory use for media indices, etc.
37 # That's why this separate geometry implementation. With Richard Kramer's
38 # phantom data the implementation requires ~12 MB to hold all the organ data
39 # in memory. Richard's ASCII phantom data was also transformed into binary
40 # format with the result that loading them is almost instantaneous.
41 #
42 # For now the implementation only provides the macro geometry, i.e., no
43 # spongiosa micro matrices. Adding the micro matrices is the next step.
44 #
45 ###############################################################################
46 */
47 
48 
55 #ifndef EGS_ND_GEOMETRY_
56 #define EGS_ND_GEOMETRY_
57 
58 #ifdef WIN32
59 
60  #ifdef BUILD_VHP_DLL
61  #define EGS_VHP_EXPORT __declspec(dllexport)
62  #else
63  #define EGS_VHP_EXPORT __declspec(dllimport)
64  #endif
65  #define EGS_VHP_LOCAL
66 
67 #else
68 
69  #ifdef HAVE_VISIBILITY
70  #define EGS_VHP_EXPORT __attribute__ ((visibility ("default")))
71  #define EGS_VHP_LOCAL __attribute__ ((visibility ("hidden")))
72  #else
73  #define EGS_VHP_EXPORT
74  #define EGS_VHP_LOCAL
75  #endif
76 
77 #endif
78 
79 
80 #include "egs_base_geometry.h"
81 #include "egs_functions.h"
82 
83 #include <iostream>
84 #include <string>
85 
86 using namespace std;
87 
88 #ifndef SKIP_DOXYGEN
93 struct EGS_VHP_LOCAL VHP_RowInfo {
94  unsigned short nbin;
95  unsigned short *pix;
96  unsigned char *org;
97  VHP_RowInfo() : nbin(0) {};
98  ~VHP_RowInfo() {
99  if (nbin > 0) {
100  delete [] pix;
101  delete [] org;
102  }
103  };
104  int getBin(int i) const {
105  int ml=0, mu=nbin;
106  while (mu - ml > 1) {
107  int mav = (ml+mu)/2;
108  if (i < pix[mav]) {
109  mu = mav;
110  }
111  else {
112  ml = mav;
113  }
114  }
115  return mu-1;
116  };
117  int getOrgan(int pixel) const {
118  return pixel < pix[0] || pixel >= pix[nbin] ? 0 : org[getBin(pixel)];
119  };
120  int getSize() const {
121  int result = sizeof(nbin) + sizeof(unsigned short *) +
122  sizeof(unsigned char *);
123  if (nbin > 0) result += nbin*sizeof(unsigned char) +
124  (nbin+1)*sizeof(unsigned short);
125  return result;
126  };
127  void readData(istream &in);
128 };
129 #endif
130 
131 #ifndef SKIP_DOXYGEN
136 struct EGS_VHP_LOCAL VHP_SliceInfo {
137  unsigned short firstr;
138  unsigned short lastr;
139  VHP_RowInfo *rinfo;
140  VHP_SliceInfo() : rinfo(0) {};
141  ~VHP_SliceInfo() {
142  if (rinfo) {
143  delete [] rinfo;
144  }
145  };
146  void getSliceDimensions(int &nx, int &ny) {
147  ny = lastr+1;
148  nx = 0;
149  int nrow = lastr - firstr + 1;
150  for (int j=0; j<nrow; ++j) {
151  int ni = rinfo[j].pix[rinfo[j].nbin];
152  if (ni > nx) {
153  nx = ni;
154  }
155  }
156  };
157  int getOrgan(int row, int pixel) const {
158  return rinfo && row >= firstr && row <= lastr ?
159  rinfo[row-firstr].getOrgan(pixel) : 0;
160  };
161  int getSize() const {
162  int result = sizeof(firstr) + sizeof(lastr) +
163  sizeof(VHP_RowInfo *);
164  if (rinfo) {
165  int nrow = lastr-firstr+1;
166  for (int j=0; j<nrow; ++j) {
167  result += rinfo[j].getSize();
168  }
169  }
170  return result;
171  };
172  void readData(istream &in);
173 };
174 #endif
175 
176 #ifndef SKIP_DOXYGEN
181 class EGS_VHP_LOCAL VHP_OrganData {
182 public:
183  VHP_OrganData(const char *fname, int slice_min=0, int slice_max=1000000);
184  ~VHP_OrganData();
185 
186  int nSlice() const {
187  return nslice;
188  };
189  int getOrgan(int islice, int i, int j) const {
190  return islice >= 0 && islice < nslice ?
191  slices[islice].getOrgan(i,j) : 0;
192  };
193  const VHP_SliceInfo &getVHP_SliceInfo(int j) const {
194  return slices[j];
195  };
196  const VHP_SliceInfo &operator[](int j) const {
197  return slices[j];
198  };
199  int getSize() const {
200  int result = sizeof(nslice) + sizeof(VHP_SliceInfo *);
201  if (nslice > 0) {
202  for (int j=0; j<nslice; ++j) {
203  result += slices[j].getSize();
204  }
205  }
206  return result;
207  };
208  void getPhantomDimensions(int &Nx, int &Ny, int &Nz) {
209  Nx = nx;
210  Ny = ny;
211  Nz = nslice;
212  };
213  void getVoxelSizes(EGS_Float &Dx, EGS_Float &Dy, EGS_Float &Dz) {
214  Dx = dx;
215  Dy = dy;
216  Dz = dz;
217  };
218  int getOrgan(int ireg) const {
219  int k = ireg/nxy;
220  if (k < 0 || k >= nslice) {
221  return 0;
222  }
223  int jj = ireg - k*nx*ny;
224  int j = jj/nx;
225  int i = jj - j*nx;
226  return getOrgan(k,j,i);
227  };
228 
229  bool isOK() const {
230  return ok;
231  };
232 
233 protected:
234  double dx, dy, dz;
235  VHP_SliceInfo *slices;
236  unsigned short nslice;
237  int nx, ny, nxy;
238  bool ok;
239 };
240 #endif
241 
242 #ifndef SKIP_DOXYGEN
247 struct EGS_VHP_LOCAL EGS_VoxelInfo {
248  unsigned short ix, iy, iz;
249 };
250 #endif
251 
252 #ifndef SKIP_DOXYGEN
257 class EGS_VHP_LOCAL EGS_VoxelGeometry {
258 
259 public:
260 
261  EGS_Float dx, dxi, xmin, xmax;
262  EGS_Float dy, dyi, ymin, ymax;
263  EGS_Float dz, dzi, zmin, zmax;
264  int nx, ny, nz, nxy;
265  int ix, iy, iz;
266 
267  static EGS_VoxelInfo *v;
268  static int nv;
269 
270  EGS_VoxelGeometry(int Nx, int Ny, int Nz,
271  EGS_Float Dx, EGS_Float Dy, EGS_Float Dz) :
272  dx(Dx), dxi(1/Dx), xmin(0), xmax(Dx*Nx),
273  dy(Dy), dyi(1/Dy), ymin(0), ymax(Dy*Ny),
274  dz(Dz), dzi(1/Dz), zmin(0), zmax(Dz*Nz),
275  nx(Nx), ny(Ny), nz(Nz), nxy(nx*ny) {};
276 
277  ~EGS_VoxelGeometry() {
278  if (v) {
279  delete [] v;
280  v = 0;
281  nv = 0;
282  }
283  };
284 
285  bool isInside(const EGS_Vector &x) {
286  return (x.x >= xmin && x.x <= xmax &&
287  x.y >= ymin && x.y <= ymax &&
288  x.z >= zmin && x.z <= zmax) ? true : false;
289  };
290 
291  void setIndeces(const EGS_Vector &x) {
292  ix = (int)(x.x*dxi);
293  if (ix >= nx) {
294  ix = nx-1;
295  }
296  iy = (int)(x.y*dyi);
297  if (iy >= ny) {
298  iy = ny-1;
299  }
300  iz = (int)(x.z*dzi);
301  if (iz >= nz) {
302  iz = nz-1;
303  }
304  };
305 
306  int isWhere(const EGS_Vector &x) {
307  if (!isInside(x)) {
308  return -1;
309  }
310  setIndeces(x);
311  return ix + iy*nx + iz*nxy;
312  };
313 
314  int isWhereFast(const EGS_Vector &x) {
315  setIndeces(x);
316  return ix + iy*nx + iz*nxy;
317  };
318 
319  int computeIntersections(int ireg, int n, const EGS_Vector &X,
320  const EGS_Vector &u, EGS_GeometryIntersections *isections) {
321  if (n < 1) {
322  return -1;
323  }
324  if (nv < n) {
325  if (nv > 0) {
326  delete [] v;
327  }
328  v = new EGS_VoxelInfo [n];
329  nv = n;
330  }
331  int ifirst = 0;
332  EGS_Float t;
333  EGS_Vector x(X);
334  if (ireg < 0) {
335  t = veryFar;
336  ireg = howfar(ireg,x,u,t);
337  if (ireg < 0) {
338  return 0;
339  }
340  isections[0].t = t;
341  isections[0].rhof = 1;
342  isections[0].ireg = -1;
343  isections[0].imed = -1;
344  ++ifirst;
345  x += u*t;
346  }
347  int iz = ireg/nxy;
348  int ir = ireg - iz*nxy;
349  int iy = ir/nx;
350  int ix = ir - iy*nx;
351  EGS_Float uxi, uyi, uzi, nextx, nexty, nextz, sx, sy, sz;
352  int dirx, diry, dirz;
353  if (u.x > 0) {
354  uxi = 1/u.x;
355  dirx = 1;
356  nextx = (dx*(ix+1)-x.x)*uxi;
357  sx = uxi*dx;
358  }
359  else if (u.x < 0) {
360  uxi = 1/u.x;
361  dirx = -1;
362  nextx = (dx*ix-x.x)*uxi;
363  sx = -uxi*dx;
364  }
365  else {
366  dirx = 1;
367  nextx = veryFar*1e3;
368  sx = veryFar*1e3;
369  }
370  if (u.y > 0) {
371  uyi = 1/u.y;
372  diry = 1;
373  nexty = (dy*(iy+1)-x.y)*uyi;
374  sy = uyi*dy;
375  }
376  else if (u.y < 0) {
377  uyi = 1/u.y;
378  diry = -1;
379  nexty = (dy*iy-x.y)*uyi;
380  sy = -uyi*dy;
381  }
382  else {
383  diry = 1;
384  nexty = veryFar*1e3;
385  sy = veryFar*1e3;
386  }
387  if (u.z > 0) {
388  uzi = 1/u.z;
389  dirz = 1;
390  nextz = (dz*(iz+1)-x.z)*uzi;
391  sz = uzi*dz;
392  }
393  else if (u.z < 0) {
394  uzi = 1/u.z;
395  dirz = -1;
396  nextz = (dz*iz-x.z)*uzi;
397  sz = -uzi*dz;
398  }
399  else {
400  dirz = 1;
401  nextz = veryFar*1e3;
402  sz = veryFar*1e3;
403  }
404  for (int j=ifirst; j<n; j++) {
405  isections[j].ireg = ireg;
406  v[j].ix = ix;
407  v[j].iy = iy;
408  v[j].iz = iz;
409  int inew;
410  if (nextx < nexty && nextx < nextz) {
411  t = nextx;
412  ix += dirx;
413  if (ix < 0 || ix >= nx) {
414  inew = -1;
415  }
416  else {
417  inew = ireg + dirx;
418  nextx += sx;
419  }
420  }
421  else if (nexty < nextz) {
422  t = nexty;
423  iy += diry;
424  if (iy < 0 || iy >= ny) {
425  inew = -1;
426  }
427  else {
428  inew = ireg + nx*diry;
429  nexty += sy;
430  }
431  }
432  else {
433  t = nextz;
434  iz += dirz;
435  if (iz < 0 || iz >= nz) {
436  inew = -1;
437  }
438  else {
439  inew = ireg + nxy*dirz;
440  nextz += sz;
441  }
442  }
443  isections[j].t = t;
444  if (inew < 0) {
445  return j+1;
446  }
447  }
448  return ireg >= 0 ? -1 : n;
449  };
450 
451  EGS_Float howfarToOutside(int ireg, const EGS_Vector &x,
452  const EGS_Vector &u) {
453  if (ireg < 0) {
454  return 0;
455  }
456  EGS_Float tx = u.x > 0 ? (xmax-x.x)/u.x :
457  u.x < 0 ? (xmin-x.x)/u.x : veryFar*1e5;
458  EGS_Float ty = u.y > 0 ? (ymax-x.y)/u.y :
459  u.y < 0 ? (ymin-x.y)/u.y : veryFar*1e5;
460  EGS_Float tz = u.z > 0 ? (zmax-x.z)/u.z :
461  u.z < 0 ? (zmin-x.z)/u.z : veryFar*1e5;
462  return tx < ty && tx < tz ? tx : ty < tz ? ty : tz;
463  };
464 
465  void setIndeces(int ireg) {
466  iz = ireg/nxy;
467  int ir = ireg - iz*nxy;
468  iy = ir/nx;
469  ix = ir - iy*nx;
470  };
471 
472  //int howfar(int ireg, const EGS_Vector &x, const EGS_Vector &u,
473  // EGS_Float &t, EGS_Vector *normal=0) {
474  // setIndeces(ireg);
475  // return howfar(ireg,ix,iy,iz,x,u,t,normal);
476  //};
477  int howfarIn(int ireg, int Ix, int Iy, int Iz, const EGS_Vector &x,
478  const EGS_Vector &u, EGS_Float &t, EGS_Vector *normal=0) {
479  ix = Ix;
480  iy = Iy;
481  iz = Iz;
482  return howfarIn(ireg,x,u,t,normal);
483  };
484 
485  int howfarIn(int ireg, const EGS_Vector &x, const EGS_Vector &u,
486  EGS_Float &t, EGS_Vector *normal=0) {
487  int ixs = ix, iys = iy;
488  int inew = ireg;
489  if (u.x > 0) {
490  EGS_Float d = (dx*(ix+1)-x.x)/u.x;
491  if (d <= t) {
492  t = d;
493  if ((++ix) < nx) {
494  inew = ireg + 1;
495  }
496  else {
497  inew = -1;
498  }
499  if (normal) {
500  *normal = EGS_Vector(-1,0,0);
501  }
502  }
503  }
504  else if (u.x < 0) {
505  EGS_Float d = (dx*ix-x.x)/u.x;
506  if (d <= t) {
507  t = d;
508  if ((--ix) >= 0) {
509  inew = ireg - 1;
510  }
511  else {
512  inew = -1;
513  }
514  if (normal) {
515  *normal = EGS_Vector(1,0,0);
516  }
517  }
518  }
519  if (u.y > 0) {
520  EGS_Float d = (dy*(iy+1)-x.y)/u.y;
521  if (d <= t) {
522  ix = ixs;
523  t = d;
524  if ((++iy) < ny) {
525  inew = ireg + nx;
526  }
527  else {
528  inew = -1;
529  }
530  if (normal) {
531  *normal = EGS_Vector(0,-1,0);
532  }
533  }
534  }
535  else if (u.y < 0) {
536  EGS_Float d = (dy*iy-x.y)/u.y;
537  if (d <= t) {
538  ix = ixs;
539  t = d;
540  if ((--iy) >= 0) {
541  inew = ireg - nx;
542  }
543  else {
544  inew = -1;
545  }
546  if (normal) {
547  *normal = EGS_Vector(0,1,0);
548  }
549  }
550  }
551  if (u.z > 0) {
552  EGS_Float d = (dz*(iz+1)-x.z)/u.z;
553  if (d <= t) {
554  ix = ixs;
555  iy = iys;
556  t = d;
557  if ((++iz) < nz) {
558  inew = ireg+nxy;
559  }
560  else {
561  inew = -1;
562  }
563  if (normal) {
564  *normal = EGS_Vector(0,0,-1);
565  }
566  }
567  }
568  else if (u.z < 0) {
569  EGS_Float d = (dz*iz-x.z)/u.z;
570  if (d <= t) {
571  ix = ixs;
572  iy = iys;
573  t = d;
574  if ((--iz) >= 0) {
575  inew = ireg-nxy;
576  }
577  else {
578  inew = -1;
579  }
580  if (normal) {
581  *normal = EGS_Vector(0,0,1);
582  }
583  }
584  }
585  return inew;
586  };
587 
588  int howfarOut(const EGS_Vector &x, const EGS_Vector &u,
589  EGS_Float &t, EGS_Vector *normal=0) {
590  EGS_Float tlong = 2*t, d;
591  if (x.x <= xmin && u.x > 0) {
592  ix = 0;
593  d = (xmin-x.x)/u.x;
594  }
595  else if (x.x >= xmax && u.x < 0) {
596  ix = nx-1;
597  d = (xmax-x.x)/u.x;
598  }
599  else {
600  d = tlong;
601  }
602  if (d <= t) {
603  EGS_Float yy = x.y + u.y*d, zz = x.z + u.z*d;
604  if (yy >= ymin && yy <= ymax && zz >= zmin && zz <= zmax) {
605  t = d;
606  iy = (int)(yy*dyi);
607  if (iy > ny-1) {
608  iy = ny-1;
609  }
610  iz = (int)(zz*dzi);
611  if (iz > nz-1) {
612  iz = nz-1;
613  }
614  if (normal) *normal = ix == 0 ? EGS_Vector(-1,0,0) :
615  EGS_Vector(1,0,0);
616  return ix + iy*nx + iz*nxy;
617  }
618  }
619  if (x.y <= ymin && u.y > 0) {
620  iy = 0;
621  d = (ymin-x.y)/u.y;
622  }
623  else if (x.y >= ymax && u.y < 0) {
624  iy = ny-1;
625  d = (ymax-x.y)/u.y;
626  }
627  else {
628  d = tlong;
629  }
630  if (d <= t) {
631  EGS_Float xx = x.x + u.x*d, zz = x.z + u.z*d;
632  if (xx >= xmin && xx <= xmax && zz >= zmin && zz <= zmax) {
633  t = d;
634  ix = (int)(xx*dxi);
635  if (ix > nx-1) {
636  ix = nx-1;
637  }
638  iz = (int)(zz*dzi);
639  if (iz > nz-1) {
640  iz = nz-1;
641  }
642  if (normal) *normal = iy == 0 ? EGS_Vector(0,-1,0) :
643  EGS_Vector(0,1,0);
644  return ix + iy*nx + iz*nxy;
645  }
646  }
647  if (x.z <= zmin && u.z > 0) {
648  iz = 0;
649  d = (zmin-x.z)/u.z;
650  }
651  else if (x.z >= zmax && u.z < 0) {
652  iz = nz-1;
653  d = (zmax-x.z)/u.z;
654  }
655  else {
656  d = tlong;
657  }
658  if (d <= t) {
659  EGS_Float xx = x.x + u.x*d, yy = x.y + u.y*d;
660  if (xx >= xmin && xx <= xmax && yy >= ymin && yy <= ymax) {
661  t = d;
662  ix = (int)(xx*dxi);
663  if (ix > nx-1) {
664  ix = nx-1;
665  }
666  iy = (int)(yy*dyi);
667  if (iy > ny-1) {
668  iy = ny-1;
669  }
670  if (normal) *normal = iz == 0 ? EGS_Vector(0,0,-1) :
671  EGS_Vector(0,0,1);
672  return ix + iy*nx + iz*nxy;
673  }
674  }
675  return -1;
676  };
677 
678  int howfar(int ireg, const EGS_Vector &x, const EGS_Vector &u,
679  EGS_Float &t, EGS_Vector *normal=0) {
680  if (ireg >= 0) {
681  setIndeces(ireg);
682  return howfarIn(ireg,x,u,t,normal);
683  }
684  return howfarOut(x,u,t,normal);
685  };
686 
687  EGS_Float hownearIn(int Ix, int Iy, int Iz, const EGS_Vector &x) {
688  ix = Ix;
689  iy = Iy;
690  iz = Iz;
691  return hownearIn(x);
692  };
693 
694  EGS_Float hownearIn(const EGS_Vector &x) {
695  EGS_Float t1,t2,tx,ty,tz;
696  t1 = x.x-dx*ix;
697  t2 = dx - t1;
698  tx = t1 < t2 ? t1 : t2;
699  t1 = x.y-dy*iy;
700  t2 = dy - t1;
701  ty = t1 < t2 ? t1 : t2;
702  t1 = x.z-dz*iz;
703  t2 = dz - t1;
704  tz = t1 < t2 ? t1 : t2;
705  return tx < ty && tx < tz ? tx : ty < tz ? ty : tz;
706  };
707 
708  EGS_Float hownear(int ireg, const EGS_Vector &x) {
709  if (ireg >= 0) {
710  setIndeces(ireg);
711  return hownearIn(x);
712  }
713  int nc = 0;
714  EGS_Float s1=0, s2=0;
715  if (x.x < xmin || x.x > xmax) {
716  EGS_Float t = x.x < xmin ? xmin-x.x : x.x-xmax;
717  nc++;
718  s1 += t;
719  s2 += t*t;
720  }
721  if (x.y < ymin || x.y > ymax) {
722  EGS_Float t = x.y < ymin ? ymin-x.y : x.y-ymax;
723  nc++;
724  s1 += t;
725  s2 += t*t;
726  }
727  if (x.z < zmin || x.z > zmax) {
728  EGS_Float t = x.z < zmin ? zmin-x.z : x.z-zmax;
729  nc++;
730  s1 += t;
731  s2 += t*t;
732  }
733  return nc == 1 ? s1 : sqrt(s2);
734  };
735 };
736 #endif
737 
738 #ifndef SKIP_DOXYGEN
743 struct VHPBox {
744  EGS_Float xmin, xmax;
745  EGS_Float ymin, ymax;
746  EGS_Float zmin, zmax;
747  bool isInside(const EGS_Vector &x) {
748  return
749  (x.x>=xmin&&x.x<=xmax&&x.y>=ymin&&ymax<=ymin&&x.z>=zmin&&x.z<=zmax);
750  };
751  bool willEnter(const EGS_Vector &x, const EGS_Vector &u, EGS_Float &t) {
752  EGS_Float tlong = 2*t, d;
753  if (x.x <= xmin && u.x > 0) {
754  d = (xmin-x.x)/u.x;
755  }
756  else if (x.x >= xmax && u.x < 0) {
757  d = (xmax-x.x)/u.x;
758  }
759  else {
760  d = tlong;
761  }
762  if (d <= t) {
763  EGS_Float yy = x.y + u.y*d, zz = x.z + u.z*d;
764  if (yy >= ymin && yy <= ymax && zz >= zmin && zz <= zmax) {
765  t = d;
766  return true;
767  }
768  }
769  if (x.y <= ymin && u.y > 0) {
770  d = (ymin-x.y)/u.y;
771  }
772  else if (x.y >= ymax && u.y < 0) {
773  d = (ymax-x.y)/u.y;
774  }
775  else {
776  d = tlong;
777  }
778  if (d <= t) {
779  EGS_Float xx = x.x + u.x*d, zz = x.z + u.z*d;
780  if (xx >= xmin && xx <= xmax && zz >= zmin && zz <= zmax) {
781  t = d;
782  return true;
783  }
784  }
785  if (x.z <= zmin && u.z > 0) {
786  d = (zmin-x.z)/u.z;
787  }
788  else if (x.z >= zmax && u.z < 0) {
789  d = (zmax-x.z)/u.z;
790  }
791  else {
792  d = tlong;
793  }
794  if (d <= t) {
795  EGS_Float xx = x.x + u.x*d, yy = x.y + u.y*d;
796  if (xx >= xmin && xx <= xmax && yy >= ymin && yy <= ymax) {
797  t = d;
798  return true;
799  }
800  }
801  return false;
802  };
803  EGS_Float howfarToOut(const EGS_Vector &x, const EGS_Vector &u) {
804  EGS_Float t = veryFar*1e5, t1;
805  if (u.x > 0) {
806  t1 = (xmax-x.x)/u.x;
807  if (t1 < t) {
808  t = t1;
809  }
810  }
811  else if (u.x < 0) {
812  t1 = (xmin-x.x)/u.x;
813  if (t1 < t) {
814  t = t1;
815  }
816  }
817  if (u.y > 0) {
818  t1 = (ymax-x.y)/u.y;
819  if (t1 < t) {
820  t = t1;
821  }
822  }
823  else if (u.y < 0) {
824  t1 = (ymin-x.y)/u.y;
825  if (t1 < t) {
826  t = t1;
827  }
828  }
829  if (u.z > 0) {
830  t1 = (zmax-x.z)/u.z;
831  if (t1 < t) {
832  t = t1;
833  }
834  }
835  else if (u.z < 0) {
836  t1 = (zmin-x.z)/u.z;
837  if (t1 < t) {
838  t = t1;
839  }
840  }
841  return t;
842  };
843 };
844 #endif
845 
846 #ifndef SKIP_DOXYGEN
851 class EGS_VHP_LOCAL EGS_MicroMatrixCluster {
852 
853 public:
854 
855  EGS_MicroMatrixCluster(EGS_Float Dx, EGS_Float Dy, EGS_Float Dz,
856  EGS_Float bsc_thickness, int tb_med, int bm_med,
857  const char *micro_file);
858 
859  ~EGS_MicroMatrixCluster();
860 
861 
862  int howfar(int imic, int ireg, int &lax, int &lay, int &laz, int &newmic,
863  const EGS_Vector &x, const EGS_Vector &u,
864  EGS_Float &t, int *newmed=0, EGS_Vector *normal=0) {
865  int inew = vg->howfar(ireg,x,u,t,normal);
866  newmic = imic;
867  if (inew < 0) {
868  int iz = imic/(mx*my);
869  int itmp = imic - iz*mx*my;
870  int iy = itmp/mx;
871  int ix = itmp - mx*iy;
872  if (vg->ix < 0) {
873  --lax;
874  if ((--ix) < 0) {
875  ix = mx-1;
876  }
877  }
878  else if (vg->ix >= vg->nx) {
879  ++lax;
880  if ((++ix) >= mx) {
881  ix = 0;
882  }
883  }
884  else if (vg->iy < 0) {
885  --lay;
886  if ((--iy) < 0) {
887  iy = my-1;
888  }
889  }
890  else if (vg->iy >= vg->ny) {
891  ++lay;
892  if ((++iy) >= my) {
893  iy = 0;
894  }
895  }
896  else if (vg->iz < 0) {
897  --laz;
898  if ((--iz) < 0) {
899  iz = mz-1;
900  }
901  }
902  else if (vg->iz >= vg->nz) {
903  ++laz;
904  if ((++iz) >= mz) {
905  iz = 0;
906  }
907  }
908  newmic = ix + iy*mx + iz*mx*my;
909  }
910  if (newmed) {
911  *newmed = micros[newmic][inew] ? med_bm : med_tb;
912  }
913  return inew;
914  };
915 
916  EGS_Float hownear(int ireg, const EGS_Vector &x) {
917  return vg->hownear(ireg,x);
918  };
919 
920  int medium(int imic, int ireg) const {
921  return micros[imic][ireg] ? med_bm : med_tb;
922  };
923 
924  int medium(int ireg) const {
925  int imic = ireg/nreg;
926  ireg -= imic*nreg;
927  return medium(imic,ireg);
928  };
929 
930  int isWhere(int &imic, const EGS_Vector &x) {
931  int ix = (int)(x.x*dxi), iy = (int)(x.y*dyi), iz = (int)(x.z*dzi);
932  EGS_Vector x1(x - EGS_Vector(dx*ix, dy*iy, dz*iz));
933  ix = ix%mx;
934  iy = iy%my;
935  iz = iz%mz;
936  imic = ix + iy*mx + iz*mz;
937  return vg->isWhere(x1);
938  };
939 
940  void getIndeces(int ireg, int &imic, int &ix, int &iy, int &iz) {
941  imic = ireg/nreg;
942  ireg -= imic*nreg;
943  iz = ireg/nxy;
944  iy = (ireg - iz*nxy)/nx;
945  ix = ireg - iz*nxy - iy*nx;
946  };
947 
948  int isWhere(int jx, int jy, int jz, const EGS_Vector &x1, int *newmed=0) {
949  int ix = jx%mx, iy = jy%my, iz = jz%mz;
950  int imic = ix + iy*mx + iz*mz;
951  int iloc = vg->isWhereFast(x1);
952  if (newmed) {
953  *newmed = medium(imic,iloc);
954  }
955  return imic*nreg + iloc;
956  };
957 
958  bool isValid() const {
959  return (vg && micros);
960  };
961 
962  void getStepFractions(int imic, int ireg,
963  const EGS_Vector &xi, const EGS_Vector &xf,
964  EGS_Float &bsc_step, EGS_Float &bm_step) {
965  bsc_step = 0;
966  bm_step = 0;
967  int mict = micros[imic][ireg];
968  if (mict) {
969  // a BM voxel => check step fractions
970  if (mict == 1) {
971  bsc_step = 1; // all BSC
972  }
973  else if (mict == 2) {
974  bm_step = 1; // all BM
975  }
976  else {
977  mict -= 2;
978  // get positions in bm_boxes co-ordinate system
979  int jz = ireg/nxy;
980  int jy = (ireg-jz*nxy)/nx;
981  int jx = ireg-jz*nxy-jy*nx;
982  EGS_Vector xcorner(vg->dx*jx,vg->dy*jy,vg->dz*jz);
983  EGS_Vector Xi(xi-xcorner), Xf(xf-xcorner);
984  bool i_inside = bm_boxes[mict].isInside(Xi),
985  f_inside = bm_boxes[mict].isInside(Xf);
986  // if both positions are inside the BM box, the entire step
987  // is in the BM box
988  if (i_inside && f_inside) {
989  bm_step = 1;
990  }
991  else {
992  // we have a fraction in BM and remainder in BSC
993  EGS_Vector u(Xf-Xi);
994  EGS_Float t = u.length();
995  EGS_Float ti = 1/t;
996  u *= ti;
997  if (!i_inside) {
998  // initial position outside of BM box => find entry
999  EGS_Float tt = t;
1000  if (!bm_boxes[mict].willEnter(Xi,u,tt)) {
1001  // line Xi->Xf does not intersect BM box
1002  bsc_step = 1;
1003  return;
1004  }
1005  // move to entry point
1006  bsc_step = tt;
1007  Xi += u*tt;
1008  }
1009  if (f_inside) {
1010  // final position inside BM box. This implies initial
1011  // position was outside. Because we already computed
1012  // distance to BM box entry, the remainder goes to BM
1013  bsc_step *= ti;
1014  bm_step = 1 - bsc_step;
1015  return;
1016  }
1017  // final position outside. this means that either the
1018  // initial position was inside or that we moved to the
1019  // entry point => BM fraction is distance to outside
1020  bm_step = ti*bm_boxes[mict].howfarToOut(Xi,u);
1021  bsc_step = 1 - bm_step;
1022  }
1023  }
1024  }
1025  };
1026 
1027  void getVolumes(int imic, EGS_Float &tb, EGS_Float &bm, EGS_Float &bsc) {
1028  tb = v_tb[imic];
1029  bm = v_bm[imic];
1030  bsc = v_bsc[imic];
1031  };
1032 
1033 
1034  // ********************************************************************
1035 
1036  int nx, ny, nz, nxy, nreg;
1037  int mx, my, mz;
1038  int mxy, nmic;
1039  EGS_Float dx, dy, dz;
1040  EGS_Float dxi, dyi, dzi;
1041  EGS_Float t_bsc;
1042 
1043  int med_tb;
1044  int med_bm;
1045 
1053  unsigned char **micros;
1054  EGS_VoxelGeometry *vg;
1055  EGS_Float *v_tb;
1056  EGS_Float *v_bm;
1057  EGS_Float *v_bsc;
1058 
1059  static VHPBox bm_boxes[64];
1060  static bool bm_boxes_initialized;
1061 
1062 };
1063 #endif
1064 
1065 class EGS_VHP_EXPORT EGS_TestMicro : public EGS_BaseGeometry {
1066 
1067 public:
1068 
1069  EGS_TestMicro(EGS_Float Dx, EGS_Float Dy, EGS_Float Dz, EGS_Float bsc_t,
1070  int tb_med, int bm_med,const char *micro_file, const string &Name="");
1071 
1072  ~EGS_TestMicro();
1073 
1074  int inside(const EGS_Vector &x) {
1075  return isWhere(x);
1076  };
1077  bool isInside(const EGS_Vector &x) {
1078  return vg->isInside(x);
1079  }
1080  int isWhere(const EGS_Vector &x) {
1081  int ireg = vg->isWhere(x);
1082  if (ireg < 0) {
1083  return ireg;
1084  }
1085  EGS_Vector x1(x-EGS_Vector(vg->dx*vg->ix,vg->dy*vg->iy,vg->dz*vg->iz));
1086  int ix = vg->ix%micro->mx,
1087  iy = vg->iy%micro->my,
1088  iz = vg->iz%micro->mz;
1089  int imic = ix + iy*micro->mx + iz*micro->mxy;
1090  int iloc = micro->vg->isWhere(x1);
1091  return imic*nr + iloc;
1092  };
1093 
1094  int howfar(int ireg, const EGS_Vector &x, const EGS_Vector &u,
1095  EGS_Float &t, int *newmed=0, EGS_Vector *normal=0) {
1096  if (ireg >= 0) {
1097  int voxel = ireg/nr;
1098  int ilocal = ireg - voxel*nr;
1099  vg->setIndeces(voxel);
1100  EGS_Vector x1(x-EGS_Vector(vg->dx*vg->ix,vg->dy*vg->iy,vg->dz*vg->iz));
1101  int inew = micro->vg->howfar(ilocal,x1,u,t,normal);
1102  if (inew == ilocal) {
1103  return ireg;
1104  }
1105  if (inew < 0) {
1106  if (micro->vg->ix < 0) {
1107  if ((--vg->ix) < 0) {
1108  return -1;
1109  }
1110  micro->vg->ix = micro->vg->nx-1;
1111  }
1112  else if (micro->vg->ix >= micro->vg->nx) {
1113  if ((++vg->ix) >= vg->nx) {
1114  return -1;
1115  }
1116  micro->vg->ix = 0;
1117  }
1118  else if (micro->vg->iy < 0) {
1119  if ((--vg->iy) < 0) {
1120  return -1;
1121  }
1122  micro->vg->iy = micro->vg->ny-1;
1123  }
1124  else if (micro->vg->iy >= micro->vg->ny) {
1125  if ((++vg->iy) >= vg->ny) {
1126  return -1;
1127  }
1128  micro->vg->iy = 0;
1129  }
1130  else if (micro->vg->iz < 0) {
1131  if ((--vg->iz) < 0) {
1132  return -1;
1133  }
1134  micro->vg->iz = micro->vg->nz-1;
1135  }
1136  else if (micro->vg->iz >= micro->vg->nz) {
1137  if ((++vg->iz) >= vg->nz) {
1138  return -1;
1139  }
1140  micro->vg->iz = 0;
1141  }
1142  inew = micro->vg->ix + micro->vg->iy*micro->vg->nx +
1143  micro->vg->iz*micro->vg->nxy;
1144  }
1145  voxel = vg->ix + vg->iy*vg->nx + vg->iz*vg->nxy;
1146  if (newmed) {
1147  *newmed = micro->medium(voxel,inew);
1148  }
1149  return voxel*nr + inew;
1150  }
1151  int voxel = vg->howfar(ireg,x,u,t,normal);
1152  if (voxel < 0) {
1153  return voxel;
1154  }
1155  EGS_Vector x1(x + u*t);
1156  x1 -= EGS_Vector(vg->dx*vg->ix,vg->dy*vg->iy,vg->dz*vg->iz);
1157  int ilocal = micro->vg->isWhereFast(x1);
1158  if (ilocal < 0) {
1159  egsFatal("x=(%g,%g,%g) x+u*t=(%g,%g,%g) x1=(%g,%g,%g)\n",x.x,x.y,x.z,(x+u*t).x,(x+u*t).y,(x+u*t).z,x1.x,x1.y,x1.z);
1160  }
1161  if (newmed) {
1162  *newmed = micro->medium(voxel,ilocal);
1163  }
1164  return voxel*nr + ilocal;
1165  };
1166 
1167  EGS_Float howfarToOutside(int ireg, const EGS_Vector &x,
1168  const EGS_Vector &u) {
1169  return vg->howfarToOutside(ireg,x,u);
1170  };
1171 
1172  EGS_Float hownear(int ireg, const EGS_Vector &x) {
1173  if (ireg < 0) {
1174  return vg->hownear(ireg,x);
1175  }
1176  int voxel = ireg/nr;
1177  int ilocal = ireg - voxel*nr;
1178  vg->setIndeces(voxel);
1179  EGS_Vector x1(x-EGS_Vector(vg->dx*vg->ix,vg->dy*vg->iy,vg->dz*vg->iz));
1180  return micro->vg->hownear(ilocal,x1);
1181  };
1182 
1183  int medium(int ireg) const {
1184  int voxel = ireg/nr;
1185  int ilocal = ireg - voxel*nr;
1186  return micro->medium(voxel,ilocal);
1187  };
1188 
1189  int getMaxStep() const {
1190  return micro->mx*micro->nx + micro->my*micro->ny +
1191  micro->mz*micro->nz + 1;
1192  };
1193 
1194  const string &getType() const {
1195  return type;
1196  };
1197 
1198 protected:
1199 
1200  EGS_VoxelGeometry *vg;
1201  EGS_MicroMatrixCluster *micro;
1202  int nr;
1203 
1204  void setMedia(EGS_Input *inp, int nmed, const int *med_ind);
1205 
1206  static string type;
1207 };
1208 
1304 class EGS_VHP_EXPORT EGS_VHPGeometry : public EGS_BaseGeometry {
1305 
1306 public:
1307 
1308  EGS_VHPGeometry(const char *phantom_file, const char *media_file,
1309  int slice_min=0, int slice_max=1000000, const string &Name="");
1310 
1311  ~EGS_VHPGeometry();
1312 
1313  bool isInside(const EGS_Vector &x) {
1314  return vg->isInside(x);
1315  };
1316 
1317  int isWhere(const EGS_Vector &x) {
1318  int ireg = vg->isWhere(x);
1319  if (!nmicro) {
1320  return ireg;
1321  }
1322  int iorg = organs->getOrgan(vg->iz,vg->iy,vg->ix);
1323  int mict = organ_micro[iorg];
1324  if (mict < 0) {
1325  return ireg; // i.e., no micros in this macro voxel
1326  }
1327  EGS_Vector x1(x - EGS_Vector(vg->dx*vg->ix,vg->dy*vg->iy,vg->dz*vg->iz));
1328  int micr = micros[mict]->isWhere(vg->ix,vg->iy,vg->iz,x1);
1329  return nmax*mict + micr + nmacro;
1330  };
1331 
1332  int inside(const EGS_Vector &x) {
1333  return isWhere(x);
1334  };
1335 
1336  int medium(int ireg) const {
1337  if (ireg < nmacro) {
1338  return organ_media[organs->getOrgan(ireg)];
1339  }
1340  ireg -= nmacro;
1341  int mict = ireg/nmax;
1342  int iloc = ireg - mict*nmax;
1343  return micros[mict]->medium(iloc);
1344  };
1345 
1346  int computeIntersections(int ireg, int n, const EGS_Vector &X,
1347  const EGS_Vector &u, EGS_GeometryIntersections *isections) {
1348  // fix the following
1349  if (nmicro > 0) return EGS_BaseGeometry::computeIntersections(ireg,
1350  n,X,u,isections);
1351  if (n < 1) {
1352  return -1;
1353  }
1354  int result = vg->computeIntersections(ireg,n,X,u,isections);
1355  if (!result) {
1356  return result; // no intersections
1357  }
1358  int nloop = result > 0 ? result : n;
1359  int first = isections[0].ireg >= 0 ? 0 : 1;
1360  EGS_VoxelInfo *v = vg->v;
1361  for (int j=first; j<nloop; ++j) isections[j].imed =
1362  organ_media[organs->getOrgan(v[j].iz,v[j].iy,v[j].ix)];
1363  return result;
1364  };
1365 
1366  EGS_Float howfarToOutside(int ireg, const EGS_Vector &x,
1367  const EGS_Vector &u) {
1368  return vg->howfarToOutside(ireg,x,u);
1369  };
1370 
1371  int howfar(int ireg, const EGS_Vector &x, const EGS_Vector &u,
1372  EGS_Float &t, int *newmed=0, EGS_Vector *normal=0) {
1373  if (nmicro < 1) { // no micro matrices
1374  int inew = vg->howfar(ireg,x,u,t,normal);
1375  if (newmed && inew >= 0 && inew != ireg) {
1376  *newmed = organ_media[organs->getOrgan(vg->iz,vg->iy,vg->ix)];
1377  }
1378  return inew;
1379  }
1380  if (ireg >= 0 && ireg < nmacro) {
1381  // in macro matrix
1382  int inew = vg->howfar(ireg,x,u,t,normal);
1383  if (t < 0) {
1384  if (t < -boundaryTolerance) {
1385  egsWarning("negative step %g in macro voxel\n",t);
1386  egsWarning("x=(%g,%g,%g) u=(%g,%g,%g)\n",x.x,x.y,x.z,
1387  u.x,u.y,u.z);
1388  egsWarning("ix=%d iy=%d iz=%d\n",vg->ix,vg->iy,vg->iz);
1389  egsFatal("quitting\n");
1390  }
1391  t = boundaryTolerance;
1392  }
1393  if (inew < 0 || inew == ireg) {
1394  return inew;
1395  }
1396  // i.e., exits macro geometry or stays in same macro voxel
1397  int iorg = organs->getOrgan(vg->iz,vg->iy,vg->ix);
1398  int mict = organ_micro[iorg];
1399  if (mict < 0) { // no micro matrix in the new macro voxel
1400  if (newmed) {
1401  *newmed = organ_media[iorg];
1402  }
1403  return inew;
1404  }
1405  // if here, enters micro matrix.
1406  EGS_Vector x1(x + u*t -
1407  EGS_Vector(vg->dx*vg->ix,vg->dy*vg->iy,vg->dz*vg->iz));
1408  int micr = micros[mict]->isWhere(vg->ix,vg->iy,vg->iz,x1,newmed);
1409  return nmax*mict + micr + nmacro;
1410  }
1411  if (ireg >= 0) {
1412  // in a micro matrix.
1413  int ireg1 = ireg - nmacro;
1414  int mict = ireg1/nmax;
1415  int iloc = ireg1 - mict*nmax;
1416  //return micros[mict]->medium(iloc);
1417  int imic, ix, iy, iz;
1418  EGS_MicroMatrixCluster *micro = micros[mict];
1419  micro->getIndeces(iloc,imic,ix,iy,iz);
1420  int ilocal = iloc - imic*micros[mict]->nreg;
1421  // now we have to figure out the macro voxel indeces
1422  // in principle one could simply use vg->isWhere(x), but
1423  // this is bound to get us in trouble with roundoff errors.
1424  // we therefore treat the special case where x is in an edge
1425  // voxel of the micro matrix in a special way
1426  int lax, lay, laz;
1427  EGS_Float aux = x.x*vg->dxi;
1428  if (ix > 0 && ix < micros[mict]->vg->nx-1) {
1429  lax = (int) aux;
1430  }
1431  else if (ix == 0) {
1432  lax = (int)(aux+0.5);
1433  }
1434  else {
1435  lax = (int)(aux-0.5);
1436  }
1437  aux = x.y*vg->dyi;
1438  if (iy > 0 && iy < micros[mict]->vg->ny-1) {
1439  lay = (int) aux;
1440  }
1441  else if (iy == 0) {
1442  lay = (int)(aux+0.5);
1443  }
1444  else {
1445  lay = (int)(aux-0.5);
1446  }
1447  aux = x.z*vg->dzi;
1448  if (iz > 0 && iz < micros[mict]->vg->nz-1) {
1449  laz = (int) aux;
1450  }
1451  else if (iz == 0) {
1452  laz = (int)(aux+0.5);
1453  }
1454  else {
1455  laz = (int)(aux-0.5);
1456  }
1457  EGS_Vector x1(x - EGS_Vector(vg->dx*lax,vg->dy*lay,vg->dz*laz));
1458  int inew = micro->vg->howfarIn(ilocal,ix,iy,iz,x1,u,t,normal);
1459  if (t < 0) {
1460  if (t < -boundaryTolerance) {
1461  egsWarning("negative step %g in micro matrix\n",t);
1462  egsWarning("x=(%g,%g,%g) u=(%g,%g,%g)\n",x.x,x.y,x.z,
1463  u.x,u.y,u.z);
1464  egsWarning("ix=%d iy=%d iz=%d\n",ix,iy,iz);
1465  egsWarning("lax=%d lay=%d laz=%d\n",lax,lay,laz);
1466  egsWarning("x1=(%g,%g,%g)\n",x1.x,x1.y,x1.z);
1467  egsFatal("quitting\n");
1468  }
1469  t = 0;
1470  }
1471  if (inew == ilocal) {
1472  return ireg;
1473  }
1474  if (inew < 0) {
1475  // exiting current micro matrix => entering new
1476  // macro voxel.
1477  // 1. find which one and set new micro region
1478  if (micro->vg->ix < 0) {
1479  if ((--lax) < 0) {
1480  return -1;
1481  }
1482  micro->vg->ix = micro->vg->nx-1;
1483  }
1484  else if (micro->vg->ix >= micro->vg->nx) {
1485  if ((++lax) >= vg->nx) {
1486  return -1;
1487  }
1488  micro->vg->ix = 0;
1489  }
1490  else if (micro->vg->iy < 0) {
1491  if ((--lay) < 0) {
1492  return -1;
1493  }
1494  micro->vg->iy = micro->vg->ny-1;
1495  }
1496  else if (micro->vg->iy >= micro->vg->ny) {
1497  if ((++lay) >= vg->ny) {
1498  return -1;
1499  }
1500  micro->vg->iy = 0;
1501  }
1502  else if (micro->vg->iz < 0) {
1503  if ((--laz) < 0) {
1504  return -1;
1505  }
1506  micro->vg->iz = micro->vg->nz-1;
1507  }
1508  else if (micro->vg->iz >= micro->vg->nz) {
1509  if ((++laz) >= vg->nz) {
1510  return -1;
1511  }
1512  micro->vg->iz = 0;
1513  }
1514  // 2. micro matrix in new macro voxel?
1515  int iorg = organs->getOrgan(laz,lay,lax);
1516  mict = organ_micro[iorg];
1517  if (mict < 0) { // no.
1518  if (newmed) {
1519  *newmed = organ_media[iorg];
1520  }
1521  return lax + lay*vg->nx + laz*vg->nxy;
1522  }
1523  inew = micro->vg->ix + micro->vg->iy*micro->vg->nx +
1524  micro->vg->iz*micro->vg->nxy;
1525  int lix = lax%micros[mict]->mx;
1526  int liy = lay%micros[mict]->my;
1527  int liz = laz%micros[mict]->mz;
1528  imic = lix + micros[mict]->mx*liy + micros[mict]->mxy*liz;
1529  }
1530  if (newmed) {
1531  *newmed = micros[mict]->medium(imic,inew);
1532  }
1533  return nmax*mict + inew + nmacro;
1534  }
1535  // outside
1536  int imac = vg->howfarOut(x,u,t,normal);
1537  if (imac < 0) {
1538  return imac;
1539  }
1540  int iorg = organs->getOrgan(vg->iz,vg->iy,vg->ix);
1541  int mict = organ_micro[iorg];
1542  if (mict < 0) {
1543  if (newmed) {
1544  *newmed = organ_media[iorg];
1545  }
1546  return imac;
1547  }
1548  EGS_Vector x1(x + u*t);
1549  x1 -= EGS_Vector(vg->dx*vg->ix,vg->dy*vg->iy,vg->dz*vg->iz);
1550  int ilocal = micros[mict]->vg->isWhereFast(x1);
1551  int ix=vg->ix%micros[mict]->mx,
1552  iy=vg->iy%micros[mict]->my,
1553  iz=vg->iz%micros[mict]->mz;
1554  int imic = ix + iy*micros[mict]->mx + iz*micros[mict]->mxy;
1555  if (newmed) {
1556  *newmed = micros[mict]->medium(imic,ilocal);
1557  }
1558  return nmax*mict + imic*micros[mict]->nreg + ilocal + nmacro;
1559  };
1560 
1561  EGS_Float hownear(int ireg, const EGS_Vector &x) {
1562  if (ireg < nmacro) {
1563  return vg->hownear(ireg,x);
1564  }
1565  ireg -= nmacro;
1566  int mict = ireg/nmax;
1567  int iloc = ireg - mict*nmax;
1568  int imic, ix, iy, iz;
1569  EGS_MicroMatrixCluster *micro = micros[mict];
1570  micro->getIndeces(iloc,imic,ix,iy,iz);
1571  // now we have to figure out the macro voxel indeces
1572  // in principle one could simply use vg->isWhere(x), but
1573  // this is bound to get us in trouble with roundoff errors.
1574  // we therefore treat the special case where x is in an edge
1575  // voxel of the micro matrix in a special way
1576  int lax, lay, laz;
1577  EGS_Float aux = x.x*vg->dxi;
1578  if (ix > 0 && ix < micro->vg->nx-1) {
1579  lax = (int) aux;
1580  }
1581  else if (ix == 0) {
1582  lax = (int)(aux+0.5);
1583  }
1584  else {
1585  lax = (int)(aux-0.5);
1586  }
1587  aux = x.y*vg->dyi;
1588  if (iy > 0 && iy < micro->vg->ny-1) {
1589  lay = (int) aux;
1590  }
1591  else if (iy == 0) {
1592  lay = (int)(aux+0.5);
1593  }
1594  else {
1595  lay = (int)(aux-0.5);
1596  }
1597  aux = x.z*vg->dzi;
1598  if (iz > 0 && iz < micro->vg->nz-1) {
1599  laz = (int) aux;
1600  }
1601  else if (iz == 0) {
1602  laz = (int)(aux+0.5);
1603  }
1604  else {
1605  laz = (int)(aux-0.5);
1606  }
1607  EGS_Vector x1(x - EGS_Vector(vg->dx*lax,vg->dy*lay,vg->dz*laz));
1608  return micro->vg->hownearIn(ix,iy,iz,x1);
1609  };
1610 
1611  int getMaxStep() const {
1612  if (!nmicro) {
1613  return vg->nx+vg->ny+vg->nz+1;
1614  }
1615  int nmx=0, nmy=0, nmz=0;
1616  for (int mict=0; mict<nmicro; ++mict) {
1617  if (micros[mict]->nx > nmx) {
1618  nmx = micros[mict]->nx;
1619  }
1620  if (micros[mict]->ny > nmy) {
1621  nmy = micros[mict]->ny;
1622  }
1623  if (micros[mict]->nx > nmz) {
1624  nmx = micros[mict]->nz;
1625  }
1626  }
1627  return vg->nx*nmx+vg->ny*nmy+vg->nz*nmz+1;
1628  };
1629 
1630  bool isOK() const {
1631  return (vg && organs);
1632  };
1633 
1634  int getNx() const {
1635  return vg->nx;
1636  };
1637  int getNy() const {
1638  return vg->ny;
1639  };
1640  int getNz() const {
1641  return vg->nz;
1642  };
1643 
1644  const string &getType() const {
1645  return type;
1646  };
1647 
1648  void printInfo() const;
1649 
1650  //static EGS_VHPGeometry* createGeometry(EGS_Input *);
1651  void setMicros(EGS_Input *);
1652 
1653 protected:
1654 
1655  EGS_VoxelGeometry *vg;
1656  VHP_OrganData *organs;
1657 
1658  int organ_media[256];
1659  int organ_micro[256];
1660  string organ_names[256];
1661  EGS_MicroMatrixCluster **micros;
1662  int nmicro;
1663  int nmax, nmacro;
1664 
1665  static string type;
1666 
1667  void setup();
1668 
1669  void setMedia(EGS_Input *inp, int nmed, const int *med_ind);
1670 
1671 };
1672 
1673 #endif
Base geometry class. Every geometry class must be derived from EGS_BaseGeometry.
virtual int computeIntersections(int ireg, int n, const EGS_Vector &x, const EGS_Vector &u, EGS_GeometryIntersections *isections)
Calculates intersection distances to region boundaries.
A class for storing information in a tree-like structure of key-value pairs. This class is used throu...
Definition: egs_input.h:182
int nr
number of voxels in a micro-matrix
A Voxelized Human Phantom (VHP) geometry.
A class representing 3D vectors.
Definition: egs_vector.h:57
EGS_Float y
y-component
Definition: egs_vector.h:62
EGS_Float z
z-component
Definition: egs_vector.h:63
EGS_Float x
x-component
Definition: egs_vector.h:61
EGS_BaseGeometry class header file.
Global egspp functions header file.
EGS_InfoFunction EGS_EXPORT egsFatal
Always use this function for reporting fatal errors.
EGS_InfoFunction EGS_EXPORT egsWarning
Always use this function for reporting warnings.
const EGS_Float veryFar
A very large float.
EGS_I32 ireg
region index
EGS_Float t
distance to next region boundary
EGS_Float rhof
relative mass density in that region
EGS_I32 imed
medium index