1 /*****************************************************************************
2 ** This is part of the CTSim program
3 ** Copyright (C) 1983-2000 Kevin Rosenberg
6 ** Generate raysum projections from phantom object
9 ** 1984 -- Final DOS version
10 ** 1999 -- First UNIX version
12 ** $Id: phm2rs.cpp,v 1.2 2000/06/07 10:12:05 kevin Exp $
13 ** $Log: phm2rs.cpp,v $
14 ** Revision 1.2 2000/06/07 10:12:05 kevin
15 ** Upgraded from MPI to MPI++
17 ** Revision 1.1 2000/06/07 02:29:05 kevin
18 ** Initial C++ versions
20 ** Revision 1.12 2000/05/24 22:50:04 kevin
21 ** Added support for new SGP library
23 ** Revision 1.11 2000/05/16 04:33:59 kevin
24 ** Improved option processing
26 ** Revision 1.9 2000/05/08 20:02:32 kevin
29 ** Revision 1.8 2000/05/04 18:16:34 kevin
30 ** renamed filter definitions
32 ** Revision 1.7 2000/05/03 08:49:50 kevin
35 ** Revision 1.6 2000/04/30 18:23:53 kevin
38 ** Revision 1.5 2000/04/30 10:13:27 kevin
41 ** Revision 1.4 2000/04/30 04:06:13 kevin
42 ** Update Raysum i/o routines
43 ** Fix MPI bug in ctrec (scatter_raysum) that referenced rs_global
45 ** Revision 1.3 2000/04/29 23:24:56 kevin
46 ** *** empty log message ***
48 ** Revision 1.2 2000/04/28 13:50:45 kevin
49 ** Removed Makefile Makefile.in that are automatically generated by autoconf
51 ** Revision 1.1.1.1 2000/04/28 13:02:44 kevin
52 ** Initial CVS import for first public release
56 ** This program is free software; you can redistribute it and/or modify
57 ** it under the terms of the GNU General Public License (version 2) as
58 ** published by the Free Software Foundation.
60 ** This program is distributed in the hope that it will be useful,
61 ** but WITHOUT ANY WARRANTY; without even the implied warranty of
62 ** MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
63 ** GNU General Public License for more details.
65 ** You should have received a copy of the GNU General Public License
66 ** along with this program; if not, write to the Free Software
67 ** Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
68 ******************************************************************************/
73 enum { O_PHANTOM, O_DESC, O_NRAY, O_ROTANGLE, O_PHMFILE,
74 O_TRACE, O_VERBOSE, O_HELP, O_DEBUG, O_VERSION };
76 static struct option phm2rs_options[] =
78 {"phantom", 1, 0, O_PHANTOM},
79 {"phmfile", 1, 0, O_PHMFILE},
80 {"desc", 1, 0, O_DESC},
81 {"nray", 1, 0, O_NRAY},
82 {"rotangle", 1, 0, O_ROTANGLE},
83 {"trace", 1, 0, O_TRACE},
84 {"verbose", 0, 0, O_VERBOSE},
85 {"help", 0, 0, O_HELP},
86 {"debug", 0, 0, O_DEBUG},
87 {"version", 0, 0, O_VERSION},
93 phm2rs_usage (const char *program)
96 if (mpi_ct.my_rank == 0)
99 fprintf(stdout,"usage: %s outfile ndet nview [--phantom phantom-name] [--phmfile filename] [OPTIONS]\n", kbasename(program));
100 fprintf(stdout,"Calculate raysums (projections) through phantom object, either\n");
101 fprintf(stdout,"a predefined --phantom or a --phmfile\n");
102 fprintf(stdout,"\n");
103 fprintf(stdout," outfile Name of output file for raysums\n");
104 fprintf(stdout," ndet Number of detectors\n");
105 fprintf(stdout," nview Number of rotated views\n");
106 fprintf(stdout," --phantom Phantom to use for projection\n");
107 fprintf(stdout," herman Herman head phantom\n");
108 fprintf(stdout," rowland Rowland head phantom\n");
109 fprintf(stdout," browland Bordered Rowland head phantom\n");
110 fprintf(stdout," unitpulse Unit pulse phantom\n");
111 fprintf(stdout," --phmfile Get Phantom from phantom file\n");
112 fprintf(stdout," --desc Description of raysum\n");
113 fprintf(stdout," --nray Number of rays per detector (default = 1)\n");
114 fprintf(stdout," --rotangle Degrees to rotate view through, multiple of PI (default = 1)\n");
115 fprintf(stdout," --trace Trace level to use\n");
116 fprintf(stdout," none No tracing (default)\n");
117 fprintf(stdout," text Trace text level\n");
118 fprintf(stdout," phm Trace phantom image\n");
119 fprintf(stdout," rays Trace rays\n");
120 fprintf(stdout," plot Trace plot\n");
121 fprintf(stdout," clipping Trace clipping\n");
122 fprintf(stdout," --verbose Verbose mode\n");
123 fprintf(stdout," --debug Debug mode\n");
124 fprintf(stdout," --version Print version\n");
125 fprintf(stdout," --help Print this help message\n");
128 mpi_ct.comm.Abort(1);
133 void mpi_gather_rs (RAYSUM *rs_global, RAYSUM *rs_local, const int opt_debug);
137 phm2rs_main (const int argc, char *const argv[])
141 RAYSUM *rs_global = NULL;
142 char str[256], *opt_outfile = NULL, opt_desc[MAXREMARK+1];
143 char opt_phmfilename[256];
144 char *endptr, *endstr;
145 int opt_ndet, opt_nview;
147 int opt_trace = 0, opt_phmnum = -1;
150 double opt_rotangle = 1;
151 double time_start = 0, time_end = 0;
154 time_start = td_current_sec();
158 char **mpi_argv = (char **) argv;
159 double mpi_t1 = 0, mpi_t2 = 0, mpi_t, mpi_t_g;
161 MPI::Init (mpi_argc, mpi_argv);
162 mpi_ct.comm = MPI::COMM_WORLD.Dup ();
163 mpi_ct.nproc = mpi_ct.comm.Get_size();
164 mpi_ct.my_rank = mpi_ct.comm.Get_rank();
166 if (mpi_ct.nproc > CT_MPI_MAX_PROCESS) {
167 sys_error(ERR_FATAL, "Number of mpi processes (%d) exceeds max processes (%d)",
168 mpi_ct.nproc, CT_MPI_MAX_PROCESS);
171 time_start = MPI::Wtime();
175 if (mpi_ct.my_rank == 0) {
177 strcpy(opt_desc, "");
178 strcpy(opt_phmfilename, "");
180 int c = getopt_long(argc, argv, "", phm2rs_options, NULL);
189 if ((opt_phmnum = opt_set_phantom (optarg)) < 0) {
190 phm2rs_usage(argv[0]);
193 phm = phm_create (opt_phmnum);
197 if (mpi_ct.my_rank == 0)
198 fprintf(stderr, "Can not read phantom from file in MPI mode\n");
201 strncpy(opt_phmfilename, optarg, sizeof(opt_phmfilename));
202 phm = phm_create_from_file(opt_phmfilename);
212 if ((opt_trace = opt_set_trace(optarg)) < 0) {
213 phm2rs_usage(argv[0]);
218 strncpy(opt_desc, optarg, sizeof(opt_desc));
221 opt_rotangle = strtod(optarg, &endptr);
222 endstr = optarg + strlen(optarg);
223 if (endptr != endstr) {
224 fprintf(stderr,"Error setting --rotangle to %s\n", optarg);
225 phm2rs_usage(argv[0]);
230 opt_nray = strtol(optarg, &endptr, 10);
231 endstr = optarg + strlen(optarg);
232 if (endptr != endstr) {
233 fprintf(stderr,"Error setting --nray to %s\n", optarg);
234 phm2rs_usage(argv[0]);
240 fprintf(stdout, "Version %s\n", VERSION);
242 fprintf(stderr, "Unknown version number");
247 phm2rs_usage(argv[0]);
250 phm2rs_usage(argv[0]);
256 fprintf(stderr, "No phantom defined\n");
257 phm2rs_usage(argv[0]);
260 if (optind + 3 != argc) {
261 phm2rs_usage(argv[0]);
265 opt_outfile = argv[optind];
266 opt_ndet = strtol(argv[optind+1], &endptr, 10);
267 endstr = argv[optind+1] + strlen(argv[optind+1]);
268 if (endptr != endstr) {
269 fprintf(stderr,"Error setting --ndet to %s\n", argv[optind+1]);
270 phm2rs_usage(argv[0]);
273 opt_nview = strtol(argv[optind+2], &endptr, 10);
274 endstr = argv[optind+2] + strlen(argv[optind+2]);
275 if (endptr != endstr) {
276 fprintf(stderr,"Error setting --nview to %s\n", argv[optind+2]);
277 phm2rs_usage(argv[0]);
281 snprintf(str, sizeof(str),
282 "Raysum_Collect: NDet=%d, Nview=%d, NRay=%d, RotAngle=%.2f, ",
283 opt_ndet, opt_nview, opt_nray, opt_rotangle);
284 if (opt_phmfilename[0]) {
285 strncat(str, "Phantom=", sizeof(str));
286 strncat(str, opt_phmfilename, sizeof(str));
287 } else if (opt_phmnum != -1) {
288 strncat(str, "Phantom=", sizeof(str));
289 strncat(str, name_of_phantom(opt_phmnum), sizeof(str));
292 strncat(str, ": ", sizeof(str));
293 strncat(str, opt_desc, sizeof(str));
295 strncpy(opt_desc, str, sizeof(opt_desc));
301 mpi_ct.comm.Barrier ();
302 mpi_ct.comm.Bcast (&opt_rotangle, 1, MPI::DOUBLE, 0);
303 mpi_ct.comm.Bcast (&opt_nview, 1, MPI::INT, 0);
304 mpi_ct.comm.Bcast (&opt_ndet, 1, MPI::INT, 0);
305 mpi_ct.comm.Bcast (&opt_nray, 1, MPI::INT, 0);
306 mpi_ct.comm.Bcast (&opt_phmnum, 1, MPI::INT, 0);
307 mpi_ct.comm.Bcast (&opt_verbose, 1, MPI::INT, 0);
308 mpi_ct.comm.Bcast (&opt_debug, 1, MPI::INT, 0);
309 mpi_ct.comm.Bcast (&opt_trace, 1, MPI::INT, 0);
311 if (mpi_ct.my_rank > 0 && opt_phmnum >= 0)
312 phm = phm_create (opt_phmnum);
316 det = detector_create (phm, DETECTOR_PARALLEL, opt_ndet, opt_nview, opt_nray, opt_rotangle);
319 mpi_ct_calc_work_units(opt_nview);
320 if (mpi_ct.my_rank == 0) {
321 rs_global = raysum_create_from_det (opt_outfile, det);
322 raysum_alloc_views(rs_global);
325 rs_local = raysum_create_from_det (NULL, det);
326 rs_local->nview = mpi_ct.local_work_units[mpi_ct.my_rank];
328 printf("rs_local->nview = %d (process %d)\n", rs_local->nview, mpi_ct.my_rank);
331 mpi_t1 = MPI::Wtime();
332 raysum_collect (rs_local, det, phm, mpi_ct.start_work_unit[mpi_ct.my_rank], opt_trace, FALSE);
334 mpi_t2 = MPI::Wtime();
335 mpi_t = mpi_t2 - mpi_t1;
336 mpi_ct.comm.Reduce(&mpi_t, &mpi_t_g, 1, MPI::DOUBLE, MPI::MAX, 0);
337 if (mpi_ct.my_rank == 0)
338 printf("Time to collect rs = %f secs, Max = %f secs\n", mpi_t, mpi_t_g);
342 mpi_t1 = MPI::Wtime();
343 mpi_gather_rs(rs_global, rs_local, opt_debug);
345 mpi_t2 = MPI::Wtime();
346 mpi_t = mpi_t2 - mpi_t1;
347 mpi_ct.comm.Reduce(&mpi_t, &mpi_t_g, 1, MPI::DOUBLE, MPI::MAX, 0);
348 if (mpi_ct.my_rank == 0)
349 printf("Time to gather rs = %f secs, Max = %f secs\n", mpi_t, mpi_t_g);
352 rs_global = raysum_create_from_det (opt_outfile, det);
353 raysum_collect (rs_global, det, phm, 0, opt_trace, FALSE);
357 time_end = td_current_sec();
359 time_end = MPI::Wtime();
360 if (mpi_ct.my_rank == 0)
363 rs_global->calctime = time_end - time_start;
364 strncpy (rs_global->remark, opt_desc, sizeof(rs_global->remark));
366 fprintf(stdout, "Remark: %s\n", rs_global->remark);
367 fprintf(stdout, "Time active = %.2f secs\n", rs_global->calctime);
372 if (mpi_ct.my_rank == 0)
375 raysum_write (rs_global);
376 raysum_close (rs_global);
390 * Gather's raysums from all processes in rs_local in rs_global in process 0
394 void mpi_gather_rs (RAYSUM *rs_global, RAYSUM *rs_local, const int opt_debug)
400 end_work_unit = mpi_ct.local_work_units[mpi_ct.my_rank] - 1;
401 for (iw = 0; iw <= end_work_unit; iw++) {
402 mpi_ct.comm.Send(&rs_local->view[iw]->view_angle, 1, MPI::DOUBLE, 0, 0);
403 mpi_ct.comm.Send(&rs_local->view[iw]->ndet, 1, MPI::INT, 0, 0);
404 mpi_ct.comm.Send(rs_local->view[iw]->detval, rs_local->ndet, MPI::FLOAT, 0, 0);
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, "Recv rs data from process %d\n", iproc);
416 for (iw = mpi_ct.start_work_unit[iproc]; iw <= end_work_unit; iw++) {
419 mpi_ct.comm.Recv(&rs_global->view[iw]->view_angle, 1, MPI::DOUBLE, iproc, 0, status);
420 mpi_ct.comm.Recv(&rs_global->view[iw]->ndet, 1, MPI::INT, iproc, 0, status);
421 mpi_ct.comm.Recv(rs_global->view[iw]->detval, rs_global->ndet, MPI::FLOAT, iproc, 0, status);
432 main (const int argc, char *const argv[])
434 return (phm2rs_main(argc, argv));