3c9138724fecc79cb2879a4ba56af41a281b64b3
[ctsim.git] / src / ctrec.c
1 /*****************************************************************************
2 **  This is part of the CTSim program
3 **  Copyright (C) 1983-2000 Kevin Rosenberg
4 **
5 **  $Id: ctrec.c,v 1.2 2000/04/29 23:24:56 kevin Exp $
6 **  $Log: ctrec.c,v $
7 **  Revision 1.2  2000/04/29 23:24:56  kevin
8 **  *** empty log message ***
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
29 /* FILE
30  *   ctrec.c            Reconstruct an image from raysums
31  *
32  * DATE
33  *   Aug 84
34  *   Jul 99  -- Converted to ANSI C
35  *              Added MPI parallel processing
36  */
37
38 #include "ct.h"
39
40 #define O_INTERP       1
41 #define O_FILTER       2
42 #define O_FILTER_PARAM 3
43 #define O_BACKPROJ     4
44 #define O_VERBOSE      5
45 #define O_TRACE        6
46 #define O_HELP         7
47 #define O_DEBUG        8
48 #define O_VERSION      9
49
50 static struct option my_options[] =
51 {
52   {"interp", 1, 0, O_INTERP},
53   {"filter", 1, 0, O_FILTER},
54   {"filter-param", 1, 0, O_FILTER_PARAM},
55   {"backproj", 1, 0, O_BACKPROJ},
56   {"trace", 1, 0, O_TRACE},
57   {"debug", 0, 0, O_DEBUG},
58   {"verbose", 0, 0, O_VERBOSE},
59   {"help", 0, 0, O_HELP},
60   {"version", 0, 0, O_VERSION},
61   {0, 0, 0, 0}
62 };
63
64 void 
65 usage (const char *program)
66 {
67   fprintf(stdout,"usage: %s raysum-file image-file nx-image ny-image [OPTIONS]\n", kbasename(program));
68   fprintf(stdout,"Image reconstruction from raysum projections\n");
69   fprintf(stdout,"\n");
70   fprintf(stdout,"   raysum-file     Input raysum file\n");
71   fprintf(stdout,"   image-file      Output image file in SDF2D format\n");
72   fprintf(stdout,"   nx-image        Number of columns in output image\n");
73   fprintf(stdout,"   ny-image        Number of rows in output image\n");
74   fprintf(stdout,"   --interp        Interpolation method during backprojection\n");
75   fprintf(stdout,"       nearest     Nearest neighbor interpolation\n");
76   fprintf(stdout,"       linear      Linear interpolation\n");
77   fprintf(stdout,"       bspline     B-spline interpolation\n");
78   fprintf(stdout,"    --filter       Filter name\n");
79   fprintf(stdout,"       abs_bandlimit Abs * Bandlimiting (default)\n");
80   fprintf(stdout,"       abs_sinc      Abs * Sinc\n");
81   fprintf(stdout,"       abs_cos       Abs * Cosine\n");
82   fprintf(stdout,"       abs_hamming   Abs * Hamming\n");
83   fprintf(stdout,"       shepp         Shepp-Logan\n");
84   fprintf(stdout,"       bandlimit     Bandlimiting\n");
85   fprintf(stdout,"       sinc          Sinc\n");
86   fprintf(stdout,"       cos           Cosine\n");
87   fprintf(stdout,"       triangle      Triangle\n");
88   fprintf(stdout,"       hamming       Hamming\n");
89   fprintf(stdout,"    --backproj     Backprojection Method\n");
90   fprintf(stdout,"       trig        Trigometric functions at every point\n");
91   fprintf(stdout,"       table       Trigometric functions with precalculated table\n");
92   fprintf(stdout,"       diff        Difference method\n");
93   fprintf(stdout,"       diff2       Optimized difference method (default)\n");
94   fprintf(stdout,"       idiff2      Optimized difference method with integer math\n");
95   fprintf(stdout,"    --filter-param Alpha level for Hamming filter\n");
96   fprintf(stdout,"    --trace        Set tracing to level\n");
97   fprintf(stdout,"         none      No tracing (default)\n");
98   fprintf(stdout,"         status    Status tracing\n");
99   fprintf(stdout,"         pic       Trace picture\n");
100   fprintf(stdout,"         rays      Trace allrays\n");
101   fprintf(stdout,"         plot      Trace plotting\n");
102   fprintf(stdout,"         clipping  Trace clipping\n");
103   fprintf(stdout,"    --verbose      Turn on verbose mode\n");
104   fprintf(stdout,"    --debug        Turn on debug mode\n");
105   fprintf(stdout,"    --version      Print version\n");
106   fprintf(stdout,"    --help         Print this help message\n");
107   exit(1);  
108 }
109
110
111 #ifdef MPI_CT
112 static void mpi_scatter_rs (RAYSUM *rs_global, RAYSUM *rs_local, const int debug);
113 #endif
114 static void print_raysum_info(const RAYSUM *rs);
115
116 int 
117 main (const int argc, char *const argv[])
118 {
119   IMAGE *im_global = NULL;
120   RAYSUM *rs_global = NULL;
121   char *rs_name, *im_filename = NULL;
122   char remark[MAXREMARK];
123   char filt_name[80];
124   int nx, ny;
125   double time_start = 0, time_end = 0;
126   char *endptr;
127   int opt_verbose = 0;
128   int opt_debug = 0;
129   int opt_trace = TRACE_NONE;
130   double opt_filter_param = 1;
131   int opt_filter = W_A_BANDLIMIT;
132   int opt_interp = I_LINEAR;
133   int opt_interp_param = 1;
134   int opt_backproj = O_BPROJ_DIFF2;
135 #ifdef MPI_CT
136   IMAGE *im_local;
137   RAYSUM *rs_local;
138   int mpi_argc = argc;
139   char **mpi_argv = (char **) argv;
140   int mpi_nview, mpi_ndet;
141   double mpi_detinc, mpi_rotinc, mpi_piclen;
142   double mpi_t1, mpi_t2, mpi_t, mpi_t_g;
143   int ix;
144
145   MPI_Init(&mpi_argc, &mpi_argv);
146   MPI_Comm_dup (MPI_COMM_WORLD, &mpi_ct.comm);
147   MPI_Comm_size(mpi_ct.comm, &mpi_ct.nproc);
148   MPI_Comm_rank(mpi_ct.comm, &mpi_ct.my_rank);
149
150   if (mpi_ct.nproc > MPI_MAX_PROCESS) {
151     sys_error(ERR_FATAL, "Number of mpi processes (%d) exceeds max processes (%d)",
152               mpi_ct.nproc, MPI_MAX_PROCESS);
153     exit(1);
154   }
155 #endif
156
157 #ifdef MPI_CT
158   time_start = MPI_Wtime();
159 #else
160   time_start = td_current_sec();
161 #endif
162
163 #ifdef MPI_CT
164   if (mpi_ct.my_rank == 0) {
165 #endif
166     while (1) {
167       int c = getopt_long(argc, argv, "", my_options, NULL);
168       char *endptr = NULL;
169       
170       if (c == -1)
171         break;
172       
173       switch (c)
174         {
175         case O_INTERP:
176           opt_interp = opt_set_interpolation(optarg, argv[0]);
177           break;
178         case O_FILTER:
179           opt_filter = opt_set_filter(optarg, argv[0]);
180           break;
181         case O_BACKPROJ:
182           opt_backproj = opt_set_backproj(optarg, argv[0]);
183           break;
184         case O_FILTER_PARAM:
185           opt_filter_param = strtod(optarg, &endptr);
186           if (endptr != optarg + strlen(optarg)) {
187             usage(argv[0]);
188             exit(1);
189           }
190           break;
191         case O_VERBOSE:
192           opt_verbose = 1;
193           break;
194         case O_DEBUG:
195           opt_debug = 1;
196           break;
197         case O_TRACE:
198           opt_trace = opt_set_trace(optarg, argv[0]);
199           break;
200         case O_VERSION:
201 #ifdef VERSION
202           fprintf(stdout, "Version %s\n", VERSION);
203 #else
204           fprintf(stderr, "Unknown version number");
205 #endif
206           exit(0);
207         case O_HELP:
208         case '?':
209           usage(argv[0]);
210           exit(0);
211         default:
212           usage(argv[0]);
213           exit(1);
214         }
215     }
216   
217     if (optind + 4 != argc) {
218       usage(argv[0]);
219       exit(1);
220     }
221
222     rs_name = argv[optind];
223   
224     im_filename = argv[optind + 1];
225   
226     nx = strtol(argv[optind + 2], &endptr, 10);
227     ny = strtol(argv[optind + 3], &endptr, 10);
228   
229     if (opt_filter == W_G_HAMMING || opt_filter == W_AG_HAMMING)
230       sprintf (filt_name, "%s: alpha = %.2f",
231                filter_name_of (opt_filter), opt_filter_param); 
232     else
233       sprintf (filt_name, "%s", filter_name_of (opt_filter));
234   
235     sprintf (remark,
236              "Reconstruct: %dx%d, %s, %s, %s",
237              nx, ny, filt_name, interp_name_of (opt_interp), name_of_backproj(opt_backproj));
238   
239     if (opt_verbose)
240       fprintf (stdout, "%s\n", remark);
241 #ifdef MPI_CT
242   }
243 #endif
244
245 #ifdef MPI_CT
246   if (mpi_ct.my_rank == 0) {
247     rs_global = raysum_open (rs_name);
248     raysum_read (rs_global);
249     print_raysum_info(rs_global);
250     if (opt_verbose) {
251       printf ("Number of detectors: %d\n", rs_global->ndet);
252       printf ("    Number of views: %d\n", rs_global->nview);
253       printf ("             Remark: %s\n", rs_global->remark);
254     }
255   
256     mpi_ndet = rs_global->ndet;
257     mpi_nview = rs_global->nview;
258     mpi_detinc = rs_global->det_inc;
259     mpi_piclen = rs_global->piclen;
260     mpi_rotinc = rs_global->rot_inc;
261   }
262
263   mpi_t1 = MPI_Wtime();
264   MPI_Bcast(&opt_verbose, 1, MPI_INT, 0, mpi_ct.comm);
265   MPI_Bcast(&opt_debug, 1, MPI_INT, 0, mpi_ct.comm);
266   MPI_Bcast(&opt_trace, 1, MPI_INT, 0, mpi_ct.comm);
267   MPI_Bcast(&opt_filter, 1, MPI_INT, 0, mpi_ct.comm);
268   MPI_Bcast(&opt_interp, 1, MPI_INT, 0, mpi_ct.comm);
269   MPI_Bcast(&opt_filter_param, 1, MPI_DOUBLE, 0, mpi_ct.comm);
270   MPI_Bcast(&opt_interp_param, 1, MPI_INT, 0, mpi_ct.comm);
271   MPI_Bcast(&opt_backproj, 1, MPI_INT, 0, mpi_ct.comm);
272   MPI_Bcast(&mpi_ndet, 1, MPI_INT, 0, mpi_ct.comm);
273   MPI_Bcast(&mpi_nview, 1, MPI_INT, 0, mpi_ct.comm);
274   MPI_Bcast(&mpi_detinc, 1, MPI_DOUBLE, 0, mpi_ct.comm);
275   MPI_Bcast(&mpi_piclen, 1, MPI_DOUBLE, 0, mpi_ct.comm);
276   MPI_Bcast(&mpi_rotinc, 1, MPI_DOUBLE, 0, mpi_ct.comm);
277   MPI_Bcast(&nx, 1, MPI_INT, 0, mpi_ct.comm);
278   MPI_Bcast(&ny, 1, MPI_INT, 0, mpi_ct.comm);
279   if (opt_verbose) {
280     mpi_t2 = MPI_Wtime();
281     mpi_t = mpi_t2 - mpi_t1;
282     MPI_Reduce(&mpi_t, &mpi_t_g, 1, MPI_DOUBLE, MPI_MAX, 0, mpi_ct.comm);
283     if (mpi_ct.my_rank == 0)
284       printf("Time to Bcast vars = %f secs, Max time = %f\n", mpi_t, mpi_t_g);
285   }
286
287   mpi_ct_calc_work_units(mpi_nview);
288
289   if (opt_debug) {
290     fprintf(stdout, "Calc'd local work units process %d: nviews=%d, local_work_units=%d, start_work_units=%d\n",
291             mpi_ct.my_rank, mpi_nview, mpi_ct.local_work_units[mpi_ct.my_rank], mpi_ct.start_work_unit[mpi_ct.my_rank]);
292     MPI_Barrier(mpi_ct.comm);
293   }
294
295   rs_local = raysum_create (NULL, mpi_ct.local_work_units[mpi_ct.my_rank], mpi_ndet);
296   if (opt_debug) {
297     fprintf(stdout, "Created rs_local %lx for process %d: local views=%4d, local dets=%4d\n", (unsigned long int) rs_local, mpi_ct.my_rank, mpi_ct.local_work_units[mpi_ct.my_rank], mpi_ndet);
298     MPI_Barrier(mpi_ct.comm);
299   }
300
301   rs_local->ndet = mpi_ndet;
302   rs_local->nview = mpi_nview;
303   rs_local->det_inc = mpi_detinc;
304   rs_local->piclen = mpi_piclen;
305   rs_local->rot_inc = mpi_rotinc;
306
307   if (opt_verbose)
308     mpi_t1 = MPI_Wtime();
309   mpi_scatter_rs(rs_global, rs_local, opt_debug);
310   if (opt_verbose) {
311     mpi_t2 = MPI_Wtime();
312     mpi_t = mpi_t2 - mpi_t1;
313     MPI_Reduce(&mpi_t, &mpi_t_g, 1, MPI_DOUBLE, MPI_MAX, 0, mpi_ct.comm);
314     if (mpi_ct.my_rank == 0)
315       printf("Time to scatter rs = %f secs, Max time = %f sec\n", mpi_t, mpi_t_g);
316   }
317
318   if (mpi_ct.my_rank == 0) {
319     im_global = image_create (im_filename, nx, ny);
320     sdf_add_label (LT_HISTORY, rs_global->remark, rs_global->calctime, im_global->dfp_2d->dfp);
321     sdf_add_empty_label (im_global->dfp_2d->dfp);
322   }
323   im_local = image_create (NULL, nx, ny);
324
325 #else
326   rs_global = raysum_open (rs_name);
327   raysum_read (rs_global);
328   if (opt_debug)
329     print_raysum_info(rs_global);
330
331   im_global = image_create (im_filename, nx, ny);
332   sdf_add_label (LT_HISTORY, rs_global->remark, rs_global->calctime, im_global->dfp_2d->dfp);
333   sdf_add_empty_label (im_global->dfp_2d->dfp);
334 #endif
335
336 #ifdef MPI_CT
337   mpi_t1 = MPI_Wtime();
338   image_reconst (im_local, rs_local, opt_filter, opt_filter_param, 
339                  opt_interp, opt_interp_param, opt_backproj, opt_trace);
340   mpi_t2 = MPI_Wtime();
341   mpi_t = mpi_t2 - mpi_t1;
342   MPI_Reduce(&mpi_t, &mpi_t_g, 1, MPI_DOUBLE, MPI_MAX, 0, mpi_ct.comm);
343   if (mpi_ct.my_rank == 0) {
344     printf("Time to reconstruct = %f, Max time = %f\n", mpi_t, mpi_t_g);
345   }
346 #else
347   image_reconst (im_global, rs_global, opt_filter, opt_filter_param, 
348                  opt_interp, opt_interp_param, opt_backproj, opt_trace);
349 #endif
350
351 #ifdef MPI_CT
352   if (mpi_ct.my_rank == 0) {
353     raysum_close (rs_global);
354   }
355
356   if (opt_verbose)
357     mpi_t1 = MPI_Wtime();
358   for (ix = 0; ix < im_local->nx; ix++) {
359     MPI_Reduce(im_local->v[ix], im_global->v[ix], 
360                im_local->ny, MPI_FLOAT, MPI_SUM, 0, mpi_ct.comm);
361   }
362   if (opt_verbose) {
363     mpi_t2 = MPI_Wtime();
364     mpi_t = mpi_t2 - mpi_t1;
365     MPI_Reduce(&mpi_t, &mpi_t_g, 1, MPI_DOUBLE, MPI_MAX, 0, mpi_ct.comm);
366     if (mpi_ct.my_rank == 0)
367       printf("Time to reduce image = %f secs, max time = %f\n", mpi_t, mpi_t_g);
368   }
369   if (mpi_ct.my_rank == 0) {
370     strncpy (im_global->remark, remark, MAXREMARK);
371     image_save (im_global);
372     time_end = MPI_Wtime();
373     im_global->calctime = time_end - time_start;
374     if (opt_verbose)
375       fprintf (stdout, "Time active = %.2f\n", im_global->calctime);
376     }
377 #else  
378   raysum_close (rs_global);
379   strncpy (im_global->remark, remark, MAXREMARK);
380   image_save (im_global);
381
382   time_end = td_current_sec();
383   im_global->calctime = time_end - time_start;
384   if (opt_verbose)
385     fprintf (stdout, "Time active = %.2f secs\n", im_global->calctime);
386 #endif
387
388 #ifdef MPI_CT
389         MPI_Finalize();
390 #endif
391
392         return (0);
393 }
394
395
396 #ifdef MPI_CT
397 static void mpi_scatter_rs (RAYSUM *rs_global, RAYSUM *rs_local, const int opt_debug)
398 {
399   int iproc;
400   int end_work_unit;
401   int iw = 0;
402
403   if (mpi_ct.my_rank == 0) {
404     for (iproc = 0; iproc < mpi_ct.nproc; iproc++) {
405       end_work_unit = mpi_ct.start_work_unit[iproc] 
406         + mpi_ct.local_work_units[iproc] 
407         - 1;
408
409       if (opt_debug)
410         fprintf(stdout, "Sending rs data to process %d\n", iproc);
411
412       for (iw = mpi_ct.start_work_unit[iproc]; iw <= end_work_unit; iw++) {
413         if (opt_debug)
414           fprintf(stdout, "Sending from process %2d to process %2d: view_angle=%8f, ndet=%5d\n", mpi_ct.my_rank, iproc, rs_global->view[iw]->view_angle, rs_global->view[iw]->ndet);
415         MPI_Send(&rs_global->view[iw]->ndet, 1, MPI_INT, iproc, 0, mpi_ct.comm);
416         MPI_Send(&rs_global->view[iw]->view_angle, 1, MPI_DOUBLE, iproc, 0, mpi_ct.comm);
417         MPI_Send(rs_global->view[iw]->detval, rs_global->ndet, MPI_FLOAT, iproc, 0, mpi_ct.comm);
418       }
419     }
420   }
421
422   if (opt_debug) {
423     fprintf(stdout, "Receiving rs data in process %d\n", mpi_ct.my_rank);
424     MPI_Barrier(mpi_ct.comm);
425   }
426
427   end_work_unit = mpi_ct.local_work_units[mpi_ct.my_rank] - 1;
428   for (iw = 0; iw <= end_work_unit; iw++) {
429     MPI_Status status;
430
431     if (opt_debug) {
432       fprintf(stdout,"Receiving into rs_local from mpi_scatter_rs, process %2d: rs_local=%lx, ", mpi_ct.my_rank, (unsigned long int) rs_local);
433       fprintf(stdout,"iw=%4d, nrot=%4d, ndet=%4d, view=%lx, detval=%lx\n", iw, rs_local->nview, rs_local->ndet, (unsigned long int) rs_local->view[iw], (unsigned long int) rs_local->view[iw]->detval);
434     }
435
436     // MPI_Recv(&rs_local->view[iw]->ndet, 1, MPI_INT, 0, 0, mpi_ct.comm, &status);
437     // MPI_Recv(&rs_local->view[iw]->view_angle, 1, MPI_DOUBLE, 0, 0, mpi_ct.comm, &status);
438     //    MPI_Recv(rs_local->view[iw]->detval, rs_global->ndet, MPI_FLOAT, 0, 0, mpi_ct.comm, &status);
439
440     if (opt_debug) {
441       fprintf(stdout, "Received view_angle=%8f, ndet=%5d\n", rs_local->view[iw]->view_angle, rs_local->view[iw]->ndet);
442     }
443   }
444   rs_local->nview = mpi_ct.local_work_units[mpi_ct.my_rank];
445   if (opt_debug) {
446     MPI_Barrier(MPI_COMM_WORLD);
447     fprintf(stdout, "Done with mpi_scatter_rs in process %2d\n", mpi_ct.my_rank);
448     exit(0);
449   }
450 }
451
452 #endif
453
454 static void print_raysum_info(const RAYSUM *rs)
455 {
456   printf ("Number of detectors: %d\n", rs->ndet);
457   printf ("    Number of views: %d\n", rs->nview);
458   printf ("             Remark: %s\n", rs->remark);
459   printf ("Piclen: %f\n", rs->piclen);
460   printf ("det_start: %f\n", rs->det_start);
461   printf ("det_inc: %f\n", rs->det_inc);
462   printf ("rot_start: %f\n", rs->rot_start);
463   printf ("rot_inc: %f\n", rs->rot_inc);
464 }