1 /*****************************************************************************
2 ** This is part of the CTSim program
3 ** Copyright (C) 1983-2000 Kevin Rosenberg
5 ** $Id: ctrec.c,v 1.12 2000/05/24 22:50:04 kevin Exp $
7 ** Revision 1.12 2000/05/24 22:50:04 kevin
8 ** Added support for new SGP library
10 ** Revision 1.11 2000/05/16 04:33:59 kevin
11 ** Improved option processing
13 ** Revision 1.10 2000/05/11 01:06:30 kevin
14 ** Changed sprintf to snprintf
16 ** Revision 1.9 2000/05/08 20:02:32 kevin
19 ** Revision 1.8 2000/05/04 18:16:34 kevin
20 ** renamed filter definitions
22 ** Revision 1.7 2000/05/03 08:49:50 kevin
25 ** Revision 1.6 2000/05/02 15:31:47 kevin
28 ** Revision 1.5 2000/04/30 11:41:06 kevin
29 ** Cleaned up debugging code
31 ** Revision 1.4 2000/04/30 10:13:27 kevin
34 ** Revision 1.3 2000/04/30 04:06:13 kevin
35 ** Update Raysum i/o routines
36 ** Fix MPI bug in ctrec (scatter_raysum) that referenced rs_global
38 ** Revision 1.2 2000/04/29 23:24:56 kevin
39 ** *** empty log message ***
41 ** Revision 1.1.1.1 2000/04/28 13:02:44 kevin
42 ** Initial CVS import for first public release
46 ** This program is free software; you can redistribute it and/or modify
47 ** it under the terms of the GNU General Public License (version 2) as
48 ** published by the Free Software Foundation.
50 ** This program is distributed in the hope that it will be useful,
51 ** but WITHOUT ANY WARRANTY; without even the implied warranty of
52 ** MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
53 ** GNU General Public License for more details.
55 ** You should have received a copy of the GNU General Public License
56 ** along with this program; if not, write to the Free Software
57 ** Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
58 ******************************************************************************/
61 * ctrec.c Reconstruct an image from raysums
65 * Jul 99 -- Converted to ANSI C
66 * Added MPI parallel processing
71 enum {O_INTERP, O_FILTER, O_FILTER_PARAM, O_BACKPROJ, O_VERBOSE, O_TRACE, O_HELP, O_DEBUG, O_VERSION};
73 static struct option my_options[] =
75 {"interp", 1, 0, O_INTERP},
76 {"filter", 1, 0, O_FILTER},
77 {"filter-param", 1, 0, O_FILTER_PARAM},
78 {"backproj", 1, 0, O_BACKPROJ},
79 {"trace", 1, 0, O_TRACE},
80 {"debug", 0, 0, O_DEBUG},
81 {"verbose", 0, 0, O_VERBOSE},
82 {"help", 0, 0, O_HELP},
83 {"version", 0, 0, O_VERSION},
88 ctrec_usage (const char *program)
90 fprintf(stdout,"usage: %s raysum-file image-file nx-image ny-image [OPTIONS]\n", kbasename(program));
91 fprintf(stdout,"Image reconstruction from raysum projections\n");
93 fprintf(stdout," raysum-file Input raysum file\n");
94 fprintf(stdout," image-file Output image file in SDF2D format\n");
95 fprintf(stdout," nx-image Number of columns in output image\n");
96 fprintf(stdout," ny-image Number of rows in output image\n");
97 fprintf(stdout," --interp Interpolation method during backprojection\n");
98 fprintf(stdout," nearest Nearest neighbor interpolation\n");
99 fprintf(stdout," linear Linear interpolation\n");
100 fprintf(stdout," bspline B-spline interpolation\n");
101 fprintf(stdout," --filter Filter name\n");
102 fprintf(stdout," abs_bandlimit Abs * Bandlimiting (default)\n");
103 fprintf(stdout," abs_sinc Abs * Sinc\n");
104 fprintf(stdout," abs_cos Abs * Cosine\n");
105 fprintf(stdout," abs_hamming Abs * Hamming\n");
106 fprintf(stdout," shepp Shepp-Logan\n");
107 fprintf(stdout," bandlimit Bandlimiting\n");
108 fprintf(stdout," sinc Sinc\n");
109 fprintf(stdout," cos Cosine\n");
110 fprintf(stdout," triangle Triangle\n");
111 fprintf(stdout," hamming Hamming\n");
112 fprintf(stdout," --backproj Backprojection Method\n");
113 fprintf(stdout," trig Trigometric functions at every point\n");
114 fprintf(stdout," table Trigometric functions with precalculated table\n");
115 fprintf(stdout," diff Difference method\n");
116 fprintf(stdout," diff2 Optimized difference method (default)\n");
117 fprintf(stdout," idiff2 Optimized difference method with integer math\n");
118 fprintf(stdout," --filter-param Alpha level for Hamming filter\n");
119 fprintf(stdout," --trace Set tracing to level\n");
120 fprintf(stdout," none No tracing (default)\n");
121 fprintf(stdout," text Text level tracing\n");
122 fprintf(stdout," phm Trace phantom\n");
123 fprintf(stdout," rays Trace allrays\n");
124 fprintf(stdout," plot Trace plotting\n");
125 fprintf(stdout," clipping Trace clipping\n");
126 fprintf(stdout," --verbose Turn on verbose mode\n");
127 fprintf(stdout," --debug Turn on debug mode\n");
128 fprintf(stdout," --version Print version\n");
129 fprintf(stdout," --help Print this help message\n");
134 static void mpi_scatter_rs (RAYSUM *rs_global, RAYSUM *rs_local, const int debug);
137 static void print_raysum_info(const RAYSUM *rs);
140 ctrec_main (const int argc, char *const argv[])
142 IMAGE *im_global = NULL;
143 RAYSUM *rs_global = NULL;
144 char *rs_name, *im_filename = NULL;
145 char remark[MAXREMARK];
148 double time_start = 0, time_end = 0;
152 int opt_trace = TRACE_NONE;
153 double opt_filter_param = 1;
154 int opt_filter = FILTER_ABS_BANDLIMIT;
155 int opt_interp = I_LINEAR;
156 int opt_interp_param = 1;
157 BackprojType opt_backproj = O_BPROJ_DIFF2;
162 char **mpi_argv = (char **) argv;
163 int mpi_nview, mpi_ndet;
164 double mpi_detinc, mpi_rotinc, mpi_phmlen;
165 double mpi_t1, mpi_t2, mpi_t, mpi_t_g;
168 MPI_Init(&mpi_argc, &mpi_argv);
169 MPI_Comm_dup (MPI_COMM_WORLD, &mpi_ct.comm);
170 MPI_Comm_size(mpi_ct.comm, &mpi_ct.nproc);
171 MPI_Comm_rank(mpi_ct.comm, &mpi_ct.my_rank);
173 if (mpi_ct.nproc > MPI_MAX_PROCESS) {
174 sys_error(ERR_FATAL, "Number of mpi processes (%d) exceeds max processes (%d)",
175 mpi_ct.nproc, MPI_MAX_PROCESS);
181 time_start = MPI_Wtime();
183 time_start = td_current_sec();
187 if (mpi_ct.my_rank == 0) {
190 int c = getopt_long(argc, argv, "", my_options, NULL);
199 if ((opt_interp = opt_set_interpolation(optarg)) < 0) {
200 ctrec_usage(argv[0]);
205 if ((opt_filter = opt_set_filter(optarg)) < 0) {
206 ctrec_usage(argv[0]);
211 if ((opt_backproj = opt_set_backproj(optarg)) < 0) {
212 ctrec_usage(argv[0]);
217 opt_filter_param = strtod(optarg, &endptr);
218 if (endptr != optarg + strlen(optarg)) {
219 ctrec_usage(argv[0]);
229 if ((opt_trace = opt_set_trace(optarg)) < 0) {
230 ctrec_usage(argv[0]);
236 fprintf(stdout, "Version %s\n", VERSION);
238 fprintf(stderr, "Unknown version number");
243 ctrec_usage(argv[0]);
246 ctrec_usage(argv[0]);
251 if (optind + 4 != argc) {
252 ctrec_usage(argv[0]);
256 rs_name = argv[optind];
258 im_filename = argv[optind + 1];
260 nx = strtol(argv[optind + 2], &endptr, 10);
261 ny = strtol(argv[optind + 3], &endptr, 10);
263 if (opt_filter == FILTER_G_HAMMING || opt_filter == FILTER_ABS_G_HAMMING)
264 snprintf (filt_name, sizeof(filt_name), "%s: alpha = %.2f",
265 name_of_filter (opt_filter), opt_filter_param);
267 snprintf (filt_name, sizeof(filt_name), "%s", name_of_filter (opt_filter));
269 snprintf (remark, sizeof(remark), "Reconstruct: %dx%d, %s, %s, %s",
270 nx, ny, filt_name, name_of_interpolation (opt_interp), name_of_backproj(opt_backproj));
273 fprintf (stdout, "%s\n", remark);
279 if (mpi_ct.my_rank == 0) {
280 rs_global = raysum_open (rs_name);
281 raysum_read (rs_global);
283 print_raysum_info(rs_global);
285 mpi_ndet = rs_global->ndet;
286 mpi_nview = rs_global->nview;
287 mpi_detinc = rs_global->det_inc;
288 mpi_phmlen = rs_global->phmlen;
289 mpi_rotinc = rs_global->rot_inc;
292 mpi_t1 = MPI_Wtime();
293 MPI_Bcast(&opt_verbose, 1, MPI_INT, 0, mpi_ct.comm);
294 MPI_Bcast(&opt_debug, 1, MPI_INT, 0, mpi_ct.comm);
295 MPI_Bcast(&opt_trace, 1, MPI_INT, 0, mpi_ct.comm);
296 MPI_Bcast(&opt_filter, 1, MPI_INT, 0, mpi_ct.comm);
297 MPI_Bcast(&opt_interp, 1, MPI_INT, 0, mpi_ct.comm);
298 MPI_Bcast(&opt_filter_param, 1, MPI_DOUBLE, 0, mpi_ct.comm);
299 MPI_Bcast(&opt_interp_param, 1, MPI_INT, 0, mpi_ct.comm);
300 MPI_Bcast(&opt_backproj, 1, MPI_INT, 0, mpi_ct.comm);
301 MPI_Bcast(&mpi_ndet, 1, MPI_INT, 0, mpi_ct.comm);
302 MPI_Bcast(&mpi_nview, 1, MPI_INT, 0, mpi_ct.comm);
303 MPI_Bcast(&mpi_detinc, 1, MPI_DOUBLE, 0, mpi_ct.comm);
304 MPI_Bcast(&mpi_phmlen, 1, MPI_DOUBLE, 0, mpi_ct.comm);
305 MPI_Bcast(&mpi_rotinc, 1, MPI_DOUBLE, 0, mpi_ct.comm);
306 MPI_Bcast(&nx, 1, MPI_INT, 0, mpi_ct.comm);
307 MPI_Bcast(&ny, 1, MPI_INT, 0, mpi_ct.comm);
309 mpi_t2 = MPI_Wtime();
310 mpi_t = mpi_t2 - mpi_t1;
311 MPI_Reduce(&mpi_t, &mpi_t_g, 1, MPI_DOUBLE, MPI_MAX, 0, mpi_ct.comm);
312 if (mpi_ct.my_rank == 0)
313 printf("Time to Bcast vars = %f secs, Max time = %f\n", mpi_t, mpi_t_g);
316 mpi_ct_calc_work_units(mpi_nview);
318 rs_local = raysum_create (NULL, mpi_ct.local_work_units[mpi_ct.my_rank], mpi_ndet);
320 rs_local->ndet = mpi_ndet;
321 rs_local->nview = mpi_nview;
322 rs_local->det_inc = mpi_detinc;
323 rs_local->phmlen = mpi_phmlen;
324 rs_local->rot_inc = mpi_rotinc;
327 mpi_t1 = MPI_Wtime();
328 mpi_scatter_rs(rs_global, rs_local, opt_debug);
330 mpi_t2 = MPI_Wtime();
331 mpi_t = mpi_t2 - mpi_t1;
332 MPI_Reduce(&mpi_t, &mpi_t_g, 1, MPI_DOUBLE, MPI_MAX, 0, mpi_ct.comm);
333 if (mpi_ct.my_rank == 0)
334 printf("Time to scatter rs = %f secs, Max time = %f sec\n", mpi_t, mpi_t_g);
337 if (mpi_ct.my_rank == 0) {
338 im_global = image_create (im_filename, nx, ny);
339 sdf_add_label (LT_HISTORY, rs_global->remark, rs_global->calctime, im_global->dfp_2d->dfp);
340 sdf_add_empty_label (im_global->dfp_2d->dfp);
342 im_local = image_create (NULL, nx, ny);
345 rs_global = raysum_open (rs_name);
346 raysum_read (rs_global);
348 print_raysum_info(rs_global);
350 im_global = image_create (im_filename, nx, ny);
351 sdf_add_label (LT_HISTORY, rs_global->remark, rs_global->calctime, im_global->dfp_2d->dfp);
352 sdf_add_empty_label (im_global->dfp_2d->dfp);
356 mpi_t1 = MPI_Wtime();
357 image_reconst (im_local, rs_local, opt_filter, opt_filter_param,
358 opt_interp, opt_interp_param, opt_backproj, opt_trace);
359 mpi_t2 = MPI_Wtime();
360 mpi_t = mpi_t2 - mpi_t1;
361 MPI_Reduce(&mpi_t, &mpi_t_g, 1, MPI_DOUBLE, MPI_MAX, 0, mpi_ct.comm);
362 if (mpi_ct.my_rank == 0 && opt_verbose)
363 printf("Time to reconstruct = %f, Max time = %f\n", mpi_t, mpi_t_g);
365 image_reconst (im_global, rs_global, opt_filter, opt_filter_param,
366 opt_interp, opt_interp_param, opt_backproj, opt_trace);
371 mpi_t1 = MPI_Wtime();
373 for (ix = 0; ix < im_local->nx; ix++) {
374 void *recvbuf = NULL;
375 if (mpi_ct.my_rank == 0)
376 recvbuf = im_global->v[ix];
378 MPI_Reduce(im_local->v[ix], recvbuf, im_local->ny, MPI_FLOAT, MPI_SUM, 0, mpi_ct.comm);
382 mpi_t2 = MPI_Wtime();
383 mpi_t = mpi_t2 - mpi_t1;
384 MPI_Reduce(&mpi_t, &mpi_t_g, 1, MPI_DOUBLE, MPI_MAX, 0, mpi_ct.comm);
385 if (mpi_ct.my_rank == 0)
386 printf("Time to reduce image = %f secs, max time = %f\n", mpi_t, mpi_t_g);
389 if (mpi_ct.my_rank == 0)
390 time_end = MPI_Wtime();
392 time_end = td_current_sec();
397 if (mpi_ct.my_rank == 0)
400 raysum_close (rs_global);
401 strncpy (im_global->remark, remark, sizeof(im_global->remark));
402 im_global->calctime = time_end - time_start;
403 image_save (im_global);
405 fprintf (stdout, "Time active = %.2f\n", im_global->calctime);
417 static void mpi_scatter_rs (RAYSUM *rs_global, RAYSUM *rs_local, const int opt_debug)
423 if (mpi_ct.my_rank == 0) {
424 for (iproc = 0; iproc < mpi_ct.nproc; iproc++) {
425 end_work_unit = mpi_ct.start_work_unit[iproc] + mpi_ct.local_work_units[iproc] - 1;
427 for (iw = mpi_ct.start_work_unit[iproc]; iw <= end_work_unit; iw++) {
428 MPI_Send(&rs_global->view[iw]->ndet, 1, MPI_INT, iproc, 0, mpi_ct.comm);
429 MPI_Send(&rs_global->view[iw]->view_angle, 1, MPI_DOUBLE, iproc, 0, mpi_ct.comm);
430 MPI_Send(rs_global->view[iw]->detval, rs_global->ndet, MPI_FLOAT, iproc, 0, mpi_ct.comm);
436 fprintf(stdout, "Receiving rs data in process %d\n", mpi_ct.my_rank);
438 end_work_unit = mpi_ct.local_work_units[mpi_ct.my_rank] - 1;
439 for (iw = 0; iw <= end_work_unit; iw++) {
442 MPI_Recv(&rs_local->view[iw]->ndet, 1, MPI_INT, 0, 0, mpi_ct.comm, &status);
443 MPI_Recv(&rs_local->view[iw]->view_angle, 1, MPI_DOUBLE, 0, 0, mpi_ct.comm, &status);
444 MPI_Recv(rs_local->view[iw]->detval, rs_local->ndet, MPI_FLOAT, 0, 0, mpi_ct.comm, &status);
446 rs_local->nview = mpi_ct.local_work_units[mpi_ct.my_rank];
448 MPI_Barrier(MPI_COMM_WORLD);
449 fprintf(stdout, "Done with mpi_scatter_rs in process %2d\n", mpi_ct.my_rank);
455 static void print_raysum_info(const RAYSUM *rs)
457 printf ("Number of detectors: %d\n", rs->ndet);
458 printf (" Number of views: %d\n", rs->nview);
459 printf (" Remark: %s\n", rs->remark);
460 printf (" phmlen: %f\n", rs->phmlen);
461 printf (" det_start: %f\n", rs->det_start);
462 printf (" det_inc: %f\n", rs->det_inc);
463 printf (" rot_start: %f\n", rs->rot_start);
464 printf (" rot_inc: %f\n", rs->rot_inc);
469 main (const int argc, char *const argv[])
471 return (ctrec_main(argc, argv));