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