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