r641: no message
[ctsim.git] / libctsim / imagefile.cpp
index 946e4db75e17ee766a924939f2237b02a13ee1c7..bab87bdeff5ab0fa4486a66f629089e85efd4d91 100644 (file)
@@ -1,15 +1,15 @@
 /*****************************************************************************
 ** FILE IDENTIFICATION
 **
-**     Name:         imagefile.cpp
-**      Purpose:      Imagefile classes
+**     Name:           imagefile.cpp
+**  Purpose:      Imagefile classes
 **     Programmer:   Kevin Rosenberg
 **     Date Started: June 2000
 **
 **  This is part of the CTSim program
-**  Copyright (C) 1983-2000 Kevin Rosenberg
+**  Copyright (c) 1983-2001 Kevin Rosenberg
 **
-**  $Id: imagefile.cpp,v 1.34 2001/01/04 21:28:41 kevin Exp $
+**  $Id: imagefile.cpp,v 1.42 2001/03/21 21:45:31 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
 ******************************************************************************/
 
 #include "ct.h"
+#ifdef HAVE_CTN_DICOM
+#include "ctndicom.h"
+#endif
+
+const double ImageFile::s_dRedGrayscaleFactor = 0.299;
+const double ImageFile::s_dGreenGrayscaleFactor = 0.587;
+const double ImageFile::s_dBlueGrayscaleFactor = 0.114;
 
-const int ImageFile::FORMAT_INVALID = -1;
-const int ImageFile::FORMAT_PGM = 0;
-const int ImageFile::FORMAT_PGMASCII = 1;
+
+const int ImageFile::EXPORT_FORMAT_INVALID = -1;
+const int ImageFile::EXPORT_FORMAT_TEXT = 0;
+const int ImageFile::EXPORT_FORMAT_PGM = 1;
+const int ImageFile::EXPORT_FORMAT_PGMASCII = 2;
 #ifdef HAVE_PNG
-const int ImageFile::FORMAT_PNG = 2;
-const int ImageFile::FORMAT_PNG16 = 3;
+const int ImageFile::EXPORT_FORMAT_PNG = 3;
+const int ImageFile::EXPORT_FORMAT_PNG16 = 4;
+#endif
+#ifdef HAVE_CTN_DICOM
+const int ImageFile::EXPORT_FORMAT_DICOM = 5;
 #endif
 
-const char* ImageFile::s_aszFormatName[] = 
+const char* ImageFile::s_aszExportFormatName[] = 
 {
+  {"text"},
   {"pgm"},
   {"pgmascii"},
 #ifdef HAVE_PNG
   {"png"},
   {"png16"},
 #endif
+#ifdef HAVE_CTN_DICOM
+  {"dicom"},
+#endif
 };
 
-const char* ImageFile::s_aszFormatTitle[] = 
+const char* ImageFile::s_aszExportFormatTitle[] = 
 {
+  {"Text"},
   {"PGM"},
   {"PGM ASCII"},
   {"PNG"},
   {"PNG 16-bit"},
+#ifdef HAVE_CTN_DICOM
+  {"Dicom"},
+#endif
+};
+const int ImageFile::s_iExportFormatCount = sizeof(s_aszExportFormatName) / sizeof(const char*);
+
+
+const int ImageFile::IMPORT_FORMAT_INVALID = -1;
+const int ImageFile::IMPORT_FORMAT_PPM = 0;
+#ifdef HAVE_PNG
+const int ImageFile::IMPORT_FORMAT_PNG = 1;
+#endif
+#ifdef HAVE_CTN_DICOM
+const int ImageFile::IMPORT_FORMAT_DICOM = 2;
+#endif
+
+
+const char* ImageFile::s_aszImportFormatName[] = 
+{
+  {"ppm"},
+#ifdef HAVE_PNG
+  {"png"},
+#endif
+#ifdef HAVE_CTN_DICOM
+  {"dicom"},
+#endif
 };
 
-const int ImageFile::s_iFormatCount = sizeof(s_aszFormatName) / sizeof(const char*);
+const char* ImageFile::s_aszImportFormatTitle[] = 
+{
+  {"PPM"},
+  {"PNG"},
+#ifdef HAVE_CTN_DICOM
+  {"Dicom"},
+#endif
+};
+const int ImageFile::s_iImportFormatCount = sizeof(s_aszImportFormatName) / sizeof(const char*);
 
 
 
@@ -90,7 +141,7 @@ ImageFile::getCenterCoordinates (unsigned int& iXCenter, unsigned int& iYCenter)
     iXCenter = m_nx / 2;
   else
     iXCenter = (m_nx - 1) / 2;
-
+  
   if (isEven (m_ny))
     iYCenter = m_ny / 2;
   else
@@ -99,14 +150,15 @@ ImageFile::getCenterCoordinates (unsigned int& iXCenter, unsigned int& iYCenter)
 
 
 void 
