r89: *** empty log message ***
[ctsim.git] / src / ctrec.cpp
1 /*****************************************************************************
2 ** FILE IDENTIFICATION
3 **
4 **   Name:          ctrec.cpp
5 **   Purpose:       Reconstruct an image from projections
6 **   Programmer:    Kevin Rosenberg
7 **   Date Started:  Aug 1984
8 **
9 **  This is part of the CTSim program
10 **  Copyright (C) 1983-2000 Kevin Rosenberg
11 **
12 **  $Id: ctrec.cpp,v 1.4 2000/06/08 16:43:10 kevin Exp $
13 **
14 **  This program is free software; you can redistribute it and/or modify
15 **  it under the terms of the GNU General Public License (version 2) as
16 **  published by the Free Software Foundation.
17 **
18 **  This program is distributed in the hope that it will be useful,
19 **  but WITHOUT ANY WARRANTY; without even the implied warranty of
20 **  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
21 **  GNU General Public License for more details.
22 **
23 **  You should have received a copy of the GNU General Public License
24 **  along with this program; if not, write to the Free Software
25 **  Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
26 ******************************************************************************/
27
28 #include "ct.h"
29
30 enum {O_INTERP, O_FILTER, O_FILTER_PARAM, O_BACKPROJ, O_VERBOSE, O_TRACE, O_HELP, O_DEBUG, O_VERSION};
31
32 static struct option my_options[] =
33 {
34   {"interp", 1, 0, O_INTERP},
35   {"filter", 1, 0, O_FILTER},
36   {"filter-param", 1, 0, O_FILTER_PARAM},
37   {"backproj", 1, 0, O_BACKPROJ},
38   {"trace", 1, 0, O_TRACE},
39   {"debug", 0, 0, O_DEBUG},
40   {"verbose", 0, 0, O_VERBOSE},
41   {"help", 0, 0, O_HELP},
42   {"version", 0, 0, O_VERSION},
43   {0, 0, 0, 0}
44 };
45
46 void 
47 ctrec_usage (const char *program)
48 {
49   fprintf(stdout,"usage: %s raysum-file image-file nx-image ny-image [OPTIONS]\n", kbasename(program));
50   fprintf(stdout,"Image reconstruction from raysum projections\n");
51   fprintf(stdout,"\n");
52   fprintf(stdout,"   raysum-file     Input raysum file\n");
53   fprintf(stdout,"   image-file      Output image file in SDF2D format\n");
54   fprintf(stdout,"   nx-image        Number of columns in output image\n");
55   fprintf(stdout,"   ny-image        Number of rows in output image\n");
56   fprintf(stdout,"   --interp        Interpolation method during backprojection\n");
57   fprintf(stdout,"       nearest     Nearest neighbor interpolation\n");
58   fprintf(stdout,"       linear      Linear interpolation\n");
59 #if HAVE_BSPLINE_INTERP
60   fprintf(stdout,"       bspline     B-spline interpolation\n");
61 #endif
62   fprintf(stdout,"    --filter       Filter name\n");
63   fprintf(stdout,"       abs_bandlimit Abs * Bandlimiting (default)\n");
64   fprintf(stdout,"       abs_sinc      Abs * Sinc\n");
65   fprintf(stdout,"       abs_cos       Abs * Cosine\n");
66   fprintf(stdout,"       abs_hamming   Abs * Hamming\n");
67   fprintf(stdout,"       shepp         Shepp-Logan\n");
68   fprintf(stdout,"       bandlimit     Bandlimiting\n");
69   fprintf(stdout,"       sinc          Sinc\n");
70   fprintf(stdout,"       cos           Cosine\n");
71   fprintf(stdout,"       triangle      Triangle\n");
72   fprintf(stdout,"       hamming       Hamming\n");
73   fprintf(stdout,"    --backproj     Backprojection Method\n");
74   fprintf(stdout,"       trig        Trigometric functions at every point\n");
75   fprintf(stdout,"       table       Trigometric functions with precalculated table\n");
76   fprintf(stdout,"       diff        Difference method\n");
77   fprintf(stdout,"       diff2       Optimized difference method (default)\n");
78   fprintf(stdout,"       idiff2      Optimized difference method with integer math\n");
79   fprintf(stdout,"    --filter-param Alpha level for Hamming filter\n");
80   fprintf(stdout,"    --trace        Set tracing to level\n");
81   fprintf(stdout,"         none      No tracing (default)\n");
82   fprintf(stdout,"         text      Text level tracing\n");
83   fprintf(stdout,"         phm       Trace phantom\n");
84   fprintf(stdout,"         rays      Trace allrays\n");
85   fprintf(stdout,"         plot      Trace plotting\n");
86   fprintf(stdout,"         clipping  Trace clipping\n");
87   fprintf(stdout,"    --verbose      Turn on verbose mode\n");
88   fprintf(stdout,"    --debug        Turn on debug mode\n");
89   fprintf(stdout,"    --version      Print version\n");
90   fprintf(stdout,"    --help         Print this help message\n");
91 }
92
93
94 #ifdef HAVE_MPI
95 static void mpi_scatter_rs (RAYSUM *rs_global, RAYSUM *rs_local, const int debug);
96 #endif
97
98 static void print_raysum_info(const RAYSUM *rs);
99
100 int 
101 ctrec_main (const int argc, char *const argv[])
102 {
103   ImageFile *im_global = NULL;
104   RAYSUM *rs_global = NULL;
105   char *rs_name, *im_filename = NULL;
106   char remark[MAXREMARK];
107   char filt_name[80];
108   double time_start = 0, time_end = 0;
109   char *endptr;
110   int opt_verbose = 0;
111   int opt_debug = 0;
112   int opt_trace = TRACE_NONE;
113   double opt_filter_param = 1;
114   FilterType opt_filter = FILTER_ABS_BANDLIMIT;
115   InterpolationType opt_interp = I_LINEAR;
116   int opt_interp_param = 1;
117   BackprojType opt_backproj = O_BPROJ_DIFF2;
118   int nx, ny;
119 #ifdef HAVE_MPI
120   ImageFile *im_local;
121   RAYSUM *rs_local;
122   int mpi_argc = argc;
123   char **mpi_argv = (char **) argv;
124   int mpi_nview, mpi_ndet;
125   double mpi_detinc, mpi_rotinc, mpi_phmlen;
126   double mpi_t1, mpi_t2, mpi_t, mpi_t_g;
127
128   MPI::Init (mpi_argc, mpi_argv);
129   mpi_ct.comm = MPI::COMM_WORLD.Dup();
130   mpi_ct.nproc = mpi_ct.comm.Get_size();
131   mpi_ct.my_rank = mpi_ct.comm.Get_rank();
132
133   if (mpi_ct.nproc > CT_MPI_MAX_PROCESS) {
134     sys_error(ERR_FATAL, "Number of mpi processes (%d) exceeds max processes (%d)",
135               mpi_ct.nproc, CT_MPI_MAX_PROCESS);
136     exit(1);
137   }
138 #endif
139
140 #ifdef HAVE_MPI
141   time_start = MPI::Wtime();
142 #else
143   time_start = td_current_sec();
144 #endif
145
146 #ifdef HAVE_MPI
147   if (mpi_ct.my_rank == 0) {
148 #endif
149     while (1) {
150       int c = getopt_long(argc, argv, "", my_options, NULL);
151       char *endptr = NULL;
152       
153       if (c == -1)
154         break;
155       
156       switch (c)
157         {
158         case O_INTERP:
159           if ((opt_interp = opt_set_interpolation(optarg)) < 0) {
160             ctrec_usage(argv[0]);
161             return (1);
162           }
163           break;
164         case O_FILTER:
165           if ((opt_filter = opt_set_filter(optarg)) < 0) {
166             ctrec_usage(argv[0]);
167             return (1);
168           }
169           break;
170         case O_BACKPROJ:
171           if ((opt_backproj = opt_set_backproj(optarg)) < 0) {
172               ctrec_usage(argv[0]);
173               return (1);
174           }
175           break;
176         case O_FILTER_PARAM:
177           opt_filter_param = strtod(optarg, &endptr);
178           if (endptr != optarg + strlen(optarg)) {
179             ctrec_usage(argv[0]);
180           }
181           break;
182         case O_VERBOSE:
183           opt_verbose = 1;
184           break;
185         case O_DEBUG:
186           opt_debug = 1;
187           break;
188         case O_TRACE:
189           if ((opt_trace = opt_set_trace(optarg)) < 0) {
190             ctrec_usage(argv[0]);
191             return (1);
192           }
193           break;
194         case O_VERSION:
195 #ifdef VERSION
196           fprintf(stdout, "Version %s\n", VERSION);
197 #else
198           fprintf(stderr, "Unknown version number");
199 #endif
200           exit(0);
201         case O_HELP:
202         case '?':
203           ctrec_usage(argv[0]);
204           return (0);
205         default:
206           ctrec_usage(argv[0]);
207           return (1);
208         }
209     }
210   
211     if (optind + 4 != argc) {
212       ctrec_usage(argv[0]);
213       return (1);
214     }
215
216     rs_name = argv[optind];
217   
218     im_filename = argv[optind + 1];
219   
220     nx = strtol(argv[optind + 2], &endptr, 10);
221     ny = strtol(argv[optind + 3], &endptr, 10);
222   
223     if (opt_filter == FILTER_G_HAMMING || opt_filter == FILTER_ABS_G_HAMMING)
224       snprintf (filt_name, sizeof(filt_name), "%s: alpha = %.2f",
225                name_of_filter (opt_filter), opt_filter_param); 
226     else
227       snprintf (filt_name, sizeof(filt_name), "%s", name_of_filter (opt_filter));
228   
229     snprintf (remark, sizeof(remark), "Reconstruct: %dx%d, %s, %s, %s",
230              nx, ny, filt_name, name_of_interpolation (opt_interp), name_of_backproj(opt_backproj));
231   
232     if (opt_verbose)
233       fprintf (stdout, "%s\n", remark);
234 #ifdef HAVE_MPI
235   }
236 #endif
237
238 #ifdef HAVE_MPI
239   if (mpi_ct.my_rank == 0) {
240     rs_global = raysum_open (rs_name);
241     raysum_read (rs_global);
242     if (opt_verbose)
243       print_raysum_info(rs_global);
244
245     mpi_ndet = rs_global->ndet;
246     mpi_nview = rs_global->nview;
247     mpi_detinc = rs_global->det_inc;
248     mpi_phmlen = rs_global->phmlen;
249     mpi_rotinc = rs_global->rot_inc;
250   }
251
252   mpi_t1 = MPI::Wtime();
253   mpi_ct.comm.Bcast (&opt_verbose, 1, MPI::INT, 0);
254   mpi_ct.comm.Bcast (&opt_debug, 1, MPI::INT, 0);
255   mpi_ct.comm.Bcast (&opt_trace, 1, MPI::INT, 0);
256   mpi_ct.comm.Bcast (&opt_filter, 1, MPI::INT, 0);
257   mpi_ct.comm.Bcast (&opt_interp, 1, MPI::INT, 0);
258   mpi_ct.comm.Bcast (&opt_filter_param, 1, MPI::DOUBLE, 0);
259   mpi_ct.comm.Bcast (&opt_interp_param, 1, MPI::INT, 0);
260   mpi_ct.comm.Bcast (&opt_backproj, 1, MPI::INT, 0);
261   mpi_ct.comm.Bcast (&mpi_ndet, 1, MPI::INT, 0);
262   mpi_ct.comm.Bcast (&mpi_nview, 1, MPI::INT, 0);
263   mpi_ct.comm.Bcast (&mpi_detinc, 1, MPI::DOUBLE, 0);
264   mpi_ct.comm.Bcast (&mpi_phmlen, 1, MPI::DOUBLE, 0);
265   mpi_ct.comm.Bcast (&mpi_rotinc, 1, MPI::DOUBLE, 0);
266   mpi_ct.comm.Bcast (&nx, 1, MPI::INT, 0);
267   mpi_ct.comm.Bcast (&ny, 1, MPI::INT, 0);
268   if (opt_verbose) {
269     mpi_t2 = MPI::Wtime();
270     mpi_t = mpi_t2 - mpi_t1;
271     mpi_ct.comm.Reduce(&mpi_t, &mpi_t_g, 1, MPI::DOUBLE, MPI::MAX, 0);
272     if (mpi_ct.my_rank == 0)
273       printf("Time to Bcast vars = %f secs, Max time = %f\n", mpi_t, mpi_t_g);
274   }
275
276   mpi_ct_calc_work_units(mpi_nview);
277
278   rs_local = raysum_create (NULL, mpi_ct.local_work_units[mpi_ct.my_rank], mpi_ndet);
279
280   rs_local->ndet = mpi_ndet;
281   rs_local->nview = mpi_nview;
282   rs_local->det_inc = mpi_detinc;
283   rs_local->phmlen = mpi_phmlen;
284   rs_local->rot_inc = mpi_rotinc;
285
286   if (opt_verbose)
287     mpi_t1 = MPI::Wtime();
288   mpi_scatter_rs(rs_global, rs_local, opt_debug);
289   if (opt_verbose) {
290     mpi_t2 = MPI::Wtime();
291     mpi_t = mpi_t2 - mpi_t1;
292     mpi_ct.comm.Reduce(&mpi_t, &mpi_t_g, 1, MPI::DOUBLE, MPI::MAX, 0);
293     if (mpi_ct.my_rank == 0)
294       printf("Time to scatter rs = %f secs, Max time = %f sec\n", mpi_t, mpi_t_g);
295   }
296
297   if (mpi_ct.my_rank == 0) {
298     im_global = new ImageFile (im_filename, nx, ny);
299     im_global->adf.fileCreate();
300   }
301
302   im_local = new ImageFile (nx, ny);
303 #else
304   rs_global = raysum_open (rs_name);
305   raysum_read (rs_global);
306   if (opt_verbose)
307     print_raysum_info(rs_global);
308
309   im_global = new ImageFile (im_filename, nx, ny);
310   im_global->adf.fileCreate();
311 #endif
312
313 #ifdef HAVE_MPI
314   mpi_t1 = MPI::Wtime();
315   proj_reconst (*im_local, rs_local, opt_filter, opt_filter_param, 
316                  opt_interp, opt_interp_param, opt_backproj, opt_trace);
317
318   mpi_t2 = MPI::Wtime();
319   mpi_t = mpi_t2 - mpi_t1;
320   mpi_ct.comm.Reduce(&mpi_t, &mpi_t_g, 1, MPI::DOUBLE, MPI::MAX, 0);
321   if (mpi_ct.my_rank == 0 && opt_verbose)
322     printf("Time to reconstruct = %f, Max time = %f\n", mpi_t, mpi_t_g);
323 #else
324   proj_reconst (*im_global, rs_global, opt_filter, opt_filter_param, 
325                  opt_interp, opt_interp_param, opt_backproj, opt_trace);
326 #endif
327
328 #ifdef HAVE_MPI
329   if (opt_verbose)
330     mpi_t1 = MPI::Wtime();
331
332   int nxLocal = im_local->adf.nx();
333   int nyLocal = im_local->adf.ny();
334   ImageFileArray vLocal = im_local->getArray();
335   ImageFileArray vGlobal = NULL;
336   if (mpi_ct.my_rank == 0)
337     vGlobal = im_global->getArray();
338
339   for (int ix = 0; ix < nxLocal; ix++) {
340     void *recvbuf = NULL;
341     if (mpi_ct.my_rank == 0)
342       recvbuf = vGlobal[ix];
343
344     mpi_ct.comm.Reduce(vLocal[ix], recvbuf, nyLocal, im_local->getMPIDataType(), MPI::SUM, 0);
345   }
346
347   if (opt_verbose) {
348     mpi_t2 = MPI::Wtime();
349     mpi_t = mpi_t2 - mpi_t1;
350     mpi_ct.comm.Reduce (&mpi_t, &mpi_t_g, 1, MPI::DOUBLE, MPI::MAX, 0);
351     if (mpi_ct.my_rank == 0)
352       printf("Time to reduce image = %f secs, max time = %f\n", mpi_t, mpi_t_g);
353   }
354
355   if (mpi_ct.my_rank == 0)
356       time_end = MPI::Wtime();
357 #else
358   time_end = td_current_sec();
359 #endif
360
361
362 #ifdef HAVE_MPI
363   if (mpi_ct.my_rank == 0)
364 #endif
365     {
366       raysum_close (rs_global);
367       double calctime = time_end - time_start;
368       im_global->adf.arrayDataWrite ();
369       im_global->adf.labelAdd (Array2dFileLabel::L_HISTORY, rs_global->remark, rs_global->calctime);
370       im_global->adf.labelAdd (Array2dFileLabel::L_HISTORY, remark, calctime);
371       im_global->adf.fileClose ();
372       if (opt_verbose)
373         cout << "Time active = " << calctime << " sec" << endl;
374     }
375
376 #ifdef HAVE_MPI
377         MPI::Finalize();
378 #endif
379
380         return (0);
381 }
382
383
384 #ifdef HAVE_MPI
385 static void mpi_scatter_rs (RAYSUM *rs_global, RAYSUM *rs_local, const int opt_debug)
386 {
387   int iproc;
388   int end_work_unit;
389   int iw = 0;
390
391   if (mpi_ct.my_rank == 0) {
392     for (iproc = 0; iproc < mpi_ct.nproc; iproc++) {
393       end_work_unit = mpi_ct.start_work_unit[iproc] + mpi_ct.local_work_units[iproc]  - 1;
394
395       for (iw = mpi_ct.start_work_unit[iproc]; iw <= end_work_unit; iw++) {
396         mpi_ct.comm.Send(&rs_global->view[iw]->ndet, 1, MPI::INT, iproc, 0);
397         mpi_ct.comm.Send(&rs_global->view[iw]->view_angle, 1, MPI::DOUBLE, iproc, 0);
398         mpi_ct.comm.Send(rs_global->view[iw]->detval, rs_global->ndet, MPI::FLOAT, iproc, 0);
399       }
400     }
401   }
402
403   if (opt_debug)
404     fprintf(stdout, "Receiving rs data in process %d\n", mpi_ct.my_rank);
405
406   end_work_unit = mpi_ct.local_work_units[mpi_ct.my_rank] - 1;
407   for (iw = 0; iw <= end_work_unit; iw++) {
408     MPI::Status status;
409
410     mpi_ct.comm.Recv(&rs_local->view[iw]->ndet, 1, MPI::INT, 0, 0, status);
411     mpi_ct.comm.Recv(&rs_local->view[iw]->view_angle, 1, MPI::DOUBLE, 0, 0, status);
412     mpi_ct.comm.Recv(rs_local->view[iw]->detval, rs_local->ndet, MPI::FLOAT, 0, 0, status);
413   }
414   rs_local->nview = mpi_ct.local_work_units[mpi_ct.my_rank];
415   if (opt_debug) {
416     mpi_ct.comm.Barrier();
417     fprintf(stdout, "Done with mpi_scatter_rs in process %2d\n", mpi_ct.my_rank);
418   }
419 }
420
421 #endif
422
423 static void print_raysum_info(const RAYSUM *rs)
424 {
425   printf ("Number of detectors: %d\n", rs->ndet);
426   printf ("    Number of views: %d\n", rs->nview);
427   printf ("             Remark: %s\n", rs->remark);
428   printf ("             phmlen: %f\n", rs->phmlen);
429   printf ("          det_start: %f\n", rs->det_start);
430   printf ("            det_inc: %f\n", rs->det_inc);
431   printf ("          rot_start: %f\n", rs->rot_start);
432   printf ("            rot_inc: %f\n", rs->rot_inc);
433 }
434
435 #ifndef NO_MAIN
436 int 
437 main (const int argc, char *const argv[])
438 {
439   return (ctrec_main(argc, argv));
440 }
441 #endif
442