r4203: *** empty log message ***
authorKevin M. Rosenberg <kevin@rosenberg.net>
Sat, 15 Mar 2003 10:29:06 +0000 (10:29 +0000)
committerKevin M. Rosenberg <kevin@rosenberg.net>
Sat, 15 Mar 2003 10:29:06 +0000 (10:29 +0000)
NEWS
configure
configure.ac
debian/changelog
debian/upload.sh
libctsim/projections.cpp
src/dialogs.h
src/views.cpp
src/views.h

diff --git a/NEWS b/NEWS
index 6f6c8d0432f9723fea22d4d792e448829e8a5949..d65c241329dbb33a4d6e56857231cc2c6948b22f 100644 (file)
--- a/NEWS
+++ b/NEWS
@@ -1,3 +1,7 @@
+Version 4.2 Released Mar 2003
+
+* Direct Fourier Reconstructions now supported
+
 Version 4.0 Released Jan 2003
 
 * Port to wxWindows 2.4
index df9d60b3dac093affb71a694673376af989ff824..26a4a05a6088875a0d74fd1ea1f40ad3102c7ded 100755 (executable)
--- a/configure
+++ b/configure
@@ -1492,7 +1492,7 @@ fi
 
 PACKAGE=ctsim
 
-VERSION=4.0.2
+VERSION=4.2.0
 
 if test "`cd $srcdir && pwd`" != "`pwd`" && test -f $srcdir/config.status; then
   { { echo "$as_me:$LINENO: error: source directory already configured; run \"make distclean\" there first" >&5
index 1f4c7b3a00a3dd63680e3776f106fef34908e43e..87ba96ab1b7e0a6ddd3cf0d7258e5e62d46a06b2 100644 (file)
@@ -5,7 +5,7 @@ dnl CDPATH=
 
 AC_INIT
 AC_CONFIG_SRCDIR([src/ctsim.cpp])
-AM_INIT_AUTOMAKE(ctsim,4.0.2)
+AM_INIT_AUTOMAKE(ctsim,4.2.0)
 AM_CONFIG_HEADER(config.h)
 
 dnl Checks for programs.
index 430b28f7f6a8a414925aeebf8f20b106c70360b0..dcb5eba97cc05ab860e3f1dcbca1374e190e9bf1 100644 (file)
@@ -1,3 +1,9 @@
+ctsim (4.2.0-1) unstable; urgency=low
+
+  * New feature: Direct inverse Fourier reconstructions
+
+ -- Kevin M. Rosenberg <kmr@debian.org>  Fri, 14 Mar 2003 23:12:06 -0700
+
 ctsim (4.1.2-3) unstable; urgency=low
 
   * Add build-depends-indep to control file
index 14e31d69bb7fe64f595acce125bc62cfc4fe8978..dac42880cdb4738e2113c2b726b2e6507330eb77 100755 (executable)
@@ -1,4 +1,4 @@
 #!/bin/bash -e
 
-dup ctsim -Uftp.med-info.com -D/home/ftp/ctsim -C"(cd /opt/apache/htdocs/ctsim.org; make install)" -su $*
+dup ctsim -Uftp.med-info.com -D/home/ftp/ctsim -C"(cd /opt/apache/htdocs/ctsim; make install)" -su $*
 
index 612ab6bf9161663dd42d6894e3bb1cfc371fc5de..db8fd0c3d60a63632acbea37bacdd08387d3b8d5 100644 (file)
@@ -8,7 +8,7 @@
 **  This is part of the CTSim program
 **  Copyright (c) 1983-2001 Kevin Rosenberg
 **
-**  $Id: projections.cpp,v 1.80 2002/06/27 03:19:23 kevin Exp $
+**  $Id: projections.cpp,v 1.81 2003/03/15 10:27:30 kevin Exp $
 **
 **  This program is free software; you can redistribute it and/or modify
 **  it under the terms of the GNU General Public License (version 2) as
@@ -25,7 +25,7 @@
 ******************************************************************************/
 
 #include "ct.h"
-#include <ctime>\r
+#include <ctime>
 
 const kuint16 Projections::m_signature = ('P'*256 + 'J');
 
@@ -1016,36 +1016,36 @@ Projections::convertFFTPolar (ImageFile& rIF, int iInterpolationID, int iZeropad
   if (m_geometry == Scanner::GEOMETRY_EQUIANGULAR || m_geometry == Scanner::GEOMETRY_EQUILINEAR)
     pProj = interpolateToParallel();
 
-  int iInterpDet = nx;
-//  int iInterpDet = pProj->m_nDet;
+  int iInterpDet = static_cast<int>(static_cast<double>(sqrt(nx*nx+ny*ny)));
   int iNumInterpDetWithZeros = ProcessSignal::addZeropadFactor (iInterpDet, iZeropad);
-
+  double dProjScale = iInterpDet / (pProj->viewDiameter() * 0.5);
   double dZeropadRatio = static_cast<double>(iNumInterpDetWithZeros) / static_cast<double>(iInterpDet);
 
   fftw_plan plan = fftw_create_plan (iNumInterpDetWithZeros, FFTW_FORWARD, FFTW_IN_PLACE | FFTW_ESTIMATE | FFTW_USE_WISDOM);
 
   fftw_complex* pcIn = new fftw_complex [iNumInterpDetWithZeros];
   std::complex<double>** ppcDetValue = new std::complex<double>* [pProj->m_nView];
-  double dInterpScale = (pProj->m_nDet-1) / static_cast<double>(iInterpDet-1) / SQRT2;
+  //double dInterpScale = (pProj->m_nDet-1) / static_cast<double>(iInterpDet-1);
+  double dInterpScale = pProj->m_nDet / static_cast<double>(iInterpDet);
   
   double dFFTScale = 1. / static_cast<double>(iInterpDet * iInterpDet);
   int iMidPoint = iInterpDet / 2;
   double dMidPoint = static_cast<double>(iInterpDet) / 2.;
   int iZerosAdded = iNumInterpDetWithZeros - iInterpDet;
 
-  // For each view, interpolate to nx length, shift to center at origin, and FFt transform
+  // For each view, interpolate, shift to center at origin, and FFT
   for (int iView = 0; iView < m_nView; iView++) {
     DetectorValue* detval = pProj->getDetectorArray(iView).detValues();
     LinearInterpolator<DetectorValue> projInterp (detval, pProj->m_nDet);
     for (int iDet = 0; iDet < iInterpDet; iDet++) {
       double dInterpPos = (m_nDet / 2.) + (iDet - dMidPoint) * dInterpScale;
-      pcIn[iDet].re = projInterp.interpolate (dInterpPos) * dInterpScale;
+      pcIn[iDet].re = projInterp.interpolate (dInterpPos) * dProjScale;
       pcIn[iDet].im = 0;
     }
 
     Fourier::shuffleFourierToNaturalOrder (pcIn, iInterpDet);
     if (iZerosAdded > 0) {
-      for (int iDet1 = iMidPoint; iDet1 < iInterpDet; iDet1++)
+       for (int iDet1 = iInterpDet -1; iDet1 >= iMidPoint; iDet1--)
         pcIn[iDet1+iZerosAdded] = pcIn[iDet1];
       for (int iDet2 = iMidPoint; iDet2 < iMidPoint + iZerosAdded; iDet2++) 
         pcIn[iDet2].re = pcIn[iDet2].im = 0;
@@ -1107,10 +1107,10 @@ Projections::calcArrayPolarCoordinates (unsigned int nx, unsigned int ny, double
   double xInc = (xMax - xMin) / nx;    // size of cells
   double yInc = (yMax - yMin) / ny;
 
+  double dDetCenter = (iNumDetWithZeros - 1) / 2.;     // index refering to L=0 projection 
   // +1 is correct for frequency data, ndet-1 is correct for projections
-  int iDetCenter = (iNumDetWithZeros - 1) / 2; // index refering to L=0 projection 
-  if (isEven (iNumDetWithZeros))
-    iDetCenter = (iNumDetWithZeros + 1) / 2;   
+  //  if (isEven (iNumDetWithZeros))
+  //    dDetCenter = (iNumDetWithZeros + 0) / 2;       
 
   // Calculates polar coordinates (view#, det#) for each point on phantom grid
   double x = xMin + xInc / 2;  // Rectang coords of center of pixel 
@@ -1128,7 +1128,7 @@ Projections::calcArrayPolarCoordinates (unsigned int nx, unsigned int ny, double
       }
       
       ppdView[ix][iy] = (phi - m_rotStart) / m_rotInc;
-      ppdDet[ix][iy] = (r / dDetInc) + iDetCenter;
+      ppdDet[ix][iy] = (r / dDetInc) + dDetCenter;
     }
   }
 }
index 79264ffbd52f1397f7f1f08e24178e2786439f5d..77145625d532cc83585c367dd08d6c040e1fbf4a 100644 (file)
@@ -9,7 +9,7 @@
 **  This is part of the CTSim program
 **  Copyright (c) 1983-2001 Kevin Rosenberg
 **
-**  $Id: dialogs.h,v 1.38 2001/09/24 11:46:17 kevin Exp $
+**  $Id: dialogs.h,v 1.39 2003/03/15 10:27:30 kevin Exp $
 **
 **  This program is free software; you can redistribute it and/or modify
 **  it under the terms of the GNU General Public License (version 2) as
@@ -393,7 +393,7 @@ class DialogGetConvertPolarParameters : public wxDialog
  public:
    DialogGetConvertPolarParameters (wxWindow* pParent, const char* const pszTitle, int iDefaultXSize = 0, 
      int iDefaultYSize = 0, int iDefaultInterpolationID = Projections::POLAR_INTERP_BILINEAR, 
-     int iDefaultZeropad = 1, int iHelpID = IDH_DLG_POLAR);
+     int iDefaultZeropad = 3, int iHelpID = IDH_DLG_POLAR);
    virtual ~DialogGetConvertPolarParameters ();
 
     unsigned int getXSize();
index e3922dafce711ee8e7c1c44a01ed10e76aa8bb50..cb7663780813ec236fda07c4c46d4a647c7e58c4 100644 (file)
@@ -9,7 +9,7 @@
 **  This is part of the CTSim program
 **  Copyright (c) 1983-2001 Kevin Rosenberg
 **
-**  $Id: views.cpp,v 1.170 2003/01/30 21:53:16 kevin Exp $
+**  $Id: views.cpp,v 1.171 2003/03/15 10:27:30 kevin Exp $
 **
 **  This program is free software; you can redistribute it and/or modify
 **  it under the terms of the GNU General Public License (version 2) as
@@ -2521,7 +2521,7 @@ ProjectionFileView::ProjectionFileView()
   m_iDefaultFilterMethod = ProcessSignal::FILTER_METHOD_CONVOLUTION;
   m_iDefaultFilterGeneration = ProcessSignal::FILTER_GENERATION_DIRECT;
 #endif
-  m_iDefaultZeropad = 1;
+  m_iDefaultZeropad = 2;
   m_iDefaultBackprojector = Backprojector::BPROJ_IDIFF;
   m_iDefaultInterpolation = Backprojector::INTERP_LINEAR;
   m_iDefaultInterpParam = 1;
@@ -2530,7 +2530,7 @@ ProjectionFileView::ProjectionFileView()
   m_iDefaultPolarNX = 256;
   m_iDefaultPolarNY = 256;
   m_iDefaultPolarInterpolation = Projections::POLAR_INTERP_BILINEAR;
-  m_iDefaultPolarZeropad = 1;
+  m_iDefaultPolarZeropad = 2;
 }
 
 ProjectionFileView::~ProjectionFileView()
@@ -3063,12 +3063,12 @@ ProjectionFileView::CreateChildFrame(wxDocument *doc, wxView *view)
   GetDocumentManager()->FileHistoryAddFilesToMenu(m_pFileMenu);
   GetDocumentManager()->FileHistoryUseMenu(m_pFileMenu);
   
-  wxMenu *convert_menu = new wxMenu;
-  convert_menu->Append (PJMENU_CONVERT_RECTANGULAR, "&Rectangular Image");
-  convert_menu->Append (PJMENU_CONVERT_POLAR, "&Polar Image...\tCtrl-L");
-  convert_menu->Append (PJMENU_CONVERT_FFT_POLAR, "FF&T->Polar Image...\tCtrl-T");
-  convert_menu->AppendSeparator();
-  convert_menu->Append (PJMENU_CONVERT_PARALLEL, "&Interpolate to Parallel");
+  m_pConvertMenu = new wxMenu;
+  m_pConvertMenu->Append (PJMENU_CONVERT_RECTANGULAR, "&Rectangular Image");
+  m_pConvertMenu->Append (PJMENU_CONVERT_POLAR, "&Polar Image...\tCtrl-L");
+  m_pConvertMenu->Append (PJMENU_CONVERT_FFT_POLAR, "FF&T->Polar Image...\tCtrl-T");
+  m_pConvertMenu->AppendSeparator();
+  m_pConvertMenu->Append (PJMENU_CONVERT_PARALLEL, "&Interpolate to Parallel");
   
   //  wxMenu* filter_menu = new wxMenu;
   //  filter_menu->Append (PJMENU_ARTIFACT_REDUCTION, "&Artifact Reduction");
@@ -3077,11 +3077,10 @@ ProjectionFileView::CreateChildFrame(wxDocument *doc, wxView *view)
   analyze_menu->Append (PJMENU_PLOT_HISTOGRAM, "&Plot Histogram");
   analyze_menu->Append (PJMENU_PLOT_TTHETA_SAMPLING, "Plot T-T&heta Sampling...\tCtrl-H");
 
-  wxMenu *reconstruct_menu = new wxMenu;
-  reconstruct_menu->Append (PJMENU_RECONSTRUCT_FBP, "&Filtered Backprojection...\tCtrl-R", "Reconstruct image using filtered backprojection");
-  reconstruct_menu->Append (PJMENU_RECONSTRUCT_FBP_REBIN, "Filtered &Backprojection (Rebin to Parallel)...\tCtrl-B", "Reconstruct image using filtered backprojection");
-  // still buggy
-  //   reconstruct_menu->Append (PJMENU_RECONSTRUCT_FOURIER, "&Fourier...\tCtrl-E", "Reconstruct image using inverse Fourier");
+  m_pReconstructMenu = new wxMenu;
+  m_pReconstructMenu->Append (PJMENU_RECONSTRUCT_FBP, "&Filtered Backprojection...\tCtrl-R", "Reconstruct image using filtered backprojection");
+  m_pReconstructMenu->Append (PJMENU_RECONSTRUCT_FBP_REBIN, "Filtered &Backprojection (Rebin to Parallel)...\tCtrl-B", "Reconstruct image using filtered backprojection");
+  m_pReconstructMenu->Append (PJMENU_RECONSTRUCT_FOURIER, "&Inverse Fourier...\tCtrl-E", "Direct inverse Fourier");
   
   wxMenu *help_menu = new wxMenu;
   help_menu->Append(MAINMENU_HELP_CONTENTS, "&Contents\tF1");
@@ -3092,10 +3091,10 @@ ProjectionFileView::CreateChildFrame(wxDocument *doc, wxView *view)
   wxMenuBar *menu_bar = new wxMenuBar;
   
   menu_bar->Append (m_pFileMenu, "&File");
-  menu_bar->Append (convert_menu, "&Convert");
+  menu_bar->Append (m_pConvertMenu, "&Convert");
   //  menu_bar->Append (filter_menu, "Fi&lter");
   menu_bar->Append (analyze_menu, "&Analyze");
-  menu_bar->Append (reconstruct_menu, "&Reconstruct");
+  menu_bar->Append (m_pReconstructMenu, "&Reconstruct");
   menu_bar->Append (help_menu, "&Help");
   
   subframe->SetMenuBar(menu_bar);  
@@ -3158,6 +3157,14 @@ ProjectionFileView::OnUpdate (wxView *WXUNUSED(sender), wxObject *WXUNUSED(hint)
   const Projections& rProj = GetDocument()->getProjections();
   const int nDet = rProj.nDet();
   const int nView = rProj.nView();
+  if (rProj.geometry() == Scanner::GEOMETRY_PARALLEL) { 
+    m_pReconstructMenu->Enable (PJMENU_RECONSTRUCT_FBP_REBIN, false);
+    m_pConvertMenu->Enable (PJMENU_CONVERT_PARALLEL, false);
+  } else {
+    m_pReconstructMenu->Enable (PJMENU_RECONSTRUCT_FBP_REBIN, true);
+    m_pConvertMenu->Enable (PJMENU_CONVERT_PARALLEL, true);
+  }
+
   if (nDet != 0 && nView != 0) {
     const DetectorArray& detarray = rProj.getDetectorArray(0);
     const DetectorValue* detval = detarray.detValues();
index 42394021f3cc270dd7ef73c76b3c152ae49804f3..0dc47c91f55f241147f44285c4a708df006a3a15 100644 (file)
@@ -9,7 +9,7 @@
 **  This is part of the CTSim program
 **  Copyright (c) 1983-2001 Kevin Rosenberg
 **
-**  $Id: views.h,v 1.60 2003/01/24 05:24:19 kevin Exp $
+**  $Id: views.h,v 1.61 2003/03/15 10:27:30 kevin Exp $
 **
 **  This program is free software; you can redistribute it and/or modify
 **  it under the terms of the GNU General Public License (version 2) as
@@ -204,6 +204,8 @@ private:
   
   ProjectionFileCanvas *m_pCanvas;
   wxMenu* m_pFileMenu;
+  wxMenu* m_pReconstructMenu;
+  wxMenu* m_pConvertMenu;
 
   int m_iDefaultNX;
   int m_iDefaultNY;
@@ -262,6 +264,7 @@ public:
   void setInitialClientSize();
 
   wxMenu* getFileMenu()  { return m_pFileMenu; }
+  wxMenu* getReconstructMenu()  { return m_pReconstructMenu; }
 
   ProjectionFileDocument* GetDocument() 
   { return dynamic_cast<ProjectionFileDocument*>(wxView::GetDocument()); }