r31: Code cleanup
[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.6 2000/05/03 08:49:49 kevin Exp $
6 **  $Log: ir.h,v $
7 **  Revision 1.6  2000/05/03 08:49:49  kevin
8 **  Code cleanup
9 **
10 **  Revision 1.5  2000/05/02 20:00:25  kevin
11 **  *** empty log message ***
12 **
13 **  Revision 1.4  2000/05/02 15:31:39  kevin
14 **  code cleaning
15 **
16 **  Revision 1.3  2000/04/29 23:24:29  kevin
17 **  *** empty log message ***
18 **
19 **  Revision 1.2  2000/04/28 14:14:16  kevin
20 **  *** empty log message ***
21 **
22 **  Revision 1.1.1.1  2000/04/28 13:02:43  kevin
23 **  Initial CVS import for first public release
24 **
25 **
26 **
27 **  This program is free software; you can redistribute it and/or modify
28 **  it under the terms of the GNU General Public License (version 2) as
29 **  published by the Free Software Foundation.
30 **
31 **  This program is distributed in the hope that it will be useful,
32 **  but WITHOUT ANY WARRANTY; without even the implied warranty of
33 **  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
34 **  GNU General Public License for more details.
35 **
36 **  You should have received a copy of the GNU General Public License
37 **  along with this program; if not, write to the Free Software
38 **  Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
39 ******************************************************************************/
40 /* FILE IDENTIFICATION
41  *
42  *      Name:         ir.h          Header File for Image Reconstruction System
43  *      Programmer:   Kevin Rosenberg
44  *      Date Started: 7-1-84
45  *      Last Change:  1-20-85
46  */
47
48 #ifndef IR_H
49 #define IR_H
50
51 #ifdef MPI_CT
52 #define MPI_MAX_PROCESS 128
53 struct mpi_ct_st
54 {
55   int my_rank;
56   int nproc;
57   int base_local_work_units;
58   int remainder_work_units;
59   int local_work_units[MPI_MAX_PROCESS];
60   int start_work_unit[MPI_MAX_PROCESS];
61   MPI_Comm comm;
62 };
63
64 extern struct mpi_ct_st mpi_ct;
65 void mpi_ct_calc_work_units(const unsigned int global_work_units);
66 #endif
67
68
69 struct histo_st {
70     int *b;                     /* Histogram array (# of elements in each bin) */
71     int nbin;                   /* Number of histogram bins */
72     double xmin, xmax, xinc;    /* Limits of histogram boundaries */
73 };
74
75 typedef struct histo_st HISTOGRAM;
76
77 /*---------------------------------------------------------------------------*/
78
79 #define POINTS_PER_CIRCLE  36
80
81 #define MAXREMARK       99
82 #define LENREMARK       (MAXREMARK+1)
83 #define IMAGE_VAL       FMTX_2D         /* use 2d floating point matrix */
84 #define IMAGE_ELEM_TYPE DT_FLOAT        /* use floats for images */
85 typedef float IMAGE_ELEM_VAL;           /* use floats for image storing */
86
87 struct image_st {
88     IMAGE_VAL v;                        /* values of voxels in matrix form */
89     SDF_2D *dfp_2d;                     /* Pointer to disk image file */
90     int nx, ny;                         /* size of voxel matrix */
91     double xmin, xmax, ymin, ymax;      /* extent of voxel matrix in phm coord */
92     char remark[LENREMARK];             /* description of voxel data */
93     float calctime;                     /* time to calculate voxels in seconds */
94 };
95
96 struct object_st {
97     int type;                      /* object type (box, ellipse, etc) */
98     double atten;                  /* X-ray attenuation coefficient */
99     double cx,cy;                  /* center of object */
100     double u,v;                    /* size of object */
101     double rot;                    /* object rotation angle (in radians) */
102     double *x, *y;                 /* ptr to array of points in obj world coord */
103     int pts;                       /* number of points in outline arrays */
104     double xmin, xmax, ymin, ymax; /* object limits */
105     double radius;                 /*   "   */
106     struct {                       /* transform matrices        */
107         GRFMTX_2D p_to_o;          /* map from phm to standard obj coords */
108         GRFMTX_2D o_to_p;          /* map from std object coords to phm coords */
109     } xform;
110     struct object_st *next;       /* pointer to next object in PHANTOM */
111 };
112
113 typedef struct object_st OBJECT;
114
115 struct phm_st {                    /* PHANTOM structure */
116     OBJECT *objlist;               /* object list */
117     int type;
118     int numobj;                    /* number of objects in PHANTOM */
119     double xmin, xmax, ymin, ymax; /* extent of objects in object coordinates */
120     double radius;                 /*     " "    */
121 };
122
123 typedef struct image_st  IMAGE;
124 typedef struct phm_st    PHANTOM;
125
126 #define P_OBJECTS       0       /* PHANTOM made of objects */
127 #define P_UNIT_PULSE    1       /* Special PHANTOM, not made of objects */
128 #define P_FILTER        9       /* defined only by a type */
129
130 /*----------------------------------------------------------------------*/
131 /*                              RAYSUM SYMBOLS                          */
132 /*----------------------------------------------------------------------*/
133
134 /* Ray sums are collected along an array of ndet detectors.  The data
135  * for these detectors is stored in the structure DETECTARRAY
136  */
137
138 #define DETECT_TYPE float
139
140 struct detarray_st {
141   DETECT_TYPE *detval;  /* Pointer to array of values recorded by detector */
142   int ndet;             /* Number of detectors in array */
143   double view_angle;    /* View angle in radians */
144 };
145
146 #define DETECTOR_PARALLEL   1
147 #define DETECTOR_EQUIANGLE  2
148 #define DETECTOR_EQUILINEAR 3
149
150 struct detect_st {
151   int geometry;                 /* Geometry of detectory */
152   int ndet;                     /* Number of detectors in array */
153   int nview;                    /* Number of rotated views */
154   int nsample;                  /* Number of rays per detector */
155   double detlen;                /* Total length of detector array */
156   double rotlen;                /* Rotation angle length in radians (norm 2PI) */
157   double det_inc;               /* Increment between centers of detectors */
158   double rot_inc;               /* Increment in rotation angle between views */
159   double radius;                /* Radius of rotation.  Distance from */
160                                 /*   center of phm to center of det */
161   double phmlen;                /* Maximum Length of PHANTOM or area of interest */
162   struct {
163     double xd1,yd1,xd2,yd2;     /* Coordinates of detector endpoints */
164     double xs1,ys1,xs2,ys2;     /* Coordinates of source endpoints */
165     double angle;               /* Starting angle */
166   } init;
167 };
168
169 struct raysum_st {
170   int fd;
171   int file_mode;
172   int header_size;
173   struct detarray_st **view;    /* Pointer to array of detarray_st pointers */
174
175   char remark[LENREMARK];       /* description of raysum data */
176   double calctime;              /* time required to calculate raysums */
177
178   int ndet;                     /* number of detectors in array */
179   int nview;                    /* number of rotated views */
180   double rot_start;             /* starting view rotation */
181   double rot_inc;               /* angle between rotations */
182   double det_start;             /* distance of beginning detector to center */
183                                 /*    of PHANTOM */
184   double det_inc;               /* increment between detectors */
185   double phmlen;                /* Length of PHANTOM edge (phm is square) */
186 };
187
188 typedef struct detarray_st DETARRAY;
189 typedef struct detect_st   DETECTOR;
190 typedef struct raysum_st   RAYSUM;
191
192 /*----------------------------------------------------------------------*/
193 /*                              USER SYMBOLS                            */
194 /*----------------------------------------------------------------------*/
195
196 /* codes for object types, passed to add_obj() */
197
198 #define O_RECTANGLE  1
199 #define O_TRIANGLE   2
200 #define O_ELLIPSE    3
201 #define O_SECTOR     4
202 #define O_SEGMENT    5
203
204 /* Codes for Coordinate Types      */
205 /* Defines coords for inside_obj() */
206
207 #define OBJ_COORD -1            /* Normalized Object Coordinates */
208 #define PHM_COORD -2            /* User's PHANTOM Coordinates */
209
210 /* Codes for Filter types */
211
212 #define W_BANDLIMIT     1       /* filter types for genfilter() */
213 #define W_SINC          2
214 #define W_G_HAMMING     3
215 #define W_COSINE        4
216 #define W_TRIANGLE      5
217
218 #define W_A_BANDLIMIT   11      /* filters times abs() of function */
219 #define W_A_SINC        12
220 #define W_AG_HAMMING    13
221 #define W_A_COSINE      14
222
223 #define W_SHEPP         21
224
225 /* function domains */
226
227 #define O_FREQ_STR      "freq"
228 #define O_SPATIAL_STR   "spatial"
229  
230 #define D_FREQ          1       /* Domain names */
231 #define D_SPATIAL       2
232
233 /* function symmetry */
234
235 #define FUNC_EVEN       1       /* function types, f[-n] = f[n] */
236 #define FUNC_ODD        2       /* f[-n] = -f[n] */
237 #define FUNC_BOTH       3       /* function has both odd & even components */
238
239 /* interpolation methods */
240
241 #define I_NEAREST       1       /* Interpolation methods */
242 #define I_LINEAR        2       /* Linear interpolation */
243 #define I_BSPLINE       3
244 #define I_1BSPLINE      3       /* 1st order B-Spline */
245 #define I_2BSPLINE      4
246 #define I_3BSPLINE      5
247
248 /* Constants for sizing PHANTOM */
249
250 #define PERCENT_PHM_SIZE_INCR   0.0     /* Fractional increase in PHANTOM */
251                                         /* limits compared to object size */
252 #define N_EXTRA_DETECTORS         4     /* Number of extra detectors */
253                                         /* widths when calculating detlen */
254
255 #define DET_PARALLEL 1
256 #define DET_FAN      2
257
258 #define O_TRACE_NONE_STR     "none"
259 #define O_TRACE_TEXT_STR     "text"
260 #define O_TRACE_PHM_STR      "phm"
261 #define O_TRACE_RAYS_STR     "rays"
262 #define O_TRACE_PLOT_STR     "plot"
263 #define O_TRACE_CLIPPING_STR "clipping"
264
265
266 #define TRACE_NONE     0                /* No tracing */
267 #define TRACE_TEXT     1                /* Minimal status */
268 #define TRACE_PHM      2                /* Show PHANTOM */
269 #define TRACE_RAYS     3                /* Show all rays */
270 #define TRACE_PLOT     4                /* Plot raysums */
271 #define TRACE_CLIPPING 5                /* Plot clipping */
272
273 #define O_PHM_HERMAN    1               /* Herman head phantom */
274 #define O_PHM_ROWLAND   2               /* Rowland head phantom */
275 #define O_PHM_BROWLAND  3               /* Bordered Rowland head phantom */
276 #define O_PHM_UNITPULSE 4               /* Unit pulse phantom */
277
278 #define O_PHM_HERMAN_STR    "herman"
279 #define O_PHM_ROWLAND_STR   "rowland"
280 #define O_PHM_BROWLAND_STR  "browland"
281 #define O_PHM_UNITPULSE_STR "unitpulse"
282
283 #define O_INTERP_NEAREST_STR  "nearest"
284 #define O_INTERP_LINEAR_STR   "linear"
285 #define O_INTERP_BSPLINE_STR  "bspline"
286
287 #define O_FILTER_ABS_BANDLIMIT_STR "abs_bandlimit"
288 #define O_FILTER_ABS_SINC_STR      "abs_sinc"
289 #define O_FILTER_ABS_COS_STR       "abs_cos"
290 #define O_FILTER_ABS_HAMMING_STR   "abs_hamming"
291 #define O_FILTER_SHEPP_STR         "shepp"
292 #define O_FILTER_BANDLIMIT_STR     "bandlimit"
293 #define O_FILTER_SINC_STR          "sinc"
294 #define O_FILTER_COS_STR           "cos"
295 #define O_FILTER_HAMMING_STR       "hamming"
296 #define O_FILTER_TRIANGLE_STR      "triangle"
297
298 #define O_BPROJ_TRIG   1 
299 #define O_BPROJ_TABLE  2
300 #define O_BPROJ_DIFF   3
301 #define O_BPROJ_DIFF2  4
302 #define O_BPROJ_IDIFF2 5
303
304 #define O_BPROJ_TRIG_STR     "trig"
305 #define O_BPROJ_TABLE_STR    "table"
306 #define O_BPROJ_DIFF_STR     "diff"
307 #define O_BPROJ_DIFF2_STR    "diff2"
308 #define O_BPROJ_IDIFF2_STR   "idiff2"
309
310 #define RS_TRACE_ROW_TITLE      1
311 #define RS_TRACE_ROW_TITLE2     2
312 #define RS_TRACE_ROW_PHANT_ID   4
313 #define RS_TRACE_ROW_CHROMATIC  7
314 #define RS_TRACE_ROW_SCATTER    8
315 #define RS_TRACE_ROW_PHOT_STAT  9
316 #define RS_TRACE_ROW_NDET       12
317 #define RS_TRACE_ROW_NVIEW      13
318 #define RS_TRACE_ROW_SAMPLES    14
319 #define RS_TRACE_ROW_CURR_VIEW  17
320 #define RS_TRACE_ROW_ATTEN      18
321
322
323 /*----------------------------------------------------------------------*/
324 /*                      GRAY SCALE STRUCTURES                           */
325 /*----------------------------------------------------------------------*/
326
327 #define GS_MAX_CELL_SIZE   4
328
329 typedef int GS_BITMASK[4][4];
330
331 struct greyscale_st {
332     int dev;                            /* Device to output to */
333     int (*dotfunc)(int x, int y, int color);                    /* Pointer to dot function for device */
334     int cur_x, cur_y;                   /* Current cell location */
335     int nxcell, nycell;                 /* size of cell in pixels */
336     int xmin, ymin;                     /* starting position of grey scale */
337     int num_color;                      /* Number of primary colors available */
338     int num_intens;                     /* Number of intensities available */
339     int max_level;                      /* gs levels range from 0 to max_level */
340     char *fg_color_tbl;                 /* Hold foreground color for each level */
341     char *bg_color_tbl;                 /* Holds background color */
342     char *level_sub_tbl;                /* Holds value to subtract for level */
343                                         /* before accessing bit mask */
344     GS_BITMASK *bm;                     /* Holds grey-scale bit mask */
345     struct greyscale_st *next_dev;      /* Pointer to next open device */
346                                         /* == NULL when no more devices */
347 };
348
349 typedef struct greyscale_st GREYSCALE;
350
351
352 /* From reconstr.c */
353 IMAGE *image_reconst (IMAGE *im, RAYSUM *rs, int filt_type, double filt_param, int interp_type, int interp_param, const int backproj_type, int ir_trace);
354
355 /* From bproj.c */
356 void backproj_init (const RAYSUM *rs, IMAGE *im, const int bproj_method);
357 int  backproj_calc (const RAYSUM *rs, IMAGE *im, const double *t, const double view_angle, 
358                     const int interp_type, const int bproj_method);
359 void backproj_term (const RAYSUM *rs, IMAGE *im, const int bproj_method);
360
361 void backproj_init_trig (const RAYSUM *rs, IMAGE *im);
362 int  backproj_calc_trig (const RAYSUM *rs, IMAGE *im, const double *t, 
363                          const double view_angle, const int interp_type);
364 void backproj_term_trig (const RAYSUM *rs, IMAGE *im);
365 void backproj_init_table (const RAYSUM *rs, IMAGE *im);
366 int  backproj_calc_table (const RAYSUM *rs, IMAGE *im, const double *t, 
367                           const double view_angle, const int interp_type);
368 void backproj_term_table (const RAYSUM *rs, IMAGE *im);
369 void backproj_init_d (const RAYSUM *rs, IMAGE *im);
370 int  backproj_calc_d (const RAYSUM *rs, IMAGE *im, const double *t, 
371                       const double view_angle, const int interp_type);
372 void backproj_term_d (const RAYSUM *rs, IMAGE *im);
373 void backproj_init_d2 (const RAYSUM *rs, IMAGE *im);
374 int  backproj_calc_d2 (const RAYSUM *rs, IMAGE *im, const double *t, 
375                        const double view_angle, const int interp_type);
376 void backproj_term_d2 (const RAYSUM *rs, IMAGE *im);
377 void backproj_init_id (const RAYSUM *rs, IMAGE *im);
378 int  backproj_calc_id (const RAYSUM *rs, IMAGE *im, const double *t, 
379                        const double view_angle, const int interp_type);
380 void backproj_term_id (const RAYSUM *rs, IMAGE *im);
381 void backproj_init_id2 (const RAYSUM *rs, IMAGE *im);
382 int  backproj_calc_id2 (const RAYSUM *rs, IMAGE *im, const double *t, 
383                         const double view_angle, const int interp_type);
384 void backproj_term_id2 (const RAYSUM *rs, IMAGE *im);
385
386 void usage (const char *program);
387 int main(const int argc, char * const argv[]);
388
389
390
391 /* bspline.c */
392 int bspline(int samples, int zoom_factor, int spline_order, double input[], double output[]);
393 /* convolve.c */
394 double convolve(const double f1[], const double f2[], const double dx, const int n, const int np, const int func_type);
395 double convolve_both(const double f1[], const double f2[], const double dx, const int n, const int np);
396 void rotate2d(double x[], double y[], int pts, double angle);
397 void xlat2d(double x[], double y[], int pts, double xoffset, double yoffset);
398 void scale2d(double x[], double y[], int pts, double xfact, double yfact);
399 int circle_pts(double theta);
400 /* dialogs.c */
401 int add_objs_kb(PHANTOM *phm);
402 PHANTOM *phm_select(void);
403 int interpolation_select(void);
404 int filter_select(double *filter_param);
405 /* filter.c */
406 double *filter_generate(int filt_type, double bw, double xmin, double xmax, int n, double param, int domain, int numint);
407 double filter_spatial_response_calc(int filt_type, double x, double bw, double param, int n);
408 double filter_spatial_response_analytic(int filt_type, double x, double bw, double param);
409 double filter_frequency_response(int filt_type, double u, double bw, double param);
410 double sinc(double x, double mult);
411 double abscos_int(double u, double w);
412 /* image.c */
413 IMAGE *image_create(const char *fname, const int nx, const int ny);
414 int image_clear(IMAGE *im);
415 int image_save(IMAGE *im);
416 IMAGE *image_load(const char *fname);
417 void image_filter_response(IMAGE *im, int domain, double bw, int filt_type, double filt_param, int opt_trace);
418 /* options.c */
419 int opt_set_trace(const char *optarg, const char *program);
420 const char *name_of_phantom(const int phmid);
421 int opt_set_phantom(const char *optarg, const char *program);
422 int opt_set_interpolation(const char *optarg, const char *program);
423 const char *name_of_interpolation(int interp_type);
424 int opt_set_filter(const char *optarg, const char *program);
425 const char *name_of_filter(const int filter);
426 int opt_set_filter_domain(const char *optarg, const char *program);
427 const char *name_of_filter_domain(const int domain);
428 int opt_set_backproj(const char *optarg, const char *program);
429 const char *name_of_backproj(const int backproj);
430 /* phm.c */
431 PHANTOM *phm_create(const int phmid);
432 PHANTOM *phm_create_from_file(const char *fname);
433 PHANTOM *phm_init(void);
434 int add_objs_kb(PHANTOM *phm);
435 int add_objs_file(PHANTOM *phm, const char *fname);
436 void addobject(PHANTOM *phm, const int type, const double cx, const double cy, 
437                const double u, const double v, const double rot, const double atten);
438 int makeobjpts(OBJECT *obj);
439 void makeobjxform(OBJECT *obj);
440 void calc_arc(double x[], double y[], const int pts, const double xcent, const double ycent, 
441               const double r, const double start, const double stop);
442 void calc_ellipse(double x[], double y[], const int pts, const double u, const double v);
443 OBJECT *alloc_obj(void);
444 void phm_print(PHANTOM *phm);
445 #if HAVE_INTERACTIVE_GRAPHICS
446 void phm_show(const PHANTOM *phm);
447 void phm_draw(const PHANTOM *phm);
448 #endif
449 /* phm2image.c */
450 void phm_to_image(const PHANTOM *phm, IMAGE *im, const int col_start, const int col_count,
451                   const int nsample, const int trace);
452 int inside_obj(const OBJECT *obj, double x, double y, const int coord_type);
453 /* phmstd.c */
454 void herm_head(PHANTOM *phm);
455 void row_head(PHANTOM *phm);
456 void row_bord_head(PHANTOM *phm);
457 /* ray.c */
458 DETECTOR *detect_create(const PHANTOM *phm, int ndet, int nview, int nsample, const double rot_anglen);
459 void detect_free(DETECTOR *det);
460 double calc_rsum(const PHANTOM *phm, const double x1, const double y1, const double x2, const double y2);
461 double calc_objsum(const OBJECT *obj, const double x1, const double y1, const double x2, const double y2);
462 int clipobj(const OBJECT *obj, double *x1, double *y1, double *x2, double *y2);
463 /* raycollect.c */
464 int raysum_collect(RAYSUM *rs, const DETECTOR *det, const PHANTOM *phm, const int start_view, 
465                    const int trace, const int unit_pulse);
466 void rayview(const PHANTOM *phm, DETARRAY *darray, const DETECTOR *det, 
467              const double xd1, const double yd1, const double xd2, const double yd2, 
468              const double xs1, const double ys1, const double xs2, const double ys2, const int unit_pulse);
469 void rs_trace_showprm (const char *label, const char *fmt, int row, int color, ...);
470 /* rayio.c */
471 RAYSUM *raysum_create(const char *fname, const int nview, const int ndet);
472 RAYSUM *raysum_create_from_det(const char *fname, const DETECTOR *det);
473 RAYSUM *raysum_open(const char *filename);
474 void raysum_alloc_views(RAYSUM *rs);
475 void raysum_free(RAYSUM *rs);
476 int raysum_is_open(RAYSUM *rs);
477 int raysum_close(RAYSUM *rs);
478 int raysum_read_header(RAYSUM *rs);
479 int raysum_write_header(RAYSUM *rs);
480 int raysum_read(RAYSUM *rs);
481 int raysum_write(RAYSUM *rs);
482 DETARRAY *detarray_alloc(const int n);
483 void detarray_free(DETARRAY *darray);
484 int detarray_read(RAYSUM *rs, DETARRAY *darray, const int view_num);
485 int detarray_write(RAYSUM *rs, const DETARRAY *darray, const int view_num);
486 int raysum_print(const RAYSUM *rs);
487
488
489 #endif