-ImageFile::filterResponse (const char* const domainName, double bw, const char* const filterName, double filt_param, double dInputScale, double dOutputScale)
+ImageFile::filterResponse (const char* const domainName, double bw, const char* const filterName, 
+                           double filt_param, double dInputScale, double dOutputScale)
 {
   ImageFileArray v = getArray();
   SignalFilter filter (filterName, domainName, bw, filt_param);
   
   unsigned int iXCenter, iYCenter;
   getCenterCoordinates (iXCenter, iYCenter);
-
+  
   for (unsigned int ix = 0; ix < m_nx; ix++)
     for (unsigned int iy = 0; iy < m_ny; iy++) {
       long lD2 = ((ix - iXCenter) * (ix - iXCenter)) + ((iy - iYCenter) * (iy - iYCenter));
@@ -115,60 +167,6 @@ ImageFile::filterResponse (const char* const domainName, double bw, const char*
     }
 }
 
-int
-ImageFile::display (void) const
-{
-  double pmin, pmax;
-  
-  getMinMax (pmin, pmax);
-  
-  return (displayScaling (1, pmin, pmax));
-}
-
-int 
-ImageFile::displayScaling (const int scale, const ImageFileValue pmin, const ImageFileValue pmax) const
-{
-  int nx = m_nx;
-  int ny = m_ny;
-  ImageFileArrayConst v = getArray();
-  if (v == NULL || nx == 0 || ny == 0)
-    return 0;
-  
-#if HAVE_G2_H
-  int* pPens = new int [nx * ny * scale * scale ];
-  
-  double view_scale = 255 / (pmax - pmin);
-  int id_X11 = g2_open_X11 (nx * scale, ny * scale);
-  int grayscale[256];
-  for (int i = 0; i < 256; i++) {
-    double cval = i / 255.;
-    grayscale[i] = g2_ink (id_X11, cval, cval, cval);
-  }
-  
-  for (int iy = ny - 1; iy >= 0; iy--) {
-    int iRowPos = ((ny - 1 - iy) * scale) * (nx * scale);
-    for (int ix = 0; ix < nx; ix++) {
-      int cval = static_cast<int>((v[ix][iy] - pmin) * view_scale);
-      if (cval < 0)  
-        cval = 0;
-      else if (cval > 255) 
-        cval = 255;
-      for (int sy = 0; sy < scale; sy++)
-        for (int sx = 0; sx < scale; sx++)
-          pPens[iRowPos+(sy * nx * scale)+(sx + (ix * scale))] = grayscale[cval];
-    }
-  }
-  
-  g2_image (id_X11, 0., 0., nx * scale, ny * scale, pPens);
-  
-  delete pPens;
-  return (id_X11);
-#else
-  return 0;
-#endif
-}
-
-
 
 // ImageFile::comparativeStatistics    Calculate comparative stats
 //
@@ -599,7 +597,7 @@ bool
 ImageFile::log (ImageFile& result) const
 {
   if (m_nx != result.nx() || m_ny != result.ny()) {
-    sys_error (ERR_WARNING, "Difference sizes of images [ImageFile::invertPixelValues]");
+    sys_error (ERR_WARNING, "Difference sizes of images [ImageFile::log]");
     return false;
   }
   
@@ -621,8 +619,12 @@ ImageFile::log (ImageFile& result) const
         std::complex<double> cResult = std::log (cLHS);
         vResult[ix][iy] = cResult.real();
         vResultImag[ix][iy] = cResult.imag();
-      } else
-        vResult[ix][iy] = ::log (vLHS[ix][iy]);
+      } else {
+        if (vLHS[ix][iy] > 0)
+          vResult[ix][iy] = ::log (vLHS[ix][iy]);
+        else
+          vResult[ix][iy] = 0;
+      }
     }
   }
   
@@ -684,38 +686,22 @@ ImageFile::scaleImage (ImageFile& result) const
   ImageFileArray vResult = result.getArray();
   ImageFileArray vResultImag = result.getImaginaryArray();
   
+  BilinearInterpolator<ImageFileValue> realInterp (vReal, nx, ny);
+  BilinearInterpolator<ImageFileValue> imagInterp (vImag, nx, ny);
+
   for (unsigned int ix = 0; ix < newNX; ix++) {
     for (unsigned int iy = 0; iy < newNY; iy++) {
       double dXPos = ix / dXScale;
       double dYPos = iy / dYScale;
-      unsigned int scaleNX = static_cast<unsigned int> (dXPos);
-      unsigned int scaleNY = static_cast<unsigned int> (dYPos);
-      double dXFrac = dXPos - scaleNX;
-      double dYFrac = dYPos - scaleNY;
-      if (scaleNX >= nx - 1 || scaleNY >= ny - 1) {
-        vResult[ix][iy] = vReal[scaleNX][scaleNY];
-        if (result.isComplex()) {
-          if (isComplex())
-            vResultImag[ix][iy] = vImag[scaleNX][scaleNY];
-          else
-            vResultImag[ix][iy] = 0;
-        }
-      } else {
-        vResult[ix][iy] = vReal[scaleNX][scaleNY] + 
-          dXFrac * (vReal[scaleNX+1][scaleNY] - vReal[scaleNX][scaleNY]) + 
-          dYFrac * (vReal[scaleNX][scaleNY+1] - vReal[scaleNX][scaleNY]);
-        if (result.isComplex()) {
-          if (isComplex())
-            vResultImag[ix][iy] = vImag[scaleNX][scaleNY] + 
-            dXFrac * (vImag[scaleNX+1][scaleNY] - vImag[scaleNX][scaleNY]) + 
-            dYFrac * (vImag[scaleNX][scaleNY+1] - vImag[scaleNX][scaleNY]);
-          else
-            vResultImag[ix][iy] = 0;
-        }
-      }
+      vResult[ix][iy] = realInterp.interpolate (dXPos, dYPos);
+      if (result.isComplex())
+        if (isComplex())
+          vResultImag[ix][iy] = imagInterp.interpolate (dXPos, dYPos);
+        else
+          vResultImag[ix][iy] = 0;
     }
   }
