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