cf8b6a09ff80573e945533e2101b6d7c1bf0b597
[ctsim.git] / include / ir.h
1 /*****************************************************************************
2 **  This is part of the CTSim program
3 **  Copyright (C) 1983-2000 Kevin Rosenberg
4 **
5 **  $Id: ir.h,v 1.16 2000/05/24 22:48:17 kevin Exp $
6 **  $Log: ir.h,v $
7 **  Revision 1.16  2000/05/24 22:48:17  kevin
8 **  First functional version of SDF library for X-window
9 **
10 **  Revision 1.15  2000/05/16 04:33:17  kevin
11 **  Updated documentation
12 **
13 **  Revision 1.14  2000/05/11 14:07:00  kevin
14 **  Added support for Windows NT
15 **
16 **  Revision 1.13  2000/05/08 20:00:48  kevin
17 **  ANSI C changes
18 **
19 **  Revision 1.12  2000/05/07 12:46:19  kevin
20 **  made c++ compatible
21 **
22 **  Revision 1.11  2000/05/05 02:37:31  kevin
23 **  renamed phmelm to pelm
24 **
25 **  Revision 1.10  2000/05/04 18:16:34  kevin
26 **  renamed filter definitions
27 **
28 **  Revision 1.9  2000/05/04 04:29:18  kevin
29 **  *** empty log message ***
30 **
31 **  Revision 1.8  2000/05/04 04:25:55  kevin
32 **  Renamed phantom and phantom-element functions/variables
33 **
34 **  Revision 1.7  2000/05/03 19:51:41  kevin
35 **  function renaming for phantoms and phantom elements
36 **
37 **  Revision 1.6  2000/05/03 08:49:49  kevin
38 **  Code cleanup
39 **
40 **  Revision 1.5  2000/05/02 20:00:25  kevin
41 **  *** empty log message ***
42 **
43 **  Revision 1.4  2000/05/02 15:31:39  kevin
44 **  code cleaning
45 **
46 **  Revision 1.3  2000/04/29 23:24:29  kevin
47 **  *** empty log message ***
48 **
49 **  Revision 1.2  2000/04/28 14:14:16  kevin
50 **  *** empty log message ***
51 **
52 **  Revision 1.1.1.1  2000/04/28 13:02:43  kevin
53 **  Initial CVS import for first public release
54 **
55 **
56 **
57 **  This program is free software; you can redistribute it and/or modify
58 **  it under the terms of the GNU General Public License (version 2) as
59 **  published by the Free Software Foundation.
60 **
61 **  This program is distributed in the hope that it will be useful,
62 **  but WITHOUT ANY WARRANTY; without even the implied warranty of
63 **  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
64 **  GNU General Public License for more details.
65 **
66 **  You should have received a copy of the GNU General Public License
67 **  along with this program; if not, write to the Free Software
68 **  Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
69 ******************************************************************************/
70
71 /*      ir.h        Header File for Image Reconstruction System
72  *      Programmer:   Kevin Rosenberg
73  *      Date Started: 7-1-84
74  */
75
76 #ifndef IR_H
77 #define IR_H
78
79 #ifdef __cplusplus
80 extern "C" {
81 #endif /* __cplusplus */
82
83 #ifdef MPI_CT
84 #define MPI_MAX_PROCESS 128 
85 struct mpi_ct_st
86 {
87   int my_rank;
88   int nproc;
89   int base_local_work_units;
90   int remainder_work_units;
91   int local_work_units[MPI_MAX_PROCESS];
92   int start_work_unit[MPI_MAX_PROCESS];
93   MPI_Comm comm;
94 };
95
96 extern struct mpi_ct_st mpi_ct;
97 void mpi_ct_calc_work_units(const unsigned int global_work_units);
98 #endif
99
100
101 struct histo_st {
102     int *b;                     /* Histogram array (# of elements in each bin) */
103     int nbin;                   /* Number of histogram bins */
104     double xmin, xmax, xinc;    /* Limits of histogram boundaries */
105 };
106
107 typedef struct histo_st HISTOGRAM;
108
109 /*---------------------------------------------------------------------------*/
110
111 static const int POINTS_PER_CIRCLE=90;
112 #define MAXREMARK  99
113
114 typedef FMTX_2D         IMAGE_ARRAY;    /* use 2d floating point matrix */
115 typedef float           IMAGE_ELEM_VAL; /* use floats for image storing */
116
117 struct image_st {
118     IMAGE_ARRAY v;                      /* values of voxels in matrix form */
119     SDF_2D *dfp_2d;                     /* Pointer to disk image file */
120     int nx, ny;                         /* size of voxel matrix */
121     double xmin, xmax, ymin, ymax;      /* extent of voxel matrix in phm coord */
122     char remark[MAXREMARK+1];           /* description of voxel data */
123     float calctime;                     /* time to calculate voxels in seconds */
124 };
125 typedef struct image_st  IMAGE;
126
127 typedef enum {
128   RECTANGLE,
129   TRIANGLE,
130   ELLIPSE,
131   SECTOR,
132   SEGMENT
133 } PElmType;
134
135 struct pelm_st {
136     PElmType type;           /* pelm type (box, ellipse, etc) */
137     double atten;            /* X-ray attenuation coefficient */
138     double cx,cy;            /* center of pelm */
139     double u,v;              /* size of pelm */
140     double rot;              /* pelm rotation angle (in radians) */
141     double *x, *y;           /* ptr to array of points in obj world coord */
142     int pts;                 /* number of points in outline arrays */
143     double xmin, xmax, ymin, ymax; /* pelm limits */
144     double radius;           /*   "   */
145     GRFMTX_2D p_to_o;        /* map from phantom to standard pelm coords */
146     GRFMTX_2D o_to_p;        /* map from std pelm coords to phantom coords */
147     struct pelm_st *next;    /* pointer to next pelm in phantom */
148 };
149 typedef struct pelm_st PELM;
150
151 typedef enum {
152   P_PELMS,        /* Phantom made of Pelms */
153   P_UNIT_PULSE,   /* Special phantom, not made of pelms */
154   P_FILTER        /* defined only by this type */
155 } PhmType;
156
157 struct phm_st {                    /* Phantom structure */
158     PELM *pelm_list;               /* pelm linked-list */
159     PhmType type;
160     int n_pelm;                    /* number of pelms in phantom */
161     double xmin, xmax, ymin, ymax; /* extent of pelms in pelm coordinates */
162     double radius;                 /*      " "    */
163 };
164 typedef struct phm_st    PHANTOM;
165
166
167 /*----------------------------------------------------------------------*/
168 /*                              RAYSUM SYMBOLS                          */
169 /*----------------------------------------------------------------------*/
170
171 /* Ray sums are collected along an array of ndet detectors.  The data
172  * for these detectors is stored in the structure DETECTARRAY
173  */
174
175 typedef float DETECT_TYPE;
176
177 struct detarray_st {
178   DETECT_TYPE *detval;  /* Pointer to array of values recorded by detector */
179   int ndet;             /* Number of detectors in array */
180   double view_angle;    /* View angle in radians */
181 };
182 typedef struct detarray_st DETARRAY;
183
184
185 typedef enum {
186     DETECTOR_PARALLEL,
187     DETECTOR_EQUIANGLE,
188     DETECTOR_EQUILINEAR
189 } ScannerGeometry;
190   
191 struct detector_st {
192   ScannerGeometry geometry;     /* Geometry of detectory */
193   int ndet;                     /* Number of detectors in array */
194   int nview;                    /* Number of rotated views */
195   int nsample;                  /* Number of rays per detector */
196   double detlen;                /* Total length of detector array */
197   double rotlen;                /* Rotation angle length in radians (norm 2PI) */
198   double det_inc;               /* Increment between centers of detectors */
199   double rot_inc;               /* Increment in rotation angle between views */
200   double radius;                /* Radius of rotation.  Distance from */
201                                 /*   center of phm to center of det */
202   double phmlen;                /* Maximum Length of phantom or area of interest */
203   struct {
204     double xd1,yd1,xd2,yd2;     /* Coordinates of detector endpoints */
205     double xs1,ys1,xs2,ys2;     /* Coordinates of source endpoints */
206     double angle;               /* Starting angle */
207   } init;
208 };
209 typedef struct detector_st DETECTOR;
210
211 struct raysum_st {
212   int fd;                       /* Disk file descriptor */
213   int file_mode;                /* Current file mode (read or write) */
214   int header_size;              /* Size of disk file header */
215   int geometry;                 /* Geometry of scanner */
216   struct detarray_st **view;    /* Pointer to array of detarray_st pointers */
217
218   char remark[MAXREMARK+1];     /* description of raysum data */
219   double calctime;              /* time required to calculate raysums */
220
221   int ndet;                     /* number of detectors in array */
222   int nview;                    /* number of rotated views */
223   double rot_start;             /* starting view rotation */
224   double rot_inc;               /* angle between rotations */
225   double det_start;             /* distance of beginning detector to center */
226                                 /*    of PHANTOM */
227   double det_inc;               /* increment between detectors */
228   double phmlen;                /* Length of PHANTOM edge (phm is square) */
229 };
230 typedef struct raysum_st   RAYSUM;
231
232 /*----------------------------------------------------------------------*/
233 /*                              USER SYMBOLS                            */
234 /*----------------------------------------------------------------------*/
235
236 /* Codes for Coordinate Types      */
237 /* Defines coords for pelm_is_point_inside() */
238
239 typedef enum {
240   PELM_COORD,         /* Normalized Pelm Coordinates */
241   PHM_COORD           /* User's phantom Coordinates */
242 } CoordType;
243
244 /* Codes for Filter types */
245
246 typedef enum {  /* filter types for filter_generate() */
247   FILTER_BANDLIMIT, 
248   FILTER_SINC,
249   FILTER_G_HAMMING,
250   FILTER_COSINE,
251   FILTER_TRIANGLE,
252   FILTER_ABS_BANDLIMIT,         /* filter times |x| */
253   FILTER_ABS_SINC, 
254   FILTER_ABS_G_HAMMING,
255   FILTER_ABS_COSINE,
256   FILTER_SHEPP
257 } FilterType;
258
259 /* function domains */
260
261 static const char D_FREQ_STR[]=    "freq";
262 static const char D_SPATIAL_STR[]= "spatial";
263  
264 typedef enum {
265   D_FREQ,
266   D_SPATIAL 
267 } DomainType;
268
269 typedef enum {
270     FUNC_EVEN,    /* function types, f[-n] = f[n] */
271     FUNC_ODD,     /* f[-n] = -f[n] */
272     FUNC_BOTH    /* function has both odd & even components */
273 } FunctionSymmetry;
274
275 /* interpolation methods */
276 typedef enum {     /* Interpolation methods */
277   I_NEAREST,       /* Nearest neighbor */
278   I_LINEAR,        /* Linear interpolation */
279   I_BSPLINE,
280   I_1BSPLINE,      /* 1st order B-Spline */
281   I_2BSPLINE,
282   I_3BSPLINE
283 } InterpolationType;
284
285 /* Constants for sizing PHANTOM */
286
287 static const double PERCENT_PHM_SIZE_INCR=0.0;  /* Fractional increase in phantom limits compared to pelm size */
288 static const int N_EXTRA_DETECTORS=4;           /* Number of extra detectors widths when calculating detlen */
289
290 static const char O_TRACE_NONE_STR[]=     "none";
291 static const char O_TRACE_TEXT_STR[]=     "text";
292 static const char O_TRACE_PHM_STR[]=      "phm";
293 static const char O_TRACE_RAYS_STR[]=     "rays";
294 static const char O_TRACE_PLOT_STR[]=     "plot";
295 static const char O_TRACE_CLIPPING_STR[]= "clipping";
296
297 enum {
298   TRACE_NONE,           /* No tracing */
299   TRACE_TEXT,           /* Minimal status */
300   TRACE_PHM,            /* Show phantom */
301   TRACE_RAYS,           /* Show all rays */
302   TRACE_PLOT,           /* Plot raysums */
303   TRACE_CLIPPING        /* Plot clipping */
304 };
305
306 typedef enum {
307   O_PHM_HERMAN,               /* Herman head phantom */
308   O_PHM_ROWLAND,              /* Rowland head phantom */
309   O_PHM_BROWLAND,             /* Bordered Rowland head phantom */
310   O_PHM_UNITPULSE             /* Unit pulse phantom */
311 } PhantomType;
312
313 static const char O_PHM_HERMAN_STR[]=    "herman";
314 static const char O_PHM_ROWLAND_STR[]=   "rowland";
315 static const char O_PHM_BROWLAND_STR[]=  "browland";
316 static const char O_PHM_UNITPULSE_STR[]= "unitpulse";
317
318 static const char O_INTERP_NEAREST_STR[]=  "nearest";
319 static const char O_INTERP_LINEAR_STR[]=   "linear";
320 static const char O_INTERP_BSPLINE_STR[]=  "bspline";
321
322 static const char O_FILTER_ABS_BANDLIMIT_STR[]= "abs_bandlimit";
323 static const char O_FILTER_ABS_SINC_STR[]=      "abs_sinc";
324 static const char O_FILTER_ABS_COS_STR[]=       "abs_cos";
325 static const char O_FILTER_ABS_HAMMING_STR[]=   "abs_hamming";
326 static const char O_FILTER_SHEPP_STR[]=         "shepp";
327 static const char O_FILTER_BANDLIMIT_STR[]=     "bandlimit";
328 static const char O_FILTER_SINC_STR[]=          "sinc";
329 static const char O_FILTER_COS_STR[]=           "cos";
330 static const char O_FILTER_HAMMING_STR[]=       "hamming";
331 static const char O_FILTER_TRIANGLE_STR[]=      "triangle";
332
333 typedef  enum {
334   O_BPROJ_TRIG,
335   O_BPROJ_TABLE,
336   O_BPROJ_DIFF,
337   O_BPROJ_DIFF2,
338   O_BPROJ_IDIFF2
339 } BackprojType;
340
341 static const char O_BPROJ_TRIG_STR[]=     "trig";
342 static const char O_BPROJ_TABLE_STR[]=    "table";
343 static const char O_BPROJ_DIFF_STR[]=     "diff";
344 static const char O_BPROJ_DIFF2_STR[]=    "diff2";
345 static const char O_BPROJ_IDIFF2_STR[]=   "idiff2";
346
347 const static int RAYSUM_TRACE_ROW_TITLE=1;
348 const static int RAYSUM_TRACE_ROW_TITLE2=2;
349 const static int RAYSUM_TRACE_ROW_PHANT_ID=4;
350 const static int RAYSUM_TRACE_ROW_CHROMATIC=7;
351 const static int RAYSUM_TRACE_ROW_SCATTER=8;
352 const static int RAYSUM_TRACE_ROW_PHOT_STAT=9;
353 const static int RAYSUM_TRACE_ROW_NDET=12;
354 const static int RAYSUM_TRACE_ROW_NVIEW=13;
355 const static int RAYSUM_TRACE_ROW_SAMPLES=14;
356 const static int RAYSUM_TRACE_ROW_CURR_VIEW=17;
357 const static int RAYSUM_TRACE_ROW_ATTEN=18;
358
359
360 /*************************************************************************
361  *  FUNCTION DECLARATIONS
362  ************************************************************************/
363  
364 void usage (const char *program);
365
366 /* From reconstr.c */
367 IMAGE *image_reconst (IMAGE *im, RAYSUM *rs, int filt_type, double filt_param, const InterpolationType interp_type, int interp_param, const BackprojType backproj_type, int const ir_trace);
368
369 /* From bproj.c */
370 void backproj_init (const RAYSUM *rs, IMAGE *im, const BackprojType bproj_method);
371 int  backproj_calc (const RAYSUM *rs, IMAGE *im, const double *t, const double view_angle, 
372                     const int interp_type, const int bproj_method);
373 void backproj_term (const RAYSUM *rs, IMAGE *im, const int bproj_method);
374
375 void backproj_init_trig (const RAYSUM *rs, IMAGE *im);
376 int  backproj_calc_trig (const RAYSUM *rs, IMAGE *im, const double *t, 
377                          const double view_angle, const int interp_type);
378 void backproj_term_trig (const RAYSUM *rs, IMAGE *im);
379 void backproj_init_table (const RAYSUM *rs, IMAGE *im);
380 int  backproj_calc_table (const RAYSUM *rs, IMAGE *im, const double *t, 
381                           const double view_angle, const int interp_type);
382 void backproj_term_table (const RAYSUM *rs, IMAGE *im);
383 void backproj_init_d (const RAYSUM *rs, IMAGE *im);
384 int  backproj_calc_d (const RAYSUM *rs, IMAGE *im, const double *t, 
385                       const double view_angle, const int interp_type);
386 void backproj_term_d (const RAYSUM *rs, IMAGE *im);
387 void backproj_init_d2 (const RAYSUM *rs, IMAGE *im);
388 int  backproj_calc_d2 (const RAYSUM *rs, IMAGE *im, const double *t, 
389                        const double view_angle, const int interp_type);
390 void backproj_term_d2 (const RAYSUM *rs, IMAGE *im);
391 void backproj_init_id (const RAYSUM *rs, IMAGE *im);
392 int  backproj_calc_id (const RAYSUM *rs, IMAGE *im, const double *t, 
393                        const double view_angle, const int interp_type);
394 void backproj_term_id (const RAYSUM *rs, IMAGE *im);
395 void backproj_init_id2 (const RAYSUM *rs, IMAGE *im);
396 int  backproj_calc_id2 (const RAYSUM *rs, IMAGE *im, const double *t, 
397                         const double view_angle, const int interp_type);
398 void backproj_term_id2 (const RAYSUM *rs, IMAGE *im);
399
400 /* bspline.c */
401 int bspline(int samples, int zoom_factor, int spline_order, double input[], double output[]);
402
403 /* convolve.c */
404 double convolve(const double f1[], const double f2[], const double dx, const int n, const int np, const FunctionSymmetry func_type);
405 double convolve_both(const double f1[], const double f2[], const double dx, const int n, const int np);
406
407 /* dialogs.c */
408 int phm_add_pelm_kb(PHANTOM *phm);
409 PHANTOM *phm_select(void);
410 int interpolation_select(void);
411 int filter_select(double *filter_param);
412
413 /* filter.c */
414 double *filter_generate(const FilterType filt_type, double bw, double xmin, double xmax, int n, double param, const DomainType domain, int numint);
415 double filter_spatial_response_calc(int filt_type, double x, double bw, double param, int n);
416 double filter_spatial_response_analytic(int filt_type, double x, double bw, double param);
417 double filter_frequency_response(int filt_type, double u, double bw, double param);
418 double sinc(double x, double mult);
419 double integral_abscos(double u, double w);
420
421 /* image.c */
422 IMAGE *image_create(const char *fname, const int nx, const int ny);
423 int image_clear(IMAGE *im);
424 int image_save(IMAGE *im);
425 IMAGE *image_load(const char *fname);
426 void image_filter_response(IMAGE *im, const DomainType domain, double bw, const FilterType filt_type, double filt_param, const int opt_trace);
427 int image_display (const IMAGE *im);
428 int image_display_scale (const IMAGE *im, const int scale, const double pmin, const double pmax);
429
430 /* options.c */
431 int opt_set_trace(const char *optarg);
432 const char *name_of_phantom(const int phmid);
433 int opt_set_phantom(const char *optarg);
434 int opt_set_interpolation(const char *optarg);
435 const char *name_of_interpolation(int interp_type);
436 int opt_set_filter(const char *optarg);
437 const char *name_of_filter(const int filter);
438 DomainType opt_set_filter_domain(const char *optarg);
439 const char *name_of_filter_domain(const DomainType domain);
440 int opt_set_backproj(const char *optarg);
441 const char *name_of_backproj(const BackprojType backproj);
442
443 /* phm.c */
444 PHANTOM *phm_create(const int phmid);
445 PHANTOM *phm_create_from_file(const char *fname);
446 PHANTOM *phm_init(void);
447 void phm_free (PHANTOM *phm);
448 int phm_add_pelm_file(PHANTOM *phm, const char *fname);
449 void phm_add_pelm (PHANTOM *phm, const char *type, const double cx, const double cy, 
450                const double u, const double v, const double rot, const double atten);
451 int pelm_make_points(PELM *obj);
452 void pelm_make_xform (PELM *obj);
453 PELM *pelm_alloc(void);
454 void calc_arc(double x[], double y[], const int pts, const double xcent, const double ycent, 
455               const double r, const double start, const double stop);
456 void calc_ellipse(double x[], double y[], const int pts, const double u, const double v);
457 int circle_pts(double theta);
458 void phm_print(PHANTOM *phm);
459 #if HAVE_SGP
460 void phm_show(const PHANTOM *phm);
461 void phm_draw(const PHANTOM *phm);
462 #endif
463
464 /* phm2image.c */
465 void phm_to_image(const PHANTOM *phm, IMAGE *im, const int col_start, const int col_count,
466                   const int nsample, const int trace);
467 int pelm_is_point_inside(PELM *obj, const double x, const double y, const CoordType coord_type);
468
469 /* phmstd.c */
470 void phm_std_herman (PHANTOM *phm);
471 void phm_std_rowland (PHANTOM *phm);
472 void phm_std_rowland_bordered (PHANTOM *phm);
473
474 /* raycollect.c */
475 int raysum_collect(RAYSUM *rs, const DETECTOR *det, const PHANTOM *phm, const int start_view, const int trace, const int unit_pulse);
476 void rayview(const PHANTOM *phm, DETARRAY *darray, const DETECTOR *det, const double xd1, const double yd1, const double xd2, const double yd2, const double xs1, const double ys1, const double xs2, const double ys2, const int unit_pulse);
477 double phm_ray_attenuation (const PHANTOM *phm, const double x1, const double y1, const double x2, const double y2);
478 double pelm_ray_attenuation (PELM *pelm, const double x1, const double y1, const double x2, const double y2);
479 int pelm_clip_line (const PELM *pelm, double *x1, double *y1, double *x2, double *y2);
480 void raysum_trace_show_param (const char *label, const char *fmt, int row, int color, ...);
481
482 /* scanner.c */
483 DETECTOR *detector_create(const PHANTOM *phm, int geometry, int ndet, int nview, int nsample, const double rot_anglen);
484 void detector_free(DETECTOR *det);
485
486 /* rayio.c */
487 RAYSUM *raysum_create(const char *fname, const int nview, const int ndet);
488 RAYSUM *raysum_create_from_det(const char *fname, const DETECTOR *det);
489 RAYSUM *raysum_open(const char *filename);
490 void raysum_alloc_views(RAYSUM *rs);
491 void raysum_free(RAYSUM *rs);
492 int raysum_is_open(RAYSUM *rs);
493 int raysum_close(RAYSUM *rs);
494 int raysum_read_header(RAYSUM *rs);
495 int raysum_write_header(RAYSUM *rs);
496 int raysum_read(RAYSUM *rs);
497 int raysum_write(RAYSUM *rs);
498 DETARRAY *detarray_alloc(const int n);
499 void detarray_free(DETARRAY *darray);
500 int detarray_read(RAYSUM *rs, DETARRAY *darray, const int view_num);
501 int detarray_write(RAYSUM *rs, const DETARRAY *darray, const int view_num);
502 int raysum_print(const RAYSUM *rs);
503
504 #ifdef __cplusplus
505
506 #endif /* __cplusplus */
507
508 #endif