1 /*****************************************************************************
2 ** This is part of the CTSim program
3 ** Copyright (C) 1983-2000 Kevin Rosenberg
5 ** $Id: ctrec.c,v 1.9 2000/05/08 20:02:32 kevin Exp $
7 ** Revision 1.9 2000/05/08 20:02:32 kevin
10 ** Revision 1.8 2000/05/04 18:16:34 kevin
11 ** renamed filter definitions
13 ** Revision 1.7 2000/05/03 08:49:50 kevin
16 ** Revision 1.6 2000/05/02 15:31:47 kevin
19 ** Revision 1.5 2000/04/30 11:41:06 kevin
20 ** Cleaned up debugging code
22 ** Revision 1.4 2000/04/30 10:13:27 kevin
25 ** Revision 1.3 2000/04/30 04:06:13 kevin
26 ** Update Raysum i/o routines
27 ** Fix MPI bug in ctrec (scatter_raysum) that referenced rs_global
29 ** Revision 1.2 2000/04/29 23:24:56 kevin
30 ** *** empty log message ***
32 ** Revision 1.1.1.1 2000/04/28 13:02:44 kevin
33 ** Initial CVS import for first public release
37 ** This program is free software; you can redistribute it and/or modify
38 ** it under the terms of the GNU General Public License (version 2) as
39 ** published by the Free Software Foundation.
41 ** This program is distributed in the hope that it will be useful,
42 ** but WITHOUT ANY WARRANTY; without even the implied warranty of
43 ** MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
44 ** GNU General Public License for more details.
46 ** You should have received a copy of the GNU General Public License
47 ** along with this program; if not, write to the Free Software
48 ** Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
49 ******************************************************************************/
52 * ctrec.c Reconstruct an image from raysums
56 * Jul 99 -- Converted to ANSI C
57 * Added MPI parallel processing
62 enum {O_INTERP, O_FILTER, O_FILTER_PARAM, O_BACKPROJ, O_VERBOSE, O_TRACE, O_HELP, O_DEBUG, O_VERSION};
64 static struct option my_options[] =
66 {"interp", 1, 0, O_INTERP},
67 {"filter", 1, 0, O_FILTER},
68 {"filter-param", 1, 0, O_FILTER_PARAM},
69 {"backproj", 1, 0, O_BACKPROJ},
70 {"trace", 1, 0, O_TRACE},
71 {"debug", 0, 0, O_DEBUG},
72 {"verbose", 0, 0, O_VERBOSE},
73 {"help", 0, 0, O_HELP},
74 {"version", 0, 0, O_VERSION},
79 usage (const char *program)
81 fprintf(stdout,"usage: %s raysum-file image-file nx-image ny-image [OPTIONS]\n", kbasename(program));
82 fprintf(stdout,"Image reconstruction from raysum projections\n");
84 fprintf(stdout," raysum-file Input raysum file\n");
85 fprintf(stdout," image-file Output image file in SDF2D format\n");
86 fprintf(stdout," nx-image Number of columns in output image\n");
87 fprintf(stdout," ny-image Number of rows in output image\n");
88 fprintf(stdout," --interp Interpolation method during backprojection\n");
89 fprintf(stdout," nearest Nearest neighbor interpolation\n");
90 fprintf(stdout," linear Linear interpolation\n");
91 fprintf(stdout," bspline B-spline interpolation\n");
92 fprintf(stdout," --filter Filter name\n");
93 fprintf(stdout," abs_bandlimit Abs * Bandlimiting (default)\n");
94 fprintf(stdout," abs_sinc Abs * Sinc\n");
95 fprintf(stdout," abs_cos Abs * Cosine\n");
96 fprintf(stdout," abs_hamming Abs * Hamming\n");
97 fprintf(stdout," shepp Shepp-Logan\n");
98 fprintf(stdout," bandlimit Bandlimiting\n");
99 fprintf(stdout," sinc Sinc\n");
100 fprintf(stdout," cos Cosine\n");
101 fprintf(stdout," triangle Triangle\n");
102 fprintf(stdout," hamming Hamming\n");
103 fprintf(stdout," --backproj Backprojection Method\n");
104 fprintf(stdout," trig Trigometric functions at every point\n");
105 fprintf(stdout," table Trigometric functions with precalculated table\n");
106 fprintf(stdout," diff Difference method\n");
107 fprintf(stdout," diff2 Optimized difference method (default)\n");
108 fprintf(stdout," idiff2 Optimized difference method with integer math\n");
109 fprintf(stdout," --filter-param Alpha level for Hamming filter\n");
110 fprintf(stdout," --trace Set tracing to level\n");
111 fprintf(stdout," none No tracing (default)\n");
112 fprintf(stdout," status Status tracing\n");
113 fprintf(stdout," phm Trace phantom\n");
114 fprintf(stdout," rays Trace allrays\n");
115 fprintf(stdout," plot Trace plotting\n");
116 fprintf(stdout," clipping Trace clipping\n");
117 fprintf(stdout," --verbose Turn on verbose mode\n");
118 fprintf(stdout," --debug Turn on debug mode\n");
119 fprintf(stdout," --version Print version\n");
120 fprintf(stdout," --help Print this help message\n");
126 static void mpi_scatter_rs (RAYSUM *rs_global, RAYSUM *rs_local, const int debug);
128 static void print_raysum_info(const RAYSUM *rs);
131 main (const int argc, char *const argv[])
133 IMAGE *im_global = NULL;
134 RAYSUM *rs_global = NULL;
135 char *rs_name, *im_filename = NULL;
136 char remark[MAXREMARK];
139 double time_start = 0, time_end = 0;
143 int opt_trace = TRACE_NONE;
144 double opt_filter_param = 1;
145 int opt_filter = FILTER_ABS_BANDLIMIT;
146 int opt_interp = I_LINEAR;
147 int opt_interp_param = 1;
148 BackprojType opt_backproj = O_BPROJ_DIFF2;
153 char **mpi_argv = (char **) argv;
154 int mpi_nview, mpi_ndet;
155 double mpi_detinc, mpi_rotinc, mpi_phmlen;
156 double mpi_t1, mpi_t2, mpi_t, mpi_t_g;
159 MPI_Init(&mpi_argc, &mpi_argv);
160 MPI_Comm_dup (MPI_COMM_WORLD, &mpi_ct.comm);
161 MPI_Comm_size(mpi_ct.comm, &mpi_ct.nproc);
162 MPI_Comm_rank(mpi_ct.comm, &mpi_ct.my_rank);
164 if (mpi_ct.nproc > MPI_MAX_PROCESS) {
165 sys_error(ERR_FATAL, "Number of mpi processes (%d) exceeds max processes (%d)",
166 mpi_ct.nproc, MPI_MAX_PROCESS);
172 time_start = MPI_Wtime();
174 time_start = td_current_sec();
178 if (mpi_ct.my_rank == 0) {
181 int c = getopt_long(argc, argv, "", my_options, NULL);
190 opt_interp = opt_set_interpolation(optarg, argv[0]);
193 opt_filter = opt_set_filter(optarg, argv[0]);
196 opt_backproj = opt_set_backproj(optarg, argv[0]);
199 opt_filter_param = strtod(optarg, &endptr);
200 if (endptr != optarg + strlen(optarg)) {
212 opt_trace = opt_set_trace(optarg, argv[0]);
216 fprintf(stdout, "Version %s\n", VERSION);
218 fprintf(stderr, "Unknown version number");
231 if (optind + 4 != argc) {
236 rs_name = argv[optind];
238 im_filename = argv[optind + 1];
240 nx = strtol(argv[optind + 2], &endptr, 10);
241 ny = strtol(argv[optind + 3], &endptr, 10);
243 if (opt_filter == FILTER_G_HAMMING || opt_filter == FILTER_ABS_G_HAMMING)
244 sprintf (filt_name, "%s: alpha = %.2f",
245 name_of_filter (opt_filter), opt_filter_param);
247 sprintf (filt_name, "%s", name_of_filter (opt_filter));
249 sprintf (remark, "Reconstruct: %dx%d, %s, %s, %s",
250 nx, ny, filt_name, name_of_interpolation (opt_interp), name_of_backproj(opt_backproj));
253 fprintf (stdout, "%s\n", remark);
259 if (mpi_ct.my_rank == 0) {
260 rs_global = raysum_open (rs_name);
261 raysum_read (rs_global);
263 print_raysum_info(rs_global);
265 mpi_ndet = rs_global->ndet;
266 mpi_nview = rs_global->nview;
267 mpi_detinc = rs_global->det_inc;
268 mpi_phmlen = rs_global->phmlen;
269 mpi_rotinc = rs_global->rot_inc;
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_phmlen, 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);
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);
296 mpi_ct_calc_work_units(mpi_nview);
298 rs_local = raysum_create (NULL, mpi_ct.local_work_units[mpi_ct.my_rank], mpi_ndet);
300 rs_local->ndet = mpi_ndet;
301 rs_local->nview = mpi_nview;
302 rs_local->det_inc = mpi_detinc;
303 rs_local->phmlen = mpi_phmlen;
304 rs_local->rot_inc = mpi_rotinc;
307 mpi_t1 = MPI_Wtime();
308 mpi_scatter_rs(rs_global, rs_local, opt_debug);
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);
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);
322 im_local = image_create (NULL, nx, ny);
325 rs_global = raysum_open (rs_name);
326 raysum_read (rs_global);
328 print_raysum_info(rs_global);
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);
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);
345 image_reconst (im_global, rs_global, opt_filter, opt_filter_param,
346 opt_interp, opt_interp_param, opt_backproj, opt_trace);
350 if (mpi_ct.my_rank == 0) {
351 raysum_close (rs_global);
355 mpi_t1 = MPI_Wtime();
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];
362 MPI_Reduce(im_local->v[ix], recvbuf, im_local->ny, MPI_FLOAT, MPI_SUM, 0, mpi_ct.comm);
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);
372 if (mpi_ct.my_rank == 0) {
373 strncpy (im_global->remark, remark, sizeof(im_global->remark));
374 time_end = MPI_Wtime();
375 im_global->calctime = time_end - time_start;
376 image_save (im_global);
378 fprintf (stdout, "Time active = %.2f\n", im_global->calctime);
381 raysum_close (rs_global);
382 strncpy (im_global->remark, remark, sizeof(im_global->remark));
383 time_end = td_current_sec();
384 im_global->calctime = time_end - time_start;
385 image_save (im_global);
387 fprintf (stdout, "Time active = %.2f secs\n", im_global->calctime);
399 static void mpi_scatter_rs (RAYSUM *rs_global, RAYSUM *rs_local, const int opt_debug)
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;
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);
418 fprintf(stdout, "Receiving rs data in process %d\n", mpi_ct.my_rank);
420 end_work_unit = mpi_ct.local_work_units[mpi_ct.my_rank] - 1;
421 for (iw = 0; iw <= end_work_unit; iw++) {
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);
428 rs_local->nview = mpi_ct.local_work_units[mpi_ct.my_rank];
430 MPI_Barrier(MPI_COMM_WORLD);
431 fprintf(stdout, "Done with mpi_scatter_rs in process %2d\n", mpi_ct.my_rank);
437 static void print_raysum_info(const RAYSUM *rs)
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 (" phmlen: %f\n", rs->phmlen);
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);