-  
+
   return true;
 }
 
@@ -740,7 +726,7 @@ ImageFile::fft (ImageFile& result) const
   
   unsigned int ix, iy;
   unsigned int iArray = 0;
-  for (ix = 0; ix < m_nx; ix++)
+  for (ix = 0; ix < m_nx; ix++) {
     for (iy = 0; iy < m_ny; iy++) {
       in[iArray].re = vReal[ix][iy];
       if (isComplex())
@@ -749,29 +735,28 @@ ImageFile::fft (ImageFile& result) const
         in[iArray].im = 0;
       iArray++;
     }
-    
-    fftwnd_plan plan = fftw2d_create_plan (m_nx, m_ny, FFTW_FORWARD, FFTW_IN_PLACE);
-    
-    fftwnd_one (plan, in, NULL);
-    
-    ImageFileArray vRealResult = result.getArray();
-    ImageFileArray vImagResult = result.getImaginaryArray();
-    iArray = 0;
-    unsigned int iScale = m_nx * m_ny;
-    for (ix = 0; ix < m_nx; ix++)
-      for (iy = 0; iy < m_ny; iy++) {
-        vRealResult[ix][iy] = in[iArray].re / iScale;
-        vImagResult[ix][iy] = in[iArray].im / iScale;
-        iArray++;
-      }
-      
-      fftwnd_destroy_plan (plan);
-      delete in;
-      
-      
-      Fourier::shuffleFourierToNaturalOrder (result);
-      
-      return true;
+  }
+  
+  fftwnd_plan plan = fftw2d_create_plan (m_nx, m_ny, FFTW_FORWARD, FFTW_IN_PLACE | FFTW_ESTIMATE | FFTW_USE_WISDOM);  
+  fftwnd_one (plan, in, NULL);
+  
+  ImageFileArray vRealResult = result.getArray();
+  ImageFileArray vImagResult = result.getImaginaryArray();
+  iArray = 0;
+  unsigned int iScale = m_nx * m_ny;
+  for (ix = 0; ix < m_nx; ix++) {
+    for (iy = 0; iy < m_ny; iy++) {
+      vRealResult[ix][iy] = in[iArray].re / iScale;
+      vImagResult[ix][iy] = in[iArray].im / iScale;
+      iArray++;
+    }
+  }
+  delete in;
+  fftwnd_destroy_plan (plan);
+  
+  Fourier::shuffleFourierToNaturalOrder (result);
+  
+  return true;
 }
 
 
@@ -793,7 +778,7 @@ ImageFile::ifft (ImageFile& result) const
   ImageFileArray vRealResult = result.getArray();
   ImageFileArray vImagResult = result.getImaginaryArray();
   unsigned int ix, iy;
