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