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