-  for (ix = 0; ix < m_nx; ix++)
+  for (ix = 0; ix < m_nx; ix++) {
     for (iy = 0; iy < m_ny; iy++) {
       vRealResult[ix][iy] = vReal[ix][iy];
       if (isComplex()) 
@@ -801,37 +786,244 @@ ImageFile::ifft (ImageFile& result) const
       else
         vImagResult[ix][iy] = 0;
     }
+  }
+  
+  Fourier::shuffleNaturalToFourierOrder (result);
+  
+  fftw_complex* in = new fftw_complex [m_nx * m_ny];
+  
+  unsigned int iArray = 0;
+  for (ix = 0; ix < m_nx; ix++) {
+    for (iy = 0; iy < m_ny; iy++) {
+      in[iArray].re = vRealResult[ix][iy];
+      in[iArray].im = vImagResult[ix][iy];
+      iArray++;
+    }
+  }
+  
+  fftwnd_plan plan = fftw2d_create_plan (m_nx, m_ny, FFTW_BACKWARD, FFTW_IN_PLACE | FFTW_ESTIMATE | FFTW_USE_WISDOM);
+  
+  fftwnd_one (plan, in, NULL);
+  
+  iArray = 0;
+  for (ix = 0; ix < m_nx; ix++) {
+    for (iy = 0; iy < m_ny; iy++) {
+      vRealResult[ix][iy] = in[iArray].re;
+      vImagResult[ix][iy] = in[iArray].im;
+      iArray++;
+    }
+  }
+  fftwnd_destroy_plan (plan);
+  
+  delete in;
+  
+  return true;
+}
+
+bool
+ImageFile::fftRows (ImageFile& result) const
+{
+  if (m_nx != result.nx() || m_ny != result.ny()) {
+    sys_error (ERR_WARNING, "Difference sizes of images [ImageFile::fftRows]");
+    return false;
+  }
+  
+  if (result.dataType() == Array2dFile::DATA_TYPE_REAL) {
+    if (! result.convertRealToComplex ())
+      return false;
+  }
+   
+  ImageFileArrayConst vReal = getArray();
+  ImageFileArrayConst vImag = getImaginaryArray();
+  
+  fftw_plan plan = fftw_create_plan (m_nx, FFTW_FORWARD, FFTW_IN_PLACE | FFTW_ESTIMATE | FFTW_USE_WISDOM);
+
+  fftw_complex* in = new fftw_complex [m_nx];
+  std::complex<double>* pcRow = new std::complex<double> [m_nx];  
+  for (int iy = 0; iy < m_ny; iy++) {
+    unsigned int ix;
+    for (ix = 0; ix < m_nx; ix++) {
+      in[ix].re = vReal[ix][iy];
+      if (isComplex())
+        in[ix].im = vImag[ix][iy];
+      else
+        in[ix].im = 0;
+    }
     
-    Fourier::shuffleNaturalToFourierOrder (result);
-    
-    fftw_complex* in = new fftw_complex [m_nx * m_ny];
+    fftw_one (plan, in, NULL);
     
-    unsigned int iArray = 0;
     for (ix = 0; ix < m_nx; ix++)
-      for (iy = 0; iy < m_ny; iy++) {
-        in[iArray].re = vRealResult[ix][iy];
-        in[iArray].im = vImagResult[ix][iy];
-        iArray++;
-      }
-      
-      fftwnd_plan plan = fftw2d_create_plan (m_nx, m_ny, FFTW_BACKWARD, FFTW_IN_PLACE);
-      
-      fftwnd_one (plan, in, NULL);
-      
-      iArray = 0;
-      for (ix = 0; ix < m_nx; ix++)
-        for (iy = 0; iy < m_ny; iy++) {
-          vRealResult[ix][iy] = in[iArray].re;
-          vImagResult[ix][iy] = in[iArray].im;
-          iArray++;
-        }
-        
-        fftwnd_destroy_plan (plan);
-        
-        delete in;
-        
-        return true;
+      pcRow[ix] = std::complex<double>(in[ix].re, in[ix].im);
+    
+    Fourier::shuffleFourierToNaturalOrder (pcRow, m_nx);
+    for (ix = 0; ix < m_nx; ix++) {
+      vReal[ix][iy] = pcRow[ix].real() / m_nx;
+      vImag[ix][iy] = pcRow[ix].imag() / m_nx;
+    }
+  }
+  delete [] pcRow;
+  
+  fftw_destroy_plan (plan);
+  delete in;
+  
+  return true;
+}     
+
+bool
+ImageFile::ifftRows (ImageFile& result) const
+{
+  if (m_nx != result.nx() || m_ny != result.ny()) {
+    sys_error (ERR_WARNING, "Difference sizes of images [ImageFile::fftRows]");
+    return false;
+  }
+  
+  if (result.dataType() == Array2dFile::DATA_TYPE_REAL) {
+    if (! result.convertRealToComplex ())
+      return false;
+  }
+  
+  fftw_complex* in = new fftw_complex [m_nx];
+  
+  ImageFileArrayConst vReal = getArray();
+  ImageFileArrayConst vImag = getImaginaryArray();
+  
+  fftw_plan plan = fftw_create_plan (m_nx, FFTW_BACKWARD, FFTW_IN_PLACE | FFTW_ESTIMATE | FFTW_USE_WISDOM);
+  std::complex<double>* pcRow = new std::complex<double> [m_nx];
+  
+  unsigned int ix, iy;
+  unsigned int iArray = 0;
+  for (iy = 0; iy < m_ny; iy++) {
+    for (ix = 0; ix < m_nx; ix++) {
+      double dImag = 0;
+      if (isComplex())
+        dImag = vImag[ix][iy];
+      pcRow[ix] = std::complex<double> (vReal[ix][iy], dImag);
+    }
+    
+    Fourier::shuffleNaturalToFourierOrder (pcRow, m_nx);
+    
+    for (ix = 0; ix < m_nx; ix++) {
+      in[ix].re = pcRow[ix].real();
+      in[ix].im = pcRow[ix].imag();
+    }
+    
+    fftw_one (plan, in, NULL);
+    
+    for (ix = 0; ix < m_nx; ix++) {
+      vReal[ix][iy] = in[ix].re;
+      vImag[ix][iy] = in[ix].im;
+    }
+  }
+  delete [] pcRow;
+  
+  fftw_destroy_plan (plan);
+  delete in;
+  
+  return true;
+}
+
+bool
+ImageFile::fftCols (ImageFile& result) const
+{
+  if (m_nx != result.nx() || m_ny != result.ny()) {
+    sys_error (ERR_WARNING, "Difference sizes of images [ImageFile::fftRows]");
+    return false;
+  }
+  
+  if (result.dataType() == Array2dFile::DATA_TYPE_REAL) {
+    if (! result.convertRealToComplex ())
+      return false;
+  }
+    
+  ImageFileArrayConst vReal = getArray();
+  ImageFileArrayConst vImag = getImaginaryArray();
+  
+  fftw_plan plan = fftw_create_plan (m_ny, FFTW_FORWARD, FFTW_IN_PLACE | FFTW_ESTIMATE | FFTW_USE_WISDOM);
+
+  std::complex<double>* pcCol = new std::complex<double> [m_ny];  
+  fftw_complex* in = new fftw_complex [m_ny];
+  for (int ix = 0; ix < m_nx; ix++) {
+    unsigned int iy;
+    for (iy = 0; iy < m_ny; iy++) {
+      in[iy].re = vReal[ix][iy];
+      if (isComplex())
+        in[iy].im = vImag[ix][iy];
+      else
+        in[iy].im = 0;
+    }
+    
+    fftw_one (plan, in, NULL);
+    
+    for (iy = 0; iy < m_ny; iy++)
+      pcCol[iy] = std::complex<double>(in[iy].re, in[iy].im);
+    
+    Fourier::shuffleFourierToNaturalOrder (pcCol, m_ny);
+    for (iy = 0; iy < m_ny; iy++) {
+      vReal[ix][iy] = pcCol[iy].real() / m_ny;
+      vImag[ix][iy] = pcCol[iy].imag() / m_ny;
+    }
+  }
+  delete [] pcCol;
+  
+  fftw_destroy_plan (plan);
+  delete in;
+  
+  return true;
+}
+
+bool
+ImageFile::ifftCols (ImageFile& result) const
+{
+  if (m_nx != result.nx() || m_ny != result.ny()) {
+    sys_error (ERR_WARNING, "Difference sizes of images [ImageFile::fftRows]");
+    return false;
+  }
+  
+  if (result.dataType() == Array2dFile::DATA_TYPE_REAL) {
+    if (! result.convertRealToComplex ())
+      return false;
+  }
+  
+  fftw_complex* in = new fftw_complex [m_ny];
+  
+  ImageFileArrayConst vReal = getArray();
+  ImageFileArrayConst vImag = getImaginaryArray();
+  
+  fftw_plan plan = fftw_create_plan (m_ny, FFTW_BACKWARD, FFTW_IN_PLACE | FFTW_ESTIMATE | FFTW_USE_WISDOM);
+  std::complex<double>* pcCol = new std::complex<double> [m_ny];
+  
+  unsigned int ix, iy;
+  unsigned int iArray = 0;
+  for (ix = 0; ix < m_nx; ix++) {
+    for (iy = 0; iy < m_ny; iy++) {
+      double dImag = 0;
+      if (isComplex())
+        dImag = vImag[ix][iy];
+      pcCol[iy] = std::complex<double> (vReal[ix][iy], dImag);
+    }
+    
+    Fourier::shuffleNaturalToFourierOrder (pcCol, m_ny);
+    
+    for (iy = 0; iy < m_ny; iy++) {
+      in[iy].re = pcCol[iy].real();
+      in[iy].im = pcCol[iy].imag();
+    }
+    
+    fftw_one (plan, in, NULL);
+    
+    for (iy = 0; iy < m_ny; iy++) {
+      vReal[ix][iy] = in[iy].re;
+      vImag[ix][iy] = in[iy].im;
+    }
+  }
+  delete [] pcCol;
+  
+  fftw_destroy_plan (plan);
+  delete in;
+  
+  return true;
 }
