r70: added imagefiles
[ctsim.git] / src / phm2rs.c
1 /*****************************************************************************
2 **  This is part of the CTSim program
3 **  Copyright (C) 1983-2000 Kevin Rosenberg
4 **
5 ** PURPOSE
6 **   Generate raysum projections from phantom object
7 **
8 ** HISTORY
9 **   1984 -- Final DOS version  
10 **   1999 -- First UNIX version
11 **
12 **  $Id: phm2rs.c,v 1.12 2000/05/24 22:50:04 kevin Exp $
13 **  $Log: phm2rs.c,v $
14 **  Revision 1.12  2000/05/24 22:50:04  kevin
15 **  Added support for new SGP library
16 **
17 **  Revision 1.11  2000/05/16 04:33:59  kevin
18 **  Improved option processing
19 **
20 **  Revision 1.9  2000/05/08 20:02:32  kevin
21 **  ANSI C changes
22 **
23 **  Revision 1.8  2000/05/04 18:16:34  kevin
24 **  renamed filter definitions
25 **
26 **  Revision 1.7  2000/05/03 08:49:50  kevin
27 **  Code cleanup
28 **
29 **  Revision 1.6  2000/04/30 18:23:53  kevin
30 **  Code cleaning
31 **
32 **  Revision 1.5  2000/04/30 10:13:27  kevin
33 **  Fixed MPI bugs
34 **
35 **  Revision 1.4  2000/04/30 04:06:13  kevin
36 **  Update Raysum i/o routines
37 **  Fix MPI bug in ctrec (scatter_raysum) that referenced rs_global
38 **
39 **  Revision 1.3  2000/04/29 23:24:56  kevin
40 **  *** empty log message ***
41 **
42 **  Revision 1.2  2000/04/28 13:50:45  kevin
43 **  Removed Makefile Makefile.in that are automatically generated by autoconf
44 **
45 **  Revision 1.1.1.1  2000/04/28 13:02:44  kevin
46 **  Initial CVS import for first public release
47 **
48 **
49 **
50 **  This program is free software; you can redistribute it and/or modify
51 **  it under the terms of the GNU General Public License (version 2) as
52 **  published by the Free Software Foundation.
53 **
54 **  This program is distributed in the hope that it will be useful,
55 **  but WITHOUT ANY WARRANTY; without even the implied warranty of
56 **  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
57 **  GNU General Public License for more details.
58 **
59 **  You should have received a copy of the GNU General Public License
60 **  along with this program; if not, write to the Free Software
61 **  Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
62 ******************************************************************************/
63
64
65 #include "ct.h"
66
67 enum { O_PHANTOM, O_DESC, O_NRAY, O_ROTANGLE, O_PHMFILE,
68        O_TRACE, O_VERBOSE, O_HELP, O_DEBUG, O_VERSION };
69
70 static struct option phm2rs_options[] = 
71 {
72   {"phantom", 1, 0, O_PHANTOM},
73   {"phmfile", 1, 0, O_PHMFILE},
74   {"desc", 1, 0, O_DESC},
75   {"nray", 1, 0, O_NRAY},
76   {"rotangle", 1, 0, O_ROTANGLE},
77   {"trace", 1, 0, O_TRACE},
78   {"verbose", 0, 0, O_VERBOSE},
79   {"help", 0, 0, O_HELP},
80   {"debug", 0, 0, O_DEBUG},
81   {"version", 0, 0, O_VERSION},
82   {0, 0, 0, 0}
83 };
84
85
86 void 
87 phm2rs_usage (const char *program)
88 {
89 #ifdef MPI_CT
90 if (mpi_ct.my_rank == 0)
91   {
92 #endif  
93   fprintf(stdout,"usage: %s outfile ndet nview [--phantom phantom-name] [--phmfile filename] [OPTIONS]\n", kbasename(program));
94   fprintf(stdout,"Calculate raysums (projections) through phantom object, either\n");
95   fprintf(stdout,"a predefined --phantom or a --phmfile\n");
96   fprintf(stdout,"\n");
97   fprintf(stdout,"     outfile      Name of output file for raysums\n");
98   fprintf(stdout,"     ndet         Number of detectors\n");
99   fprintf(stdout,"     nview        Number of rotated views\n");
100   fprintf(stdout,"     --phantom    Phantom to use for projection\n");
101   fprintf(stdout,"        herman    Herman head phantom\n");
102   fprintf(stdout,"        rowland   Rowland head phantom\n");
103   fprintf(stdout,"        browland  Bordered Rowland head phantom\n");
104   fprintf(stdout,"        unitpulse Unit pulse phantom\n");
105   fprintf(stdout,"     --phmfile    Get Phantom from phantom file\n");
106   fprintf(stdout,"     --desc       Description of raysum\n");
107   fprintf(stdout,"     --nray       Number of rays per detector (default = 1)\n");
108   fprintf(stdout,"     --rotangle   Degrees to rotate view through, multiple of PI (default = 1)\n");
109   fprintf(stdout,"     --trace      Trace level to use\n");
110   fprintf(stdout,"        none      No tracing (default)\n");
111   fprintf(stdout,"        text      Trace text level\n");
112   fprintf(stdout,"        phm       Trace phantom image\n");
113   fprintf(stdout,"        rays      Trace rays\n");
114   fprintf(stdout,"        plot      Trace plot\n");
115   fprintf(stdout,"        clipping  Trace clipping\n");
116   fprintf(stdout,"     --verbose    Verbose mode\n");
117   fprintf(stdout,"     --debug      Debug mode\n");
118   fprintf(stdout,"     --version    Print version\n");
119   fprintf(stdout,"     --help       Print this help message\n");
120 #ifdef MPI_CT
121   }
122   MPI_Abort(mpi_ct.comm, 1);
123 #endif
124 }
125
126 #ifdef MPI_CT
127 void mpi_gather_rs (RAYSUM *rs_global, RAYSUM *rs_local, const int opt_debug);
128 #endif
129
130 int 
131 phm2rs_main (const int argc, char *const argv[])
132 {
133   DETECTOR *det;
134   PHANTOM *phm = NULL;
135   RAYSUM *rs_global = NULL;
136   char str[256], *opt_outfile = NULL, opt_desc[MAXREMARK+1];
137   char opt_phmfilename[256];
138   char *endptr, *endstr;
139   int opt_ndet, opt_nview;
140   int opt_nray = 1;
141   int opt_trace = 0, opt_phmnum = -1;
142   int opt_verbose = 0;
143   int opt_debug = 0;
144   double opt_rotangle = 1;
145   double time_start = 0, time_end = 0;
146
147 #ifndef MPI_CT
148   time_start = td_current_sec();
149 #else
150   RAYSUM *rs_local;
151   int mpi_argc = argc;
152   char **mpi_argv = (char **) argv;
153   double mpi_t1 = 0, mpi_t2 = 0, mpi_t, mpi_t_g;
154
155   MPI_Init(&mpi_argc, &mpi_argv);
156   MPI_Comm_dup (MPI_COMM_WORLD, &mpi_ct.comm);
157   MPI_Comm_size(mpi_ct.comm, &mpi_ct.nproc);
158   MPI_Comm_rank(mpi_ct.comm, &mpi_ct.my_rank);
159
160   if (mpi_ct.nproc > MPI_MAX_PROCESS) {
161     sys_error(ERR_FATAL, "Number of mpi processes (%d) exceeds max processes (%d)", mpi_ct.nproc, MPI_MAX_PROCESS);
162     exit(1);
163   }
164   time_start = MPI_Wtime();
165 #endif
166   
167 #ifdef MPI_CT
168   if (mpi_ct.my_rank == 0) {
169 #endif
170     strcpy(opt_desc, "");
171     strcpy(opt_phmfilename, "");
172     while (1) {
173       int c = getopt_long(argc, argv, "", phm2rs_options, NULL);
174       char *endptr = NULL;
175       char *endstr;
176       
177       if (c == -1)
178         break;
179       
180       switch (c) {
181       case O_PHANTOM:
182         if ((opt_phmnum = opt_set_phantom (optarg)) < 0) {
183           phm2rs_usage(argv[0]);
184           return (1);
185         }
186         phm = phm_create (opt_phmnum);
187         break;
188       case O_PHMFILE:
189 #ifdef MPI_CT
190         if (mpi_ct.my_rank == 0)
191           fprintf(stderr, "Can not read phantom from file in MPI mode\n");
192         return (1);
193 #endif
194         strncpy(opt_phmfilename, optarg, sizeof(opt_phmfilename));
195         phm = phm_create_from_file(opt_phmfilename);
196         break;
197       case O_VERBOSE:
198         opt_verbose = 1;
199         break;
200       case O_DEBUG:
201         opt_debug = 1;
202         break;
203         break;
204       case O_TRACE:
205         if ((opt_trace = opt_set_trace(optarg)) < 0) {
206           phm2rs_usage(argv[0]);
207           return (1);
208         }
209         break;
210       case O_DESC:
211         strncpy(opt_desc, optarg, sizeof(opt_desc));
212         break;
213       case O_ROTANGLE:
214         opt_rotangle = strtod(optarg, &endptr);
215         endstr = optarg + strlen(optarg);
216         if (endptr != endstr) {
217           fprintf(stderr,"Error setting --rotangle to %s\n", optarg);
218           phm2rs_usage(argv[0]);
219           return (1);
220         }
221         break;
222       case O_NRAY:
223         opt_nray = strtol(optarg, &endptr, 10);
224         endstr = optarg + strlen(optarg);
225         if (endptr != endstr) {
226           fprintf(stderr,"Error setting --nray to %s\n", optarg);
227           phm2rs_usage(argv[0]);
228           return (1);
229         }
230         break;
231         case O_VERSION:
232 #ifdef VERSION
233           fprintf(stdout, "Version %s\n", VERSION);
234 #else
235           fprintf(stderr, "Unknown version number");
236 #endif
237           exit(0);
238       case O_HELP:
239       case '?':
240         phm2rs_usage(argv[0]);
241         return (0);
242       default:
243         phm2rs_usage(argv[0]);
244         return (1);
245       }
246     }
247   
248     if (phm == NULL) {
249       fprintf(stderr, "No phantom defined\n");
250       phm2rs_usage(argv[0]);
251       return (1);
252     }
253     if (optind + 3 != argc) {
254       phm2rs_usage(argv[0]);
255       return (1);
256     }
257
258     opt_outfile = argv[optind];
259     opt_ndet = strtol(argv[optind+1], &endptr, 10);
260     endstr = argv[optind+1] + strlen(argv[optind+1]);
261     if (endptr != endstr) {
262       fprintf(stderr,"Error setting --ndet to %s\n", argv[optind+1]);
263       phm2rs_usage(argv[0]);
264       return (1);
265     }
266     opt_nview = strtol(argv[optind+2], &endptr, 10);
267     endstr = argv[optind+2] + strlen(argv[optind+2]);
268     if (endptr != endstr) {
269       fprintf(stderr,"Error setting --nview to %s\n", argv[optind+2]);
270       phm2rs_usage(argv[0]);
271       return (1);
272     }
273
274     snprintf(str, sizeof(str), 
275              "Raysum_Collect: NDet=%d, Nview=%d, NRay=%d, RotAngle=%.2f, ", 
276              opt_ndet, opt_nview, opt_nray, opt_rotangle);
277     if (opt_phmfilename[0]) {
278       strncat(str, "Phantom=", sizeof(str));
279       strncat(str, opt_phmfilename, sizeof(str));
280     } else if (opt_phmnum != -1) {
281       strncat(str, "Phantom=", sizeof(str));
282       strncat(str, name_of_phantom(opt_phmnum), sizeof(str));
283     }
284     if (opt_desc[0]) {
285       strncat(str, ": ", sizeof(str));
286       strncat(str, opt_desc, sizeof(str));
287     }
288     strncpy(opt_desc, str, sizeof(opt_desc));
289 #ifdef MPI_CT
290   }
291 #endif
292
293 #ifdef MPI_CT
294   MPI_Barrier(mpi_ct.comm);
295   MPI_Bcast(&opt_rotangle, 1, MPI_DOUBLE, 0, mpi_ct.comm);
296   MPI_Bcast(&opt_nview, 1, MPI_INT, 0, mpi_ct.comm);
297   MPI_Bcast(&opt_ndet, 1, MPI_INT, 0, mpi_ct.comm);
298   MPI_Bcast(&opt_nray, 1, MPI_INT, 0, mpi_ct.comm);
299   MPI_Bcast(&opt_phmnum, 1, MPI_INT, 0, mpi_ct.comm);
300   MPI_Bcast(&opt_verbose, 1, MPI_INT, 0, mpi_ct.comm);
301   MPI_Bcast(&opt_debug, 1, MPI_INT, 0, mpi_ct.comm);
302   MPI_Bcast(&opt_trace, 1, MPI_INT, 0, mpi_ct.comm);
303
304   if (mpi_ct.my_rank > 0 && opt_phmnum >= 0)
305     phm = phm_create (opt_phmnum);
306 #endif
307
308   opt_rotangle *= PI;
309   det = detector_create (phm, DETECTOR_PARALLEL, opt_ndet, opt_nview, opt_nray, opt_rotangle);
310
311 #ifdef MPI_CT
312   mpi_ct_calc_work_units(opt_nview);
313   if (mpi_ct.my_rank == 0) {
314     rs_global = raysum_create_from_det (opt_outfile, det);
315     raysum_alloc_views(rs_global);
316   }
317   
318   rs_local = raysum_create_from_det (NULL, det);
319   rs_local->nview = mpi_ct.local_work_units[mpi_ct.my_rank];
320   if (opt_debug)
321     printf("rs_local->nview = %d (process %d)\n", rs_local->nview, mpi_ct.my_rank);
322
323   if (opt_verbose)
324     mpi_t1 = MPI_Wtime();
325   raysum_collect (rs_local, det, phm, mpi_ct.start_work_unit[mpi_ct.my_rank], opt_trace, FALSE);
326   if (opt_verbose) {
327     mpi_t2 = MPI_Wtime();
328     mpi_t = mpi_t2 - mpi_t1;
329     MPI_Reduce(&mpi_t, &mpi_t_g, 1, MPI_DOUBLE, MPI_MAX, 0, mpi_ct.comm);
330     if (mpi_ct.my_rank == 0)
331       printf("Time to collect rs = %f secs, Max = %f secs\n", mpi_t, mpi_t_g);
332   }
333
334   if (opt_verbose)
335     mpi_t1 = MPI_Wtime();
336   mpi_gather_rs(rs_global, rs_local, opt_debug);
337   if (opt_verbose) {
338     mpi_t2 = MPI_Wtime();
339     mpi_t = mpi_t2 - mpi_t1;
340     MPI_Reduce(&mpi_t, &mpi_t_g, 1, MPI_DOUBLE, MPI_MAX, 0, mpi_ct.comm);
341     if (mpi_ct.my_rank == 0)
342       printf("Time to gather rs = %f secs, Max = %f secs\n", mpi_t, mpi_t_g);
343   }
344 #else
345   rs_global = raysum_create_from_det (opt_outfile, det);
346   raysum_collect (rs_global, det, phm, 0, opt_trace, FALSE);
347 #endif
348   
349 #ifndef MPI_CT
350   time_end = td_current_sec();
351 #else
352   time_end = MPI_Wtime();
353   if (mpi_ct.my_rank == 0) 
354 #endif
355     {
356       rs_global->calctime = time_end - time_start;
357       strncpy (rs_global->remark, opt_desc, sizeof(rs_global->remark));
358       if (opt_verbose) {
359         fprintf(stdout, "Remark: %s\n", rs_global->remark);
360         fprintf(stdout, "Time active = %.2f secs\n", rs_global->calctime);
361       }
362     }
363
364 #ifdef MPI_CT
365   if (mpi_ct.my_rank == 0) 
366 #endif
367     {
368       raysum_write (rs_global);
369       raysum_close (rs_global);
370     }
371
372   phm_free (phm);
373   detector_free (det);
374
375   return (0);
376 }
377
378
379 /* FUNCTION
380  *    mpi_gather_rs
381  *
382  * SYNOPSIS
383  *    Gather's raysums from all processes in rs_local in rs_global in process 0
384  */
385
386 #ifdef MPI_CT
387 void mpi_gather_rs (RAYSUM *rs_global, RAYSUM *rs_local, const int opt_debug)
388 {
389   int iproc;
390   int end_work_unit;
391   int iw;
392
393   end_work_unit = mpi_ct.local_work_units[mpi_ct.my_rank] - 1;
394   for (iw = 0; iw <= end_work_unit; iw++) {
395     MPI_Send(&rs_local->view[iw]->view_angle, 1, MPI_DOUBLE, 0, 0, mpi_ct.comm);
396     MPI_Send(&rs_local->view[iw]->ndet, 1, MPI_INT, 0, 0, mpi_ct.comm);
397     MPI_Send(rs_local->view[iw]->detval, rs_local->ndet, MPI_FLOAT, 0, 0, mpi_ct.comm);
398   }
399
400   if (mpi_ct.my_rank == 0) {
401     for (iproc = 0; iproc < mpi_ct.nproc; iproc++) {
402       end_work_unit = mpi_ct.start_work_unit[iproc] 
403         + mpi_ct.local_work_units[iproc] 
404         - 1;
405
406       if (opt_debug)
407         fprintf(stdout, "Recv rs data from process %d\n", iproc);
408
409       for (iw = mpi_ct.start_work_unit[iproc]; iw <= end_work_unit; iw++) {
410         MPI_Status status;
411
412         MPI_Recv(&rs_global->view[iw]->view_angle, 1, MPI_DOUBLE, iproc, 0, mpi_ct.comm, &status);
413         MPI_Recv(&rs_global->view[iw]->ndet, 1, MPI_INT, iproc, 0, mpi_ct.comm, &status);
414         MPI_Recv(rs_global->view[iw]->detval, rs_global->ndet, MPI_FLOAT, iproc, 0, mpi_ct.comm, &status);
415       }
416     }
417   }
418
419 }
420 #endif
421
422
423 #ifndef NO_MAIN
424 int 
425 main (const int argc, char *const argv[])
426 {
427   return (phm2rs_main(argc, argv));
428 }
429 #endif
430