1 /*****************************************************************************
2 ** This is part of the CTSim program
3 ** Copyright (C) 1983-2000 Kevin Rosenberg
5 ** $Id: ctrec.c,v 1.3 2000/04/30 04:06:13 kevin Exp $
7 ** Revision 1.3 2000/04/30 04:06:13 kevin
8 ** Update Raysum i/o routines
9 ** Fix MPI bug in ctrec (scatter_raysum) that referenced rs_global
11 ** Revision 1.2 2000/04/29 23:24:56 kevin
12 ** *** empty log message ***
14 ** Revision 1.1.1.1 2000/04/28 13:02:44 kevin
15 ** Initial CVS import for first public release
19 ** This program is free software; you can redistribute it and/or modify
20 ** it under the terms of the GNU General Public License (version 2) as
21 ** published by the Free Software Foundation.
23 ** This program is distributed in the hope that it will be useful,
24 ** but WITHOUT ANY WARRANTY; without even the implied warranty of
25 ** MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
26 ** GNU General Public License for more details.
28 ** You should have received a copy of the GNU General Public License
29 ** along with this program; if not, write to the Free Software
30 ** Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
31 ******************************************************************************/
34 * ctrec.c Reconstruct an image from raysums
38 * Jul 99 -- Converted to ANSI C
39 * Added MPI parallel processing
46 #define O_FILTER_PARAM 3
54 static struct option my_options[] =
56 {"interp", 1, 0, O_INTERP},
57 {"filter", 1, 0, O_FILTER},
58 {"filter-param", 1, 0, O_FILTER_PARAM},
59 {"backproj", 1, 0, O_BACKPROJ},
60 {"trace", 1, 0, O_TRACE},
61 {"debug", 0, 0, O_DEBUG},
62 {"verbose", 0, 0, O_VERBOSE},
63 {"help", 0, 0, O_HELP},
64 {"version", 0, 0, O_VERSION},
69 usage (const char *program)
71 fprintf(stdout,"usage: %s raysum-file image-file nx-image ny-image [OPTIONS]\n", kbasename(program));
72 fprintf(stdout,"Image reconstruction from raysum projections\n");
74 fprintf(stdout," raysum-file Input raysum file\n");
75 fprintf(stdout," image-file Output image file in SDF2D format\n");
76 fprintf(stdout," nx-image Number of columns in output image\n");
77 fprintf(stdout," ny-image Number of rows in output image\n");
78 fprintf(stdout," --interp Interpolation method during backprojection\n");
79 fprintf(stdout," nearest Nearest neighbor interpolation\n");
80 fprintf(stdout," linear Linear interpolation\n");
81 fprintf(stdout," bspline B-spline interpolation\n");
82 fprintf(stdout," --filter Filter name\n");
83 fprintf(stdout," abs_bandlimit Abs * Bandlimiting (default)\n");
84 fprintf(stdout," abs_sinc Abs * Sinc\n");
85 fprintf(stdout," abs_cos Abs * Cosine\n");
86 fprintf(stdout," abs_hamming Abs * Hamming\n");
87 fprintf(stdout," shepp Shepp-Logan\n");
88 fprintf(stdout," bandlimit Bandlimiting\n");
89 fprintf(stdout," sinc Sinc\n");
90 fprintf(stdout," cos Cosine\n");
91 fprintf(stdout," triangle Triangle\n");
92 fprintf(stdout," hamming Hamming\n");
93 fprintf(stdout," --backproj Backprojection Method\n");
94 fprintf(stdout," trig Trigometric functions at every point\n");
95 fprintf(stdout," table Trigometric functions with precalculated table\n");
96 fprintf(stdout," diff Difference method\n");
97 fprintf(stdout," diff2 Optimized difference method (default)\n");
98 fprintf(stdout," idiff2 Optimized difference method with integer math\n");
99 fprintf(stdout," --filter-param Alpha level for Hamming filter\n");
100 fprintf(stdout," --trace Set tracing to level\n");
101 fprintf(stdout," none No tracing (default)\n");
102 fprintf(stdout," status Status tracing\n");
103 fprintf(stdout," pic Trace picture\n");
104 fprintf(stdout," rays Trace allrays\n");
105 fprintf(stdout," plot Trace plotting\n");
106 fprintf(stdout," clipping Trace clipping\n");
107 fprintf(stdout," --verbose Turn on verbose mode\n");
108 fprintf(stdout," --debug Turn on debug mode\n");
109 fprintf(stdout," --version Print version\n");
110 fprintf(stdout," --help Print this help message\n");
116 static void mpi_scatter_rs (RAYSUM *rs_global, RAYSUM *rs_local, const int debug);
118 static void print_raysum_info(const RAYSUM *rs);
121 main (const int argc, char *const argv[])
123 IMAGE *im_global = NULL;
124 RAYSUM *rs_global = NULL;
125 char *rs_name, *im_filename = NULL;
126 char remark[MAXREMARK];
129 double time_start = 0, time_end = 0;
133 int opt_trace = TRACE_NONE;
134 double opt_filter_param = 1;
135 int opt_filter = W_A_BANDLIMIT;
136 int opt_interp = I_LINEAR;
137 int opt_interp_param = 1;
138 int opt_backproj = O_BPROJ_DIFF2;
143 char **mpi_argv = (char **) argv;
144 int mpi_nview, mpi_ndet;
145 double mpi_detinc, mpi_rotinc, mpi_piclen;
146 double mpi_t1, mpi_t2, mpi_t, mpi_t_g;
149 MPI_Init(&mpi_argc, &mpi_argv);
150 MPI_Comm_dup (MPI_COMM_WORLD, &mpi_ct.comm);
151 MPI_Comm_size(mpi_ct.comm, &mpi_ct.nproc);
152 MPI_Comm_rank(mpi_ct.comm, &mpi_ct.my_rank);
154 if (mpi_ct.nproc > MPI_MAX_PROCESS) {
155 sys_error(ERR_FATAL, "Number of mpi processes (%d) exceeds max processes (%d)",
156 mpi_ct.nproc, MPI_MAX_PROCESS);
162 time_start = MPI_Wtime();
164 time_start = td_current_sec();
168 if (mpi_ct.my_rank == 0) {
171 int c = getopt_long(argc, argv, "", my_options, NULL);
180 opt_interp = opt_set_interpolation(optarg, argv[0]);
183 opt_filter = opt_set_filter(optarg, argv[0]);
186 opt_backproj = opt_set_backproj(optarg, argv[0]);
189 opt_filter_param = strtod(optarg, &endptr);
190 if (endptr != optarg + strlen(optarg)) {
202 opt_trace = opt_set_trace(optarg, argv[0]);
206 fprintf(stdout, "Version %s\n", VERSION);
208 fprintf(stderr, "Unknown version number");
221 if (optind + 4 != argc) {
226 rs_name = argv[optind];
228 im_filename = argv[optind + 1];
230 nx = strtol(argv[optind + 2], &endptr, 10);
231 ny = strtol(argv[optind + 3], &endptr, 10);
233 if (opt_filter == W_G_HAMMING || opt_filter == W_AG_HAMMING)
234 sprintf (filt_name, "%s: alpha = %.2f",
235 filter_name_of (opt_filter), opt_filter_param);
237 sprintf (filt_name, "%s", filter_name_of (opt_filter));
240 "Reconstruct: %dx%d, %s, %s, %s",
241 nx, ny, filt_name, interp_name_of (opt_interp), name_of_backproj(opt_backproj));
244 fprintf (stdout, "%s\n", remark);
250 if (mpi_ct.my_rank == 0) {
251 rs_global = raysum_open (rs_name);
252 raysum_read (rs_global);
253 print_raysum_info(rs_global);
255 printf ("Number of detectors: %d\n", rs_global->ndet);
256 printf (" Number of views: %d\n", rs_global->nview);
257 printf (" Remark: %s\n", rs_global->remark);
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;
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);
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);
291 mpi_ct_calc_work_units(mpi_nview);
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);
299 rs_local = raysum_create (NULL, mpi_ct.local_work_units[mpi_ct.my_rank], mpi_ndet);
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 MPI_Barrier(mpi_ct.comm);
305 rs_local->ndet = mpi_ndet;
306 rs_local->nview = mpi_nview;
307 rs_local->det_inc = mpi_detinc;
308 rs_local->piclen = mpi_piclen;
309 rs_local->rot_inc = mpi_rotinc;
312 mpi_t1 = MPI_Wtime();
313 mpi_scatter_rs(rs_global, rs_local, opt_debug);
315 mpi_t2 = MPI_Wtime();
316 mpi_t = mpi_t2 - mpi_t1;
317 MPI_Reduce(&mpi_t, &mpi_t_g, 1, MPI_DOUBLE, MPI_MAX, 0, mpi_ct.comm);
318 if (mpi_ct.my_rank == 0)
319 printf("Time to scatter rs = %f secs, Max time = %f sec\n", mpi_t, mpi_t_g);
322 if (mpi_ct.my_rank == 0) {
323 im_global = image_create (im_filename, nx, ny);
324 sdf_add_label (LT_HISTORY, rs_global->remark, rs_global->calctime, im_global->dfp_2d->dfp);
325 sdf_add_empty_label (im_global->dfp_2d->dfp);
327 im_local = image_create (NULL, nx, ny);
330 rs_global = raysum_open (rs_name);
331 raysum_read (rs_global);
333 print_raysum_info(rs_global);
335 im_global = image_create (im_filename, nx, ny);
336 sdf_add_label (LT_HISTORY, rs_global->remark, rs_global->calctime, im_global->dfp_2d->dfp);
337 sdf_add_empty_label (im_global->dfp_2d->dfp);
341 mpi_t1 = MPI_Wtime();
342 image_reconst (im_local, rs_local, opt_filter, opt_filter_param,
343 opt_interp, opt_interp_param, opt_backproj, opt_trace);
344 mpi_t2 = MPI_Wtime();
345 mpi_t = mpi_t2 - mpi_t1;
346 MPI_Reduce(&mpi_t, &mpi_t_g, 1, MPI_DOUBLE, MPI_MAX, 0, mpi_ct.comm);
347 if (mpi_ct.my_rank == 0) {
348 printf("Time to reconstruct = %f, Max time = %f\n", mpi_t, mpi_t_g);
351 image_reconst (im_global, rs_global, opt_filter, opt_filter_param,
352 opt_interp, opt_interp_param, opt_backproj, opt_trace);
356 if (mpi_ct.my_rank == 0) {
357 raysum_close (rs_global);
361 mpi_t1 = MPI_Wtime();
362 for (ix = 0; ix < im_local->nx; ix++) {
363 MPI_Reduce(im_local->v[ix], im_global->v[ix],
364 im_local->ny, MPI_FLOAT, MPI_SUM, 0, mpi_ct.comm);
367 mpi_t2 = MPI_Wtime();
368 mpi_t = mpi_t2 - mpi_t1;
369 MPI_Reduce(&mpi_t, &mpi_t_g, 1, MPI_DOUBLE, MPI_MAX, 0, mpi_ct.comm);
370 if (mpi_ct.my_rank == 0)
371 printf("Time to reduce image = %f secs, max time = %f\n", mpi_t, mpi_t_g);
373 if (mpi_ct.my_rank == 0) {
374 strncpy (im_global->remark, remark, MAXREMARK);
375 image_save (im_global);
376 time_end = MPI_Wtime();
377 im_global->calctime = time_end - time_start;
379 fprintf (stdout, "Time active = %.2f\n", im_global->calctime);
382 raysum_close (rs_global);
383 strncpy (im_global->remark, remark, MAXREMARK);
384 image_save (im_global);
386 time_end = td_current_sec();
387 im_global->calctime = time_end - time_start;
389 fprintf (stdout, "Time active = %.2f secs\n", im_global->calctime);
401 static void mpi_scatter_rs (RAYSUM *rs_global, RAYSUM *rs_local, const int opt_debug)
407 if (mpi_ct.my_rank == 0) {
408 for (iproc = 0; iproc < mpi_ct.nproc; iproc++) {
409 end_work_unit = mpi_ct.start_work_unit[iproc]
410 + mpi_ct.local_work_units[iproc]
414 fprintf(stdout, "Sending rs data to process %d\n", iproc);
416 for (iw = mpi_ct.start_work_unit[iproc]; iw <= end_work_unit; iw++) {
418 fprintf(stdout, "Sending from process %2d to process %2d: view_angle=%8f, ndet=%5d\n", mpi_ct.my_rank, iproc, rs_global->view[iw]->view_angle, rs_global->view[iw]->ndet);
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);
427 fprintf(stdout, "Receiving rs data in process %d\n", mpi_ct.my_rank);
428 MPI_Barrier(mpi_ct.comm);
431 end_work_unit = mpi_ct.local_work_units[mpi_ct.my_rank] - 1;
432 for (iw = 0; iw <= end_work_unit; iw++) {
436 fprintf(stdout,"Receiving into rs_local from mpi_scatter_rs, process %2d: rs_local=%lx, ", mpi_ct.my_rank, (unsigned long int) rs_local);
437 fprintf(stdout,"iw=%4d, nrot=%4d, ndet=%4d, view=%lx, detval=%lx\n", iw, rs_local->nview, rs_local->ndet, (unsigned long int) rs_local->view[iw], (unsigned long int) rs_local->view[iw]->detval);
440 MPI_Recv(&rs_local->view[iw]->ndet, 1, MPI_INT, 0, 0, mpi_ct.comm, &status);
441 MPI_Recv(&rs_local->view[iw]->view_angle, 1, MPI_DOUBLE, 0, 0, mpi_ct.comm, &status);
442 MPI_Recv(rs_local->view[iw]->detval, rs_local->ndet, MPI_FLOAT, 0, 0, mpi_ct.comm, &status);
445 fprintf(stdout, "Received view_angle=%8f, ndet=%5d\n", rs_local->view[iw]->view_angle, rs_local->view[iw]->ndet);
448 rs_local->nview = mpi_ct.local_work_units[mpi_ct.my_rank];
450 MPI_Barrier(MPI_COMM_WORLD);
451 fprintf(stdout, "Done with mpi_scatter_rs in process %2d\n", mpi_ct.my_rank);
458 static void print_raysum_info(const RAYSUM *rs)
460 printf ("Number of detectors: %d\n", rs->ndet);
461 printf (" Number of views: %d\n", rs->nview);
462 printf (" Remark: %s\n", rs->remark);
463 printf ("Piclen: %f\n", rs->piclen);
464 printf ("det_start: %f\n", rs->det_start);
465 printf ("det_inc: %f\n", rs->det_inc);
466 printf ("rot_start: %f\n", rs->rot_start);
467 printf ("rot_inc: %f\n", rs->rot_inc);