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