r31: Code cleanup
[ctsim.git] / src / phm2sdf.c
1 /*****************************************************************************
2 **  This is part of the CTSim program
3 **  Copyright (C) 1983-2000 Kevin Rosenberg
4 **
5 **  $Id: phm2sdf.c,v 1.6 2000/05/03 08:49:50 kevin Exp $
6 **  $Log: phm2sdf.c,v $
7 **  Revision 1.6  2000/05/03 08:49:50  kevin
8 **  Code cleanup
9 **
10 **  Revision 1.5  2000/04/30 19:17:54  kevin
11 **  *** empty log message ***
12 **
13 **  Revision 1.4  2000/04/30 04:06:13  kevin
14 **  Update Raysum i/o routines
15 **  Fix MPI bug in ctrec (scatter_raysum) that referenced rs_global
16 **
17 **  Revision 1.3  2000/04/28 18:19:01  kevin
18 **  removed unused files
19 **
20 **  Revision 1.2  2000/04/28 13:50:45  kevin
21 **  Removed Makefile Makefile.in that are automatically generated by autoconf
22 **
23 **  Revision 1.1.1.1  2000/04/28 13:02:44  kevin
24 **  Initial CVS import for first public release
25 **
26 **
27 **
28 **  This program is free software; you can redistribute it and/or modify
29 **  it under the terms of the GNU General Public License (version 2) as
30 **  published by the Free Software Foundation.
31 **
32 **  This program is distributed in the hope that it will be useful,
33 **  but WITHOUT ANY WARRANTY; without even the implied warranty of
34 **  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
35 **  GNU General Public License for more details.
36 **
37 **  You should have received a copy of the GNU General Public License
38 **  along with this program; if not, write to the Free Software
39 **  Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
40 ******************************************************************************/
41 /* FILE
42  *   phm2sdf.c                  Generate a SDF image from a phantom
43  *
44  * HISTORY
45  *   1984 - Final DOS verion
46  *   1999 - First UNIX version
47  */
48
49 #include "ct.h"
50
51 #define O_PHANTOM       1
52 #define O_DESC          2
53 #define O_NSAMPLE       3
54 #define O_FILTER        4
55 #define O_TRACE         5
56 #define O_VERBOSE       6
57 #define O_HELP          7
58 #define O_PHMFILE       8
59 #define O_FILTER_DOMAIN 9
60 #define O_FILTER_BW     10
61 #define O_FILTER_PARAM  11
62 #define O_DEBUG         12
63 #define O_VERSION       13
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   {"nsample", 1, 0, O_NSAMPLE},
71   {"filter", 1, 0, O_FILTER},
72   {"filter-domain", 1, 0, O_FILTER_DOMAIN},
73   {"filter-bw", 1, 0, O_FILTER_BW},
74   {"filter-param", 1, 0, O_FILTER_PARAM},
75   {"trace", 1, 0, O_TRACE},
76   {"verbose", 0, 0, O_VERBOSE},
77   {"debug", 0, 0, O_DEBUG},
78   {"help", 0, 0, O_HELP},
79   {"version", 0, 0, O_VERSION},
80   {0, 0, 0, 0}
81 };
82
83 void 
84 usage (const char *program)
85 {
86   fprintf(stdout,"usage: %s outfile nx ny [--phantom phantom-name] [--phmfile filename] [--filter filter-name] [OPTIONS]\n",kbasename(program)); 
87   fprintf(stdout,"Generate phantom image from a predefined --phantom or a --phmfile\n");
88   fprintf(stdout,"\n");
89   fprintf(stdout,"     outfile         Name of output file for image\n");
90   fprintf(stdout,"     nx              Number of pixels X-axis\n");
91   fprintf(stdout,"     ny              Number of pixels Y-axis\n");
92   fprintf(stdout,"     --phantom       Phantom to use for projection\n");
93   fprintf(stdout,"        herman       Herman head phantom\n");
94   fprintf(stdout,"        rowland      Rowland head phantom\n");
95   fprintf(stdout,"        browland     Bordered Rowland head phantom\n");
96   fprintf(stdout,"        unitpulse    Unit pulse phantom\n");
97   fprintf(stdout,"     --phmfile       Generate Phantom from phantom file\n");
98   fprintf(stdout,"     --filter        Generate Phantom from a filter function\n");
99   fprintf(stdout,"       abs_bandlimit Abs * Bandlimiting\n");
100   fprintf(stdout,"       abs_sinc      Abs * Sinc\n");
101   fprintf(stdout,"       abs_cos       Abs * Cosine\n");
102   fprintf(stdout,"       abs_hamming   Abs * Hamming\n");
103   fprintf(stdout,"       shepp         Shepp-Logan\n");
104   fprintf(stdout,"       bandlimit     Bandlimiting\n");
105   fprintf(stdout,"       sinc          Sinc\n");
106   fprintf(stdout,"       cos           Cosine\n");
107   fprintf(stdout,"       triangle      Triangle\n");
108   fprintf(stdout,"       hamming       Hamming\n");
109   fprintf(stdout,"     --filter-param  Alpha level for Hamming filter\n");
110   fprintf(stdout,"     --filter-domain Set domain of filter\n");
111   fprintf(stdout,"         spatial     Spatial domain (default)\n");
112   fprintf(stdout,"         freq        Frequency domain\n");
113   fprintf(stdout,"     --filter-bw     Filter bandwidth (default = 1)\n");
114   fprintf(stdout,"     --desc          Description of raysum\n");
115   fprintf(stdout,"     --nsample       Number of samples per axis per pixel (default = 1)\n");
116   fprintf(stdout,"     --trace         Trace level to use\n");
117   fprintf(stdout,"        none         No tracing (default)\n");
118   fprintf(stdout,"        text         Trace text level\n");
119   fprintf(stdout,"        phm          Trace phantom\n");
120   fprintf(stdout,"        rays         Trace rays\n");
121   fprintf(stdout,"        plot         Trace plot\n");
122   fprintf(stdout,"        clipping     Trace clipping\n");
123   fprintf(stdout,"     --debug         Debug mode\n");
124   fprintf(stdout,"     --verbose       Verbose mode\n");
125   fprintf(stdout,"     --version       Print version\n");
126   fprintf(stdout,"     --help          Print this help message\n");
127   exit(1);
128 }
129
130 #ifdef MPI_CT
131 void mpi_gather_image (IMAGE *im_global, IMAGE *im_local, const int opt_debug);
132 #endif
133
134 int 
135 main (const int argc, char *const argv[])
136 {
137   IMAGE *im_global = NULL;
138   PHANTOM *phm = NULL;
139   int opt_nx = 0, opt_ny = 0;
140   int opt_nsample = 1;
141   int opt_phmnum = -1;
142   int opt_filter = -1;
143   int opt_filter_domain = D_SPATIAL;
144   char *opt_outfile = NULL;
145   int opt_debug = 0;
146   char str[256];
147   char opt_desc[256], opt_phmfilename[256];
148   char *endstr, *endptr;
149   double opt_filter_param = 1;
150   double opt_filter_bw = 1.;
151   int opt_trace = TRACE_NONE;
152   int opt_verbose = 0;
153   double time_start=0, time_end=0;
154 #ifdef MPI_CT
155   IMAGE *im_local = NULL;
156   int mpi_argc = argc;
157   char **mpi_argv = (char **) argv;
158   double mpi_t1, mpi_t2, mpi_t, mpi_t_g;
159
160   MPI_Init(&mpi_argc, &mpi_argv);
161   MPI_Comm_dup (MPI_COMM_WORLD, &mpi_ct.comm);
162   MPI_Comm_size(mpi_ct.comm, &mpi_ct.nproc);
163   MPI_Comm_rank(mpi_ct.comm, &mpi_ct.my_rank);
164
165   if (mpi_ct.nproc > MPI_MAX_PROCESS) {
166     sys_error(ERR_FATAL, "Number of mpi processes (%d) exceeds max processes (%d)",
167               mpi_ct.nproc, MPI_MAX_PROCESS);
168     exit(1);
169   }
170 #endif
171
172 #ifdef MPI_CT
173   time_start = MPI_Wtime();
174 #else
175   time_start = td_current_sec();
176 #endif
177   
178 #ifdef MPI_CT
179   if (mpi_ct.my_rank == 0) {
180 #endif
181
182     strcpy(opt_desc, "");
183     strcpy(opt_phmfilename, "");
184     while (1) {
185       int c = getopt_long(argc, argv, "", my_options, NULL);
186       char *endptr = NULL;
187       char *endstr;
188       
189       if (c == -1)
190         break;
191       
192       switch (c) {
193       case O_PHANTOM:
194         opt_phmnum = opt_set_phantom(optarg, argv[0]);
195         break;
196       case O_PHMFILE:
197         strncpy(opt_phmfilename, optarg, sizeof(opt_phmfilename));
198         phm = phm_create_from_file(opt_phmfilename);
199 #ifdef MPI_CT
200         if (mpi_ct.my_rank == 0) 
201           fprintf(stderr, "Can't use phantom from file in MPI mode\n");
202         exit(1);
203 #endif
204         break;
205       case O_VERBOSE:
206         opt_verbose = 1;
207         break;
208       case O_DEBUG:
209         opt_debug = 1;
210         break;
211       case O_TRACE:
212         opt_trace = opt_set_trace(optarg, argv[0]);
213         break;
214       case O_FILTER:
215         opt_filter = opt_set_filter(optarg, argv[0]);
216         break;
217       case O_FILTER_DOMAIN:
218         opt_filter_domain = opt_set_filter_domain(optarg, argv[0]);
219         break;
220       case O_DESC:
221         strncpy(opt_desc, optarg, sizeof(opt_desc));
222         break;
223       case O_FILTER_BW:
224         opt_filter_bw = strtod(optarg, &endptr);
225         endstr = optarg + strlen(optarg);
226         if (endptr != endstr) {
227           fprintf(stderr,"Error setting --filter-bw to %s\n", optarg);
228           usage(argv[0]);
229           exit(1);
230         }
231         break;
232       case O_FILTER_PARAM:
233         opt_filter_param = strtod(optarg, &endptr);
234         endstr = optarg + strlen(optarg);
235         if (endptr != endstr) {
236           fprintf(stderr,"Error setting --filter-param to %s\n", optarg);
237           usage(argv[0]);
238           exit(1);
239         }
240         break;
241       case O_NSAMPLE:
242         opt_nsample = strtol(optarg, &endptr, 10);
243         endstr = optarg + strlen(optarg);
244         if (endptr != endstr) {
245           fprintf(stderr,"Error setting --nsample to %s\n", optarg);
246           usage(argv[0]);
247           exit(1);
248         }
249         break;
250         case O_VERSION:
251 #ifdef VERSION
252           fprintf(stdout, "Version %s\n", VERSION);
253 #else
254           fprintf(stderr, "Unknown version number");
255 #endif
256       case O_HELP:
257       case '?':
258         usage(argv[0]);
259         exit(0);
260       default:
261         usage(argv[0]);
262         exit(1);
263       }
264     }
265     
266     if (phm == NULL && opt_phmnum == -1 && opt_filter == -1) {
267       fprintf(stderr, "No phantom defined\n");
268       usage(argv[0]);
269       exit(1);
270     }
271
272 #if HAVE_INTERACTIVE_GRAPHICS
273     if (opt_trace >= TRACE_PHM)
274       show_phm(phm);
275 #endif
276
277     if (optind + 3 != argc) {
278       usage(argv[0]);
279       exit(1);
280     }
281     opt_outfile = argv[optind];
282     opt_nx = strtol(argv[optind+1], &endptr, 10);
283     endstr = argv[optind+1] + strlen(argv[optind+1]);
284     if (endptr != endstr) {
285       fprintf(stderr,"Error setting nx to %s\n", argv[optind+1]);
286       usage(argv[0]);
287       exit(1);
288     }
289     opt_ny = strtol(argv[optind+2], &endptr, 10);
290     endstr = argv[optind+2] + strlen(argv[optind+2]);
291     if (endptr != endstr) {
292       fprintf(stderr,"Error setting ny to %s\n", argv[optind+2]);
293       usage(argv[0]);
294       exit(1);
295     }
296     
297     snprintf(str, sizeof(str), "nx=%d, ny=%d, nsample=%d, ", opt_nx, opt_ny, opt_nsample);
298     if (opt_phmfilename[0]) {
299       strncat(str, "phantom=", sizeof(str));
300       strncat(str, opt_phmfilename, sizeof(str));
301     }
302     else if (opt_phmnum != -1) {
303       strncat(str, "phantom=", sizeof(str));
304       strncat(str, name_of_phantom(opt_phmnum), sizeof(str));
305     }
306     else if (opt_filter != -1) {
307       strncat(str, "filter=", sizeof(str));
308       strncat(str, name_of_filter(opt_filter), sizeof(str));
309       strncat(str, " - ", sizeof(str));
310       strncat(str, name_of_filter_domain(opt_filter_domain), sizeof(str));
311     }
312     if (opt_desc[0]) {
313       strncat(str, ": ", sizeof(str));
314       strncat(str, opt_desc, sizeof(str));
315     }
316     strncpy(opt_desc, str, sizeof(opt_desc));
317     
318     if (opt_verbose)
319       printf("Compile Phantom to Image\n\n");
320 #ifdef MPI_CT
321   }
322 #endif
323
324 #ifdef MPI_CT
325   mpi_t1 = MPI_Wtime();
326   MPI_Bcast(&opt_verbose, 1, MPI_INT, 0, mpi_ct.comm);
327   MPI_Bcast(&opt_debug, 1, MPI_INT, 0, mpi_ct.comm);
328   MPI_Bcast(&opt_trace, 1, MPI_INT, 0, mpi_ct.comm);
329   MPI_Bcast(&opt_nx, 1, MPI_INT, 0, mpi_ct.comm);
330   MPI_Bcast(&opt_ny, 1, MPI_INT, 0, mpi_ct.comm);
331   MPI_Bcast(&opt_nsample, 1, MPI_INT, 0, mpi_ct.comm);
332   MPI_Bcast(&opt_phmnum, 1, MPI_INT, 0, mpi_ct.comm);
333   MPI_Bcast(&opt_filter, 1, MPI_INT, 0, mpi_ct.comm);
334   MPI_Bcast(&opt_filter_domain, 1, MPI_INT, 0, mpi_ct.comm);
335   MPI_Bcast(&opt_filter_param, 1, MPI_DOUBLE, 0, mpi_ct.comm);
336   MPI_Bcast(&opt_filter_bw, 1, MPI_DOUBLE, 0, mpi_ct.comm);
337
338   if (opt_verbose) {
339     mpi_t2 = MPI_Wtime();
340     mpi_t = mpi_t2 - mpi_t1;
341     MPI_Reduce(&mpi_t, &mpi_t_g, 1, MPI_DOUBLE, MPI_MAX, 0, mpi_ct.comm);
342     if (mpi_ct.my_rank == 0)
343       printf("Time to Bcast vars = %f secs, Max time = %f\n", mpi_t, mpi_t_g);
344   }
345
346   mpi_ct_calc_work_units(opt_nx);
347
348   if (mpi_ct.my_rank == 0)
349     im_global = image_create (opt_outfile, opt_nx, opt_ny);
350
351   im_local = image_create(NULL, opt_nx, opt_ny);
352 #else
353   im_global = image_create (opt_outfile, opt_nx, opt_ny);
354 #endif
355   
356   if (opt_phmnum > 0)
357     phm = phm_create (opt_phmnum);
358 #ifdef MPI_CT
359   else {
360     if (mpi_ct.my_rank == 0)
361       fprintf(stderr, "phmnum < 0\n");
362     exit(1);
363   }
364 #endif
365
366 #ifdef MPI_CT
367   if (phm->type == P_UNIT_PULSE) {
368     if (mpi_ct.my_rank == 0) {
369       im_global->v[opt_nx/2][opt_ny/2] = 1.;
370       im_global->calctime = 0;
371     }
372   } else if (opt_filter != -1) {
373     if (mpi_ct.my_rank == 0) {
374       image_filter_response (im_global, opt_filter_domain, opt_filter_bw, opt_filter, opt_filter_param, opt_trace);
375       im_global->calctime = 0;
376     }
377   } else {
378     if (opt_verbose)
379       mpi_t1 = MPI_Wtime();
380     phm_to_image (phm, im_local, mpi_ct.start_work_unit[mpi_ct.my_rank],
381                   mpi_ct.local_work_units[mpi_ct.my_rank], opt_nsample, opt_trace);
382     if (opt_verbose) {
383       mpi_t2 = MPI_Wtime();
384       mpi_t = mpi_t2 - mpi_t1;
385       MPI_Reduce(&mpi_t, &mpi_t_g, 1, MPI_DOUBLE, MPI_MAX, 0, mpi_ct.comm);
386       if (mpi_ct.my_rank == 0)
387         printf("Time to compile phm = %f secs, Max time = %f\n", mpi_t, mpi_t_g);
388     }
389     mpi_gather_image(im_global, im_local, opt_debug);
390   }
391 #else
392   if (phm->type == P_UNIT_PULSE) {
393     im_global->v[opt_nx/2][opt_ny/2] = 1.;
394     im_global->calctime = 0;
395   } else if (opt_filter != -1) {
396     image_filter_response (im_global, opt_filter_domain, opt_filter_bw, opt_filter, opt_filter_param, opt_trace);
397     im_global->calctime = 0;
398   } else {
399     phm_to_image (phm, im_global, 0, opt_nx, opt_nsample, opt_trace);
400   }
401 #endif
402
403 #ifdef MPI_CT
404   time_end = MPI_Wtime();
405   if (mpi_ct.my_rank == 0) {
406     im_global->calctime = time_end - time_start;
407     strncpy (im_global->remark, opt_desc, sizeof(im_global->remark));
408   }
409
410   if (mpi_ct.my_rank == 0) {
411     image_save (im_global);
412     if (opt_verbose)
413       fprintf (stdout, "Time to compile image = %.2f sec\n\n", im_global->calctime);
414   }
415 #else
416   time_end = td_current_sec();
417   im_global->calctime = time_end - time_start;
418   strncpy (im_global->remark, opt_desc, sizeof(im_global->remark));
419   image_save (im_global);
420
421   if (opt_trace >= TRACE_PHM)
422     {
423       crt_set_cpos (1, 1);
424     }
425   
426   if (opt_verbose)
427     fprintf (stdout, "Time to compile image = %.2f sec\n\n", im_global->calctime);
428 #endif
429
430   if (opt_trace >= TRACE_PHM)
431     {
432       double dmin, dmax;
433       int nx, ny;
434
435       printf ("Enter size of cell (nx, ny): ");
436       scanf ("%d %d", &nx, &ny);
437       printf ("Enter minimum and maximum densities (min, max): ");
438       scanf ("%lf %lf", &dmin, &dmax);
439       //      image_paint (CRTDEV, im_global, 0, 0, nx, ny, dmin, dmax, 0);
440       WAITKEY();
441
442       crt_put_stra ("Finished", 8 + C_WHITE);
443       crt_set_cpos (1, 2);
444     }
445
446   return (0);
447 }
448
449
450
451 #ifdef MPI_CT
452 void mpi_gather_image (IMAGE *im_global, IMAGE *im_local, const int opt_debug)
453 {
454   int iproc;
455   int end_work_unit;
456   int iw;
457
458   end_work_unit = mpi_ct.local_work_units[mpi_ct.my_rank] - 1;
459   for (iw = 0; iw <= end_work_unit; iw++) {
460     MPI_Send(im_local->v[iw], im_local->ny, MPI_FLOAT, 0, 0, mpi_ct.comm);
461   }
462
463   if (mpi_ct.my_rank == 0) {
464     for (iproc = 0; iproc < mpi_ct.nproc; iproc++) {
465       end_work_unit = mpi_ct.start_work_unit[iproc] 
466         + mpi_ct.local_work_units[iproc] 
467         - 1;
468
469       if (opt_debug) {
470         fprintf(stdout, "Recv rs data from process %d (%d-%d)\n", iproc,
471                 mpi_ct.start_work_unit[iproc], end_work_unit);
472         fflush(stdout);
473       }
474
475
476       for (iw = mpi_ct.start_work_unit[iproc]; iw <= end_work_unit; iw++) {
477         MPI_Status status;
478
479         MPI_Recv(im_global->v[iw], im_global->ny, MPI_FLOAT, iproc, 0, mpi_ct.comm, &status);
480       }
481     }
482   }
483
484 }
485 #endif