+
 #endif // HAVE_FFTW
 
 
@@ -997,7 +1189,7 @@ ImageFile::magnitude (ImageFile& result) const
     }
     
     if (result.isComplex())
-      result.convertComplexToReal();
+      result.reallocComplexToReal();
     
     return true;
 }
@@ -1014,18 +1206,67 @@ ImageFile::phase (ImageFile& result) const
   ImageFileArray vImag = getImaginaryArray();
   ImageFileArray vRealResult = result.getArray();
   
-  for (unsigned int ix = 0; ix < m_nx; ix++)
+  for (unsigned int ix = 0; ix < m_nx; ix++) {
     for (unsigned int iy = 0; iy < m_ny; iy++) {
       if (isComplex())
         vRealResult[ix][iy] = ::atan2 (vImag[ix][iy], vReal[ix][iy]);
       else
         vRealResult[ix][iy] = 0;
     }
+  }
+  if (result.isComplex())
+    result.reallocComplexToReal();
     
-    if (result.isComplex())
-      result.convertComplexToReal();
+  return true;
+}
+
+bool
+ImageFile::real (ImageFile& result) const
+{
+  if (m_nx != result.nx() || m_ny != result.ny()) {
+    sys_error (ERR_WARNING, "Difference sizes of images [ImageFile::invertPixelValues]");
+    return false;
+  }
+  
+  ImageFileArray vReal = getArray();
+  ImageFileArray vRealResult = result.getArray();
+  
+  for (unsigned int ix = 0; ix < m_nx; ix++) {
+    for (unsigned int iy = 0; iy < m_ny; iy++) {
+        vRealResult[ix][iy] = vReal[ix][iy];
+    }
+  }
+
+  if (result.isComplex())
+    result.reallocComplexToReal();
     
-    return true;
+  return true;
+}
+
+bool
+ImageFile::imaginary (ImageFile& result) const
+{
+  if (m_nx != result.nx() || m_ny != result.ny()) {
+    sys_error (ERR_WARNING, "Difference sizes of images [ImageFile::invertPixelValues]");
+    return false;
+  }
+  
+  ImageFileArray vImag = getArray();
+  ImageFileArray vRealResult = result.getArray();
+  
+  for (unsigned int ix = 0; ix < m_nx; ix++) {
+    for (unsigned int iy = 0; iy < m_ny; iy++) {
+      if (isComplex())
+        vRealResult[ix][iy] = vImag[ix][iy];
+      else
+        vRealResult[ix][iy] = 0;
+    }
+  }
+
+  if (result.isComplex())
+      result.reallocComplexToReal();
+    
+  return true;
 }
 
 bool
@@ -1058,18 +1299,53 @@ ImageFile::square (ImageFile& result) const
         vResult[ix][iy] = vLHS[ix][iy] * vLHS[ix][iy];
     }
   }
+    
+  return true;
+}
+
+int
+ImageFile::convertExportFormatNameToID (const char* const formatName)
+{
+  int formatID = EXPORT_FORMAT_INVALID;
   
+  for (int i = 0; i < s_iExportFormatCount; i++)
+    if (strcasecmp (formatName, s_aszExportFormatName[i]) == 0) {
+      formatID = i;
+      break;
+    }
+    
+    return (formatID);
+}
+
+const char*
+ImageFile::convertExportFormatIDToName (int formatID)
+{
+  static const char *formatName = "";
   
-  return true;
+  if (formatID >= 0 && formatID < s_iExportFormatCount)
+    return (s_aszExportFormatName[formatID]);
+  
+  return (formatName);
+}
+
+const char*
+ImageFile::convertExportFormatIDToTitle (const int formatID)
+{
+  static const char *formatTitle = "";
+  
+  if (formatID >= 0 && formatID < s_iExportFormatCount)
+    return (s_aszExportFormatTitle[formatID]);
+  
+  return (formatTitle);
 }
 
 int
-ImageFile::convertFormatNameToID (const char* const formatName)
+ImageFile::convertImportFormatNameToID (const char* const formatName)
 {
-  int formatID = FORMAT_INVALID;
+  int formatID = IMPORT_FORMAT_INVALID;
   
-  for (int i = 0; i < s_iFormatCount; i++)
-    if (strcasecmp (formatName, s_aszFormatName[i]) == 0) {
+  for (int i = 0; i < s_iImportFormatCount; i++)
+    if (strcasecmp (formatName, s_aszImportFormatName[i]) == 0) {
       formatID = i;
       break;
     }
@@ -1078,49 +1354,325 @@ ImageFile::convertFormatNameToID (const char* const formatName)
 }
 
 const char*
-ImageFile::convertFormatIDToName (int formatID)
+ImageFile::convertImportFormatIDToName (int formatID)
 {
   static const char *formatName = "";
   
-  if (formatID >= 0 && formatID < s_iFormatCount)
-    return (s_aszFormatName[formatID]);
+  if (formatID >= 0 && formatID < s_iImportFormatCount)
+    return (s_aszImportFormatName[formatID]);
   
   return (formatName);
 }
 
 const char*
-ImageFile::convertFormatIDToTitle (const int formatID)
+ImageFile::convertImportFormatIDToTitle (const int formatID)
 {
   static const char *formatTitle = "";
   
-  if (formatID >= 0 && formatID < s_iFormatCount)
-    return (s_aszFormatTitle[formatID]);
+  if (formatID >= 0 && formatID < s_iImportFormatCount)
+    return (s_aszImportFormatTitle[formatID]);
   
   return (formatTitle);
 }
 
 bool
-ImageFile::exportImage (const char* const pszFormat, const char* const pszFilename, int nxcell, int nycell, double densmin, double densmax)
+ImageFile::importImage (const char* const pszFormat, const char* const pszFilename)
+{
+  int iFormatID = convertImportFormatNameToID (pszFormat);
+  
+  if (iFormatID == IMPORT_FORMAT_PPM)
+    return readImagePPM (pszFilename);
+#ifdef HAVE_PNG
+  else if (iFormatID == IMPORT_FORMAT_PNG)
+    return readImagePNG (pszFilename);
+#endif
+  
+  sys_error (ERR_SEVERE, "Invalid format %s [ImageFile::importImage]", pszFormat);
+  return false;
+}
+
+void
+ImageFile::skipSpacePPM (FILE* fp)
+{
+  int c = fgetc (fp);
+  while (isspace (c) || c == '#') {
+    if (c == '#') {   // comment until end of line
+      c = fgetc(fp);
+      while (c != 13 && c != 10)
+        c = fgetc(fp);
+    }
+    else
+      c = fgetc(fp);
+  }
+  
+  ungetc (c, fp);
+}
+
+bool
+ImageFile::readImagePPM (const char* const pszFile)
+{
+  FILE* fp = fopen (pszFile, "r");
+  if ((fp = fopen (pszFile, "r")) == NULL)
+    return false;
+  char cSignature = toupper(fgetc(fp));
+  if (cSignature != 'P') {
+    fclose(fp);
+    return false;
+  }
+  cSignature = fgetc(fp);
+  if (cSignature == '5' || cSignature == '6') { // binary modes
+    fclose(fp);
+    fp = fopen(pszFile, "rb"); // reopen in binary mode
+    fgetc(fp);
+    fgetc(fp);
+  } else if (cSignature != '2' && cSignature != '3') {
+    fclose(fp);
+    return false;
+  }
+  
+  int nRows, nCols, iMaxValue;
+  skipSpacePPM (fp); 
+  if (fscanf (fp, "%d", &nCols) != 1) {
+    fclose(fp);
+    return false;
+  }
+  skipSpacePPM (fp);
+  if (fscanf (fp, "%d", &nRows) != 1) {
+    fclose(fp);
+    return false;
+  }
+  skipSpacePPM (fp);
+  if (fscanf (fp, "%d", &iMaxValue) != 1) {
+    fclose(fp);
+    return false;
+  }
+  setArraySize (nRows, nCols);
+  
+  if (cSignature == '5' || cSignature == '6') { // binary modes
+    int c = fgetc(fp);
+    if (c == 13) {
+      c = fgetc(fp);
+      if (c != 10)  // read msdos 13-10 newline
+        ungetc(c, fp);
+    }
+  } else
+    skipSpacePPM (fp); // ascii may have comments
+  
+  bool bMonochromeImage = false;
+  double dMaxValue = iMaxValue;
+  double dMaxValue3 = iMaxValue * 3;
+
+  ImageFileArray v = getArray();
+  for (int iy = nRows - 1; iy >= 0; iy--) {
+    for (int ix = 0; ix < nCols; ix++) {
+      int iGS, iR, iG, iB;
+      double dR, dG, dB;
+      switch (cSignature) {
+      case '2':
+        if (fscanf(fp, "%d ", &iGS) != 1) {
+          fclose(fp);
+          return false;
+        }
+        v[ix][iy] = iGS / dMaxValue;
+        break;
+      case '5':
+        iGS = fgetc(fp);
+        if (iGS == EOF) {
+          fclose(fp);
+          return false;
+        }
+        v[ix][iy] = iGS / dMaxValue;
+        break;
+      case '3':
+        if (fscanf (fp, "%d %d %d ", &iR, &iG, &iB) != 3) {
+          fclose(fp);
+          return false;
+        }
+        if (ix == 0 && iy == 0 && (iR == iG && iG == iB))
+          bMonochromeImage = true;
+        if (bMonochromeImage)
+          v[ix][iy] = (iR + iG + iB) / dMaxValue3;
+        else {
+          dR = iR / dMaxValue;
+          dG = iG / dMaxValue;
+          dB = iB / dMaxValue;
+          v[ix][iy] = colorToGrayscale (dR, dG, dB);
+        }
+        break;
+      case '6':
+        iR = fgetc(fp);
+        iG = fgetc(fp);
+        iB = fgetc(fp);
+
+        if (iB == EOF) {
+          fclose(fp);
+          return false;
+        }
+        if (ix == 0 && iy == 0 && (iR == iG && iG == iB))
+          bMonochromeImage = true;
+
+        if (bMonochromeImage)
+          v[ix][iy] = (iR + iG + iB) / dMaxValue3;
+        else {
+          dR = iR / dMaxValue;
+          dG = iG / dMaxValue;
+          dB = iB / dMaxValue;
+          v[ix][iy] = colorToGrayscale (dR, dG, dB);
+        }
+        break;
+      }
+    }
+  }
+  
+  fclose(fp);
+  return true;
+}
+
+#ifdef HAVE_PNG
+bool
+ImageFile::readImagePNG (const char* const pszFile)
 {
-  int iFormatID = convertFormatNameToID (pszFormat);
-  if (iFormatID == FORMAT_INVALID) {
-    sys_error (ERR_SEVERE, "Invalid format %s [ImageFile::exportImage]", pszFormat);
+  FILE* fp = fopen(pszFile, "rb");
+  if (!fp) 
+    return false;
+  unsigned char header[8];
+  fread (header, 1, 8, fp);
+  if (png_sig_cmp (header, 0, 8)) {
+    fclose (fp);
+    return false;
+  }
+  
+  png_structp png_ptr = png_create_read_struct (PNG_LIBPNG_VER_STRING, NULL, NULL, NULL);
+  if (!png_ptr) {
+    fclose(fp);
+    return false;
+  }
+  
+  png_infop info_ptr = png_create_info_struct(png_ptr);
+  if (!info_ptr) {
+    png_destroy_read_struct(&png_ptr, (png_infopp)NULL, (png_infopp)NULL);
+    fclose(fp);
+    return false;
+  }
+  
+  png_infop end_info = png_create_info_struct(png_ptr);
+  if (!end_info) {
+    png_destroy_read_struct(&png_ptr, &info_ptr, (png_infopp)NULL);
+    fclose(fp);
+    return false;
+  }
+  
+  if (setjmp(png_ptr->jmpbuf)) {
+    png_destroy_read_struct(&png_ptr, &info_ptr, &end_info);
+    fclose(fp);
     return false;
   }
   
-  if (iFormatID == FORMAT_PGM)
+  png_init_io(png_ptr, fp);
+  png_set_sig_bytes(png_ptr, 8);
+  png_read_info(png_ptr, info_ptr);
+  
+  int width = png_get_image_width (png_ptr, info_ptr);
+  int height = png_get_image_height (png_ptr, info_ptr);
+  int bit_depth = png_get_bit_depth (png_ptr, info_ptr);
+  int color_type = png_get_color_type (png_ptr, info_ptr);
+  
+  if (color_type == PNG_COLOR_TYPE_PALETTE && bit_depth <= 8) 
+    png_set_expand(png_ptr);
+  
+  if (color_type == PNG_COLOR_TYPE_GRAY && bit_depth < 8) 
+    png_set_expand(png_ptr);
+  
+  if (bit_depth < 8)
+    png_set_packing(png_ptr);
+  
+  if (color_type & PNG_COLOR_MASK_ALPHA)
+    png_set_strip_alpha(png_ptr);
+  
+  if (bit_depth == 16)
+    png_set_swap(png_ptr); // convert to little-endian format
+  
+  png_read_update_info(png_ptr, info_ptr); // update with transformations
+  int rowbytes = png_get_rowbytes (png_ptr, info_ptr);
+  bit_depth = png_get_bit_depth (png_ptr, info_ptr);
+  color_type = png_get_color_type (png_ptr, info_ptr);
+  
+  png_bytep* row_pointers = new png_bytep [height];
+  int i;
+  for (i = 0; i < height; i++)
+    row_pointers[i] = new unsigned char [rowbytes];
+  
+  png_read_image(png_ptr, row_pointers);
+  
+  setArraySize (width, height);
+  ImageFileArray v = getArray();
+  for (int iy = 0; iy < height; iy++) {
+    for (int ix = 0; ix < width; ix++) {
+      double dV;
+      if (color_type == PNG_COLOR_TYPE_GRAY) {
+        if (bit_depth == 8)
+          dV = row_pointers[iy][ix] / 255.;
+        else if (bit_depth == 16) {
+          int iBase = ix * 2;
+          dV = (row_pointers[iy][iBase] + (row_pointers[iy][iBase+1] << 8)) / 65536.;
+        }
+      } else if (color_type == PNG_COLOR_TYPE_RGB) {
+        if (bit_depth == 8) {
+          int iBase = ix * 3;
+          double dR = row_pointers[iy][iBase] / 255.;
+          double dG = row_pointers[iy][iBase+1] / 255.;
+          double dB = row_pointers[iy][iBase+2] / 255.;
+          dV = colorToGrayscale (dR, dG, dR);
+        }
+      }
+      v[ix][height-iy-1] = dV;
+    }
+  }
+  
+  png_read_end(png_ptr, end_info);
+  png_destroy_read_struct(&png_ptr, &info_ptr, &end_info);
+  
+  for (i = 0; i < height; i++)
+    delete row_pointers[i];
+  delete row_pointers;
+  
+  fclose (fp);
+  return true;
+}
+#endif
+
+bool
+ImageFile::exportImage (const char* const pszFormat, const char* const pszFilename, int nxcell, int nycell, double densmin, double densmax)
+{
+  int iFormatID = convertExportFormatNameToID (pszFormat);
+  
+  if (iFormatID == EXPORT_FORMAT_PGM)
     return writeImagePGM (pszFilename, nxcell, nycell, densmin, densmax);
-  else if (iFormatID == FORMAT_PGMASCII)
+  else if (iFormatID == EXPORT_FORMAT_PGMASCII)
     return writeImagePGMASCII (pszFilename, nxcell, nycell, densmin, densmax);
-  else if (iFormatID == FORMAT_PNG)
+  else if (iFormatID == EXPORT_FORMAT_TEXT)
+    return writeImageText (pszFilename);
+#ifdef HAVE_PNG
+  else if (iFormatID == EXPORT_FORMAT_PNG)
     return writeImagePNG (pszFilename, 8, nxcell, nycell, densmin, densmax);
-  else if (iFormatID == FORMAT_PNG16)
+  else if (iFormatID == EXPORT_FORMAT_PNG16)
     return writeImagePNG (pszFilename, 16, nxcell, nycell, densmin, densmax);
+#endif
+#ifdef HAVE_CTN_DICOM
+  else if (iFormatID == EXPORT_FORMAT_DICOM) {
+    DicomExporter dicomExport (this);
+    bool bSuccess = dicomExport.writeFile (pszFilename);
+    if (! bSuccess) 
+      sys_error (ERR_SEVERE, dicomExport.failMessage().c_str());
+    return bSuccess;
+  }
+#endif
   
   sys_error (ERR_SEVERE, "Invalid format %s [ImageFile::exportImage]", pszFormat);
   return false;
 }
 
+
 bool
 ImageFile::writeImagePGM (const char* const outfile, int nxcell, int nycell, double densmin, double densmax)
 {
@@ -1198,6 +1750,36 @@ ImageFile::writeImagePGMASCII (const char* const outfile, int nxcell, int nycell
   return true;
 }
 
+bool
+ImageFile::writeImageText (const char* const outfile)
+{
+  FILE *fp;
+  int nx = m_nx;
+  int ny = m_ny;
+  ImageFileArray v = getArray();
+  ImageFileArray vImag = getImaginaryArray();
+  
+  if ((fp = fopen (outfile, "w")) == NULL)
+    return false;
+  
+  for (int irow = ny - 1; irow >= 0; irow--) {
+    for (int icol = 0; icol < nx; icol++) {
+      if (isComplex()) {
+        if (vImag[icol][irow] >= 0)
+          fprintf (fp, "%.9g+%.9gi ", v[icol][irow], vImag[icol][irow]);
+        else
+          fprintf (fp, "%.9g-%.9gi ", v[icol][irow], -vImag[icol][irow]);
+      } else
+        fprintf (fp, "%12.8g ", v[icol][irow]);
+    }
+    fprintf(fp, "\n");
+  }
+  
+  fclose(fp);
+  
+  return true;
+}
+
 
 #ifdef HAVE_PNG
 bool
@@ -1305,7 +1887,7 @@ ImageFile::writeImageGIF (const char* const outfile, int nxcell, int nycell, dou
   
   FILE *out;
   if ((out = fopen (outfile,"w")) == NULL) {
-    sys_error(ERR_FATAL, "Error opening output file %s for writing", outfile);
+    sys_error(ERR_SEVERE, "Error opening output file %s for writing", outfile);
     return false;
   }
   gdImageGif(gif,out);