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