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