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