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