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