r3: Initial revision
[ctsim.git] / src / ctrec.c
1 /*****************************************************************************
2 **  This is part of the CTSim program
3 **  Copyright (C) 1983-2000 Kevin Rosenberg
4 **
5 **  $Id: ctrec.c,v 1.1 2000/04/28 13:02:44 kevin Exp $
6 **  $Log: ctrec.c,v $
7 **  Revision 1.1  2000/04/28 13:02:44  kevin
8 **  Initial revision
9 **
10 **
11 **  This program is free software; you can redistribute it and/or modify
12 **  it under the terms of the GNU General Public License (version 2) as
13 **  published by the Free Software Foundation.
14 **
15 **  This program is distributed in the hope that it will be useful,
16 **  but WITHOUT ANY WARRANTY; without even the implied warranty of
17 **  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
18 **  GNU General Public License for more details.
19 **
20 **  You should have received a copy of the GNU General Public License
21 **  along with this program; if not, write to the Free Software
22 **  Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
23 ******************************************************************************/
24 /* FILE
25  *   ctrec.c            Reconstruct an image from raysums
26  *
27  * DATE
28  *   Aug 84
29  *   Jul 99  -- Converted to ANSI C
30  *              Added MPI parallel processing
31  */
32
33 #include "ct.h"
34
35 #define O_INTERP       1
36 #define O_FILTER       2
37 #define O_FILTER_PARAM 3
38 #define O_BACKPROJ     4
39 #define O_VERBOSE      5
40 #define O_TRACE        6
41 #define O_HELP         7
42 #define O_DEBUG        8
43 #define O_VERSION      9
44
45 static struct option my_options[] =
46 {
47   {"interp", 1, 0, O_INTERP},
48   {"filter", 1, 0, O_FILTER},
49   {"filter-param", 1, 0, O_FILTER_PARAM},
50   {"backproj", 1, 0, O_BACKPROJ},
51   {"trace", 1, 0, O_TRACE},
52   {"debug", 0, 0, O_DEBUG},
53   {"verbose", 0, 0, O_VERBOSE},
54   {"help", 0, 0, O_HELP},
55   {"version", 0, 0, O_VERSION},
56   {0, 0, 0, 0}
57 };
58
59 void 
60 usage (const char *program)
61 {
62   fprintf(stdout,"usage: %s raysum-file image-file nx-image ny-image [OPTIONS]\n", kbasename(program));
63   fprintf(stdout,"Image reconstruction from raysum projections\n");
64   fprintf(stdout,"\n");
65   fprintf(stdout,"   raysum-file     Input raysum file\n");
66   fprintf(stdout,"   image-file      Output image file in SDF2D format\n");
67   fprintf(stdout,"   nx-image        Number of columns in output image\n");
68   fprintf(stdout,"   ny-image        Number of rows in output image\n");
69   fprintf(stdout,"   --interp        Interpolation method during backprojection\n");
70   fprintf(stdout,"       nearest     Nearest neighbor interpolation\n");
71   fprintf(stdout,"       linear      Linear interpolation\n");
72   fprintf(stdout,"       bspline     B-spline interpolation\n");
73   fprintf(stdout,"    --filter       Filter name\n");
74   fprintf(stdout,"       abs_bandlimit Abs * Bandlimiting (default)\n");
75   fprintf(stdout,"       abs_sinc      Abs * Sinc\n");
76   fprintf(stdout,"       abs_cos       Abs * Cosine\n");
77   fprintf(stdout,"       abs_hamming   Abs * Hamming\n");
78   fprintf(stdout,"       shepp         Shepp-Logan\n");
79   fprintf(stdout,"       bandlimit     Bandlimiting\n");
80   fprintf(stdout,"       sinc          Sinc\n");
81   fprintf(stdout,"       cos           Cosine\n");
82   fprintf(stdout,"       triangle      Triangle\n");
83   fprintf(stdout,"       hamming       Hamming\n");
84   fprintf(stdout,"    --backproj     Backprojection Method\n");
85   fprintf(stdout,"       trig        Trigometric functions at every point\n");
86   fprintf(stdout,"       table       Trigometric functions with precalculated table\n");
87   fprintf(stdout,"       diff        Difference method\n");
88   fprintf(stdout,"       diff2       Optimized difference method (default)\n");
89   fprintf(stdout,"       idiff2      Optimized difference method with integer math\n");
90   fprintf(stdout,"    --filter-param Alpha level for Hamming filter\n");
91   fprintf(stdout,"    --trace        Set tracing to level\n");
92   fprintf(stdout,"         none      No tracing (default)\n");
93   fprintf(stdout,"         status    Status tracing\n");
94   fprintf(stdout,"         pic       Trace picture\n");
95   fprintf(stdout,"         rays      Trace allrays\n");
96   fprintf(stdout,"         plot      Trace plotting\n");
97   fprintf(stdout,"         clipping  Trace clipping\n");
98   fprintf(stdout,"    --verbose      Turn on verbose mode\n");
99   fprintf(stdout,"    --debug        Turn on debug mode\n");
100   fprintf(stdout,"    --version      Print version\n");
101   fprintf(stdout,"    --help         Print this help message\n");
102   exit(1);  
103 }
104
105
106 #ifdef MPI_CT
107 void mpi_scatter_rs (RAYSUM *rs_global, RAYSUM *rs_local, const int debug);
108 #endif
109 void print_raysum_info(const RAYSUM *rs);
110
111 int 
112 main (const int argc, char *const argv[])
113 {
114   IMAGE *im_global;
115   RAYSUM *rs_global;
116   char *rs_name, *im_filename;
117   char remark[MAXREMARK];
118   char filt_name[80];
119   int nx, ny;
120   double time_start, time_end;
121   char *endptr;
122   int opt_verbose = 0;
123   int opt_debug = 0;
124   int opt_trace = TRACE_NONE;
125   double opt_filter_param;
126   int opt_filter = W_A_BANDLIMIT;
127   int opt_interp = I_LINEAR;
128   int opt_interp_param = 1;
129   int opt_backproj = O_BPROJ_DIFF2;
130 #ifdef MPI_CT
131   IMAGE *im_local;
132   RAYSUM *rs_local;
133   int mpi_argc = argc;
134   char **mpi_argv = (char **) argv;
135   int mpi_nview, mpi_ndet;
136   double mpi_detinc, mpi_rotinc, mpi_piclen;
137   double mpi_t1, mpi_t2, mpi_t, mpi_t_g;
138   int ix;
139
140   MPI_Init(&mpi_argc, &mpi_argv);
141   MPI_Comm_dup (MPI_COMM_WORLD, &mpi_ct.comm);
142   MPI_Comm_size(mpi_ct.comm, &mpi_ct.nproc);
143   MPI_Comm_rank(mpi_ct.comm, &mpi_ct.my_rank);
144
145   if (mpi_ct.nproc > MPI_MAX_PROCESS) {
146     sys_error(ERR_FATAL, "Number of mpi processes (%d) exceeds max processes (%d)",
147               mpi_ct.nproc, MPI_MAX_PROCESS);
148     exit(1);
149   }
150 #endif
151
152 #ifdef MPI_CT
153   time_start = MPI_Wtime();
154 #else
155   time_start = td_current_sec();
156 #endif
157
158 #ifdef MPI_CT
159   if (mpi_ct.my_rank == 0) {
160 #endif
161     while (1) {
162       int c = getopt_long(argc, argv, "", my_options, NULL);
163       char *endptr = NULL;
164       
165       if (c == -1)
166         break;
167       
168       switch (c)
169         {
170         case O_INTERP:
171           opt_interp = opt_set_interpolation(optarg, argv[0]);
172           break;
173         case O_FILTER:
174           opt_filter = opt_set_filter(optarg, argv[0]);
175           break;
176         case O_BACKPROJ:
177           opt_backproj = opt_set_backproj(optarg, argv[0]);
178           break;
179         case O_FILTER_PARAM:
180           opt_filter_param = strtod(optarg, &endptr);
181           if (endptr != optarg + strlen(optarg)) {
182             usage(argv[0]);
183             exit(1);
184           }
185           break;
186         case O_VERBOSE:
187           opt_verbose = 1;
188           break;
189         case O_DEBUG:
190           opt_debug = 1;
191           break;
192         case O_TRACE:
193           opt_trace = opt_set_trace(optarg, argv[0]);
194           break;
195         case O_VERSION:
196 #ifdef VERSION
197           fprintf(stdout, "Version %s\n", VERSION);
198 #else
199           fprintf(stderr, "Unknown version number");
200 #endif
201           exit(0);
202         case O_HELP:
203         case '?':
204           usage(argv[0]);
205           exit(0);
206         default:
207           usage(argv[0]);
208           exit(1);
209         }
210     }
211   
212     if (optind + 4 != argc) {
213       usage(argv[0]);
214       exit(1);
215     }
216
217     rs_name = argv[optind];
218   
219     im_filename = argv[optind + 1];
220   
221     nx = strtol(argv[optind + 2], &endptr, 10);
222     ny = strtol(argv[optind + 3], &endptr, 10);
223   
224     if (opt_filter == W_G_HAMMING || opt_filter == W_AG_HAMMING)
225       sprintf (filt_name, "%s: alpha = %.2f",
226                filter_name_of (opt_filter), opt_filter_param); 
227     else
228       sprintf (filt_name, "%s", filter_name_of (opt_filter));
229   
230     sprintf (remark,
231              "Reconstruct: %dx%d, %s, %s, %s",
232              nx, ny, filt_name, interp_name_of (opt_interp), name_of_backproj(opt_backproj));
233   
234     if (opt_verbose)
235       fprintf (stdout, "%s\n", remark);
236 #ifdef MPI_CT
237   }
238 #endif
239
240 #ifdef MPI_CT
241   if (mpi_ct.my_rank == 0) {
242     if (file_exists(rs_name) == FALSE) {
243       fprintf (stderr, "File %s does not exist\n", rs_name);
244       exit(1);
245     } 
246     rs_global = raysum_open (rs_name);
247     raysum_read (rs_global);
248     print_raysum_info(rs_global);
249     if (opt_verbose) {
250       printf ("Number of detectors: %d\n", rs_global->ndet);
251       printf ("    Number of views: %d\n", rs_global->nview);
252       printf ("             Remark: %s\n", rs_global->remark);
253     }
254   
255     mpi_ndet = rs_global->ndet;
256     mpi_nview = rs_global->nview;
257     mpi_detinc = rs_global->det_inc;
258     mpi_piclen = rs_global->piclen;
259     mpi_rotinc = rs_global->rot_inc;
260   }
261
262
263   mpi_t1 = MPI_Wtime();
264   MPI_Bcast(&opt_verbose, 1, MPI_INT, 0, mpi_ct.comm);
265   MPI_Bcast(&opt_debug, 1, MPI_INT, 0, mpi_ct.comm);
266   MPI_Bcast(&opt_trace, 1, MPI_INT, 0, mpi_ct.comm);
267   MPI_Bcast(&opt_filter, 1, MPI_INT, 0, mpi_ct.comm);
268   MPI_Bcast(&opt_interp, 1, MPI_INT, 0, mpi_ct.comm);
269   MPI_Bcast(&opt_filter_param, 1, MPI_DOUBLE, 0, mpi_ct.comm);
270   MPI_Bcast(&opt_interp_param, 1, MPI_INT, 0, mpi_ct.comm);
271   MPI_Bcast(&opt_backproj, 1, MPI_INT, 0, mpi_ct.comm);
272   MPI_Bcast(&mpi_ndet, 1, MPI_INT, 0, mpi_ct.comm);
273   MPI_Bcast(&mpi_nview, 1, MPI_INT, 0, mpi_ct.comm);
274   MPI_Bcast(&mpi_detinc, 1, MPI_DOUBLE, 0, mpi_ct.comm);
275   MPI_Bcast(&mpi_piclen, 1, MPI_DOUBLE, 0, mpi_ct.comm);
276   MPI_Bcast(&mpi_rotinc, 1, MPI_DOUBLE, 0, mpi_ct.comm);
277   MPI_Bcast(&nx, 1, MPI_INT, 0, mpi_ct.comm);
278   MPI_Bcast(&ny, 1, MPI_INT, 0, mpi_ct.comm);
279   if (opt_verbose) {
280     mpi_t2 = MPI_Wtime();
281     mpi_t = mpi_t2 - mpi_t1;
282     MPI_Reduce(&mpi_t, &mpi_t_g, 1, MPI_DOUBLE, MPI_MAX, 0, mpi_ct.comm);
283     if (mpi_ct.my_rank == 0)
284       printf("Time to Bcast vars = %f secs, Max time = %f\n", mpi_t, mpi_t_g);
285   }
286
287   mpi_ct_calc_work_units(mpi_nview);
288
289   rs_local = raysum_create (NULL, mpi_ct.local_work_units[mpi_ct.my_rank], mpi_ndet, TRUE);
290   //  rs_local->ndet = mpi_ndet;
291   //  rs_local->nview = mpi_nview;
292   rs_local->det_inc = mpi_detinc;
293   rs_local->piclen = mpi_piclen;
294   rs_local->rot_inc = mpi_rotinc;
295
296   if (opt_debug)
297     fprintf(stderr, "Bcast'd rs vars, my nview=%d, local_work_units=%d, start_work_units=%d (process %d)\n",
298             mpi_nview, mpi_ct.local_work_units[mpi_ct.my_rank], mpi_ct.start_work_unit[mpi_ct.my_rank], mpi_ct.my_rank);
299
300   if (opt_verbose)
301     mpi_t1 = MPI_Wtime();
302   mpi_scatter_rs(rs_global, rs_local, opt_debug);
303   if (opt_verbose) {
304     mpi_t2 = MPI_Wtime();
305     mpi_t = mpi_t2 - mpi_t1;
306     MPI_Reduce(&mpi_t, &mpi_t_g, 1, MPI_DOUBLE, MPI_MAX, 0, mpi_ct.comm);
307     if (mpi_ct.my_rank == 0)
308       printf("Time to scatter rs = %f secs, Max time = %f sec\n", mpi_t, mpi_t_g);
309   }
310
311   if (mpi_ct.my_rank == 0) {
312     im_global = image_create (im_filename, nx, ny);
313     sdf_add_label (LT_HISTORY, rs_global->remark, rs_global->calctime, im_global->dfp_2d->dfp);
314     sdf_add_empty_label (im_global->dfp_2d->dfp);
315   }
316   im_local = image_create (NULL, nx, ny);
317
318 #else
319   rs_global = raysum_open (rs_name);
320   raysum_read (rs_global);
321   if (opt_debug)
322     print_raysum_info(rs_global);
323
324   im_global = image_create (im_filename, nx, ny);
325   sdf_add_label (LT_HISTORY, rs_global->remark, rs_global->calctime, im_global->dfp_2d->dfp);
326   sdf_add_empty_label (im_global->dfp_2d->dfp);
327 #endif
328
329 #ifdef MPI_CT
330   mpi_t1 = MPI_Wtime();
331   image_reconst (im_local, rs_local, opt_filter, opt_filter_param, 
332                  opt_interp, opt_interp_param, opt_backproj, opt_trace);
333   mpi_t2 = MPI_Wtime();
334   mpi_t = mpi_t2 - mpi_t1;
335   MPI_Reduce(&mpi_t, &mpi_t_g, 1, MPI_DOUBLE, MPI_MAX, 0, mpi_ct.comm);
336   if (mpi_ct.my_rank == 0) {
337     printf("Time to reconstruct = %f, Max time = %f\n", mpi_t, mpi_t_g);
338   }
339 #else
340   image_reconst (im_global, rs_global, opt_filter, opt_filter_param, 
341                  opt_interp, opt_interp_param, opt_backproj, opt_trace);
342 #endif
343
344 #ifdef MPI_CT
345   if (mpi_ct.my_rank == 0) {
346     raysum_close (rs_global);
347   }
348
349   if (opt_verbose)
350     mpi_t1 = MPI_Wtime();
351   for (ix = 0; ix < im_local->nx; ix++) {
352     MPI_Reduce(im_local->v[ix], im_global->v[ix], 
353                im_local->ny, MPI_FLOAT, MPI_SUM, 0, mpi_ct.comm);
354   }
355   if (opt_verbose) {
356     mpi_t2 = MPI_Wtime();
357     mpi_t = mpi_t2 - mpi_t1;
358     MPI_Reduce(&mpi_t, &mpi_t_g, 1, MPI_DOUBLE, MPI_MAX, 0, mpi_ct.comm);
359     if (mpi_ct.my_rank == 0)
360       printf("Time to reduce image = %f secs, max time = %f\n", mpi_t, mpi_t_g);
361   }
362   if (mpi_ct.my_rank == 0) {
363     strncpy (im_global->remark, remark, MAXREMARK);
364     image_save (im_global);
365     time_end = MPI_Wtime();
366     im_global->calctime = time_end - time_start;
367     if (opt_verbose)
368       fprintf (stdout, "Time active = %.2f\n", im_global->calctime);
369     }
370 #else  
371   raysum_close (rs_global);
372   strncpy (im_global->remark, remark, MAXREMARK);
373   image_save (im_global);
374
375   time_end = td_current_sec();
376   im_global->calctime = time_end - time_start;
377   if (opt_verbose)
378     fprintf (stdout, "Time active = %.2f secs\n", im_global->calctime);
379 #endif
380
381 #ifdef MPI_CT
382         MPI_Finalize();
383 #endif
384
385         return (0);
386 }
387
388
389 #ifdef MPI_CT
390 void mpi_scatter_rs (RAYSUM *rs_global, RAYSUM *rs_local, const int opt_debug)
391 {
392   int iproc;
393   int end_work_unit;
394   int iw;
395
396   if (mpi_ct.my_rank == 0) {
397     for (iproc = 0; iproc < mpi_ct.nproc; iproc++) {
398       end_work_unit = mpi_ct.start_work_unit[iproc] 
399         + mpi_ct.local_work_units[iproc] 
400         - 1;
401
402       if (opt_debug)
403         fprintf(stdout, "Sending rs data to process %d\n", iproc);
404
405       for (iw = mpi_ct.start_work_unit[iproc]; iw <= end_work_unit; iw++) {
406         MPI_Send(&rs_global->view[iw]->view_angle, 1, MPI_DOUBLE, iproc, 0, mpi_ct.comm);
407         MPI_Send(&rs_global->view[iw]->ndet, 1, MPI_INT, iproc, 0, mpi_ct.comm);
408         MPI_Send(rs_global->view[iw]->detval, rs_global->ndet, MPI_FLOAT, iproc, 0, mpi_ct.comm);
409       }
410     }
411   }
412
413   end_work_unit = mpi_ct.local_work_units[mpi_ct.my_rank] - 1;
414   for (iw = 0; iw <= end_work_unit; iw++) {
415     MPI_Status status;
416
417     MPI_Recv(&rs_local->view[iw]->view_angle, 1, MPI_DOUBLE, 0, 0, mpi_ct.comm, &status);
418     MPI_Recv(&rs_local->view[iw]->ndet, 1, MPI_INT, 0, 0, mpi_ct.comm, &status);
419     MPI_Recv(rs_local->view[iw]->detval, rs_global->ndet, MPI_FLOAT, 0, 0, mpi_ct.comm, &status);
420   }
421   rs_local->nview = mpi_ct.local_work_units[mpi_ct.my_rank];
422 }
423
424 #endif
425
426 void print_raysum_info(const RAYSUM *rs)
427 {
428   printf ("Number of detectors: %d\n", rs->ndet);
429   printf ("    Number of views: %d\n", rs->nview);
430   printf ("             Remark: %s\n", rs->remark);
431   printf ("Piclen: %f\n", rs->piclen);
432   printf ("det_start: %f\n", rs->det_start);
433   printf ("det_inc: %f\n", rs->det_inc);
434   printf ("rot_start: %f\n", rs->rot_start);
435   printf ("rot_inc: %f\n", rs->rot_inc);
436 }