X-Git-Url: http://git.kpe.io/?p=ctsim.git;a=blobdiff_plain;f=libctsim%2Fprojections.cpp;h=e9e053d96cd0d19b357fa2efe2b5125ff465a07c;hp=bdf607369f96d645a7226af85d41a90128828fe2;hb=9776c9a12ba53419d34563a5ec57c90e3d6798f4;hpb=23f5654dacb1952c15bda92c2606fae3a55e48ad diff --git a/libctsim/projections.cpp b/libctsim/projections.cpp index bdf6073..e9e053d 100644 --- a/libctsim/projections.cpp +++ b/libctsim/projections.cpp @@ -6,9 +6,9 @@ ** Date Started: Aug 84 ** ** This is part of the CTSim program -** Copyright (C) 1983-2000 Kevin Rosenberg +** Copyright (c) 1983-2001 Kevin Rosenberg ** -** $Id: projections.cpp,v 1.37 2001/01/04 21:28:41 kevin Exp $ +** $Id: projections.cpp,v 1.49 2001/02/22 18:22:40 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 @@ -33,18 +33,18 @@ const int Projections::POLAR_INTERP_NEAREST = 0; const int Projections::POLAR_INTERP_BILINEAR = 1; const int Projections::POLAR_INTERP_BICUBIC = 2; -const char* Projections::s_aszInterpName[] = +const char* const Projections::s_aszInterpName[] = { {"nearest"}, {"bilinear"}, - {"bicubic"}, +// {"bicubic"}, }; -const char* Projections::s_aszInterpTitle[] = +const char* const Projections::s_aszInterpTitle[] = { {"Nearest"}, {"Bilinear"}, - {"Bicubic"}, +// {"Bicubic"}, }; const int Projections::s_iInterpCount = sizeof(s_aszInterpName) / sizeof(char*); @@ -148,14 +148,14 @@ Projections::initFromScanner (const Scanner& scanner) deleteProjData(); init (scanner.nView(), scanner.nDet()); - m_phmLen = scanner.phmLen(); m_rotInc = scanner.rotInc(); m_detInc = scanner.detInc(); m_geometry = scanner.geometry(); - m_focalLength = scanner.focalLength(); - m_fieldOfView = scanner.fieldOfView(); + m_dFocalLength = scanner.focalLength(); + m_dViewDiameter = scanner.viewDiameter(); m_rotStart = 0; m_detStart = -(scanner.detLen() / 2); + m_dFanBeamAngle = scanner.fanBeamAngle(); } void @@ -230,10 +230,10 @@ Projections::headerWrite (fnetorderstream& fs) kfloat64 _rotInc = m_rotInc; kfloat64 _detStart = m_detStart; kfloat64 _detInc = m_detInc; - kfloat64 _phmLen = m_phmLen; - kfloat64 _fieldOfView = m_fieldOfView; - kfloat64 _focalLength = m_focalLength; - + kfloat64 _viewDiameter = m_dViewDiameter; + kfloat64 _focalLength = m_dFocalLength; + kfloat64 _fanBeamAngle = m_dFanBeamAngle; + fs.seekp(0); if (! fs) return false; @@ -248,9 +248,9 @@ Projections::headerWrite (fnetorderstream& fs) fs.writeFloat64 (_rotInc); fs.writeFloat64 (_detStart); fs.writeFloat64 (_detInc); - fs.writeFloat64 (_phmLen); + fs.writeFloat64 (_viewDiameter); fs.writeFloat64 (_focalLength); - fs.writeFloat64 (_fieldOfView); + fs.writeFloat64 (_fanBeamAngle); fs.writeInt16 (_year); fs.writeInt16 (_month); fs.writeInt16 (_day); @@ -279,7 +279,7 @@ Projections::headerRead (fnetorderstream& fs) { kuint16 _hsize, _signature, _year, _month, _day, _hour, _minute, _second, _remarksize = 0; kuint32 _nView, _nDet, _geom; - kfloat64 _calcTime, _rotStart, _rotInc, _detStart, _detInc, _phmLen, _focalLength, _fieldOfView; + kfloat64 _calcTime, _rotStart, _rotInc, _detStart, _detInc, _focalLength, _viewDiameter, _fanBeamAngle; fs.seekg(0); if (! fs) @@ -295,9 +295,9 @@ Projections::headerRead (fnetorderstream& fs) fs.readFloat64 (_rotInc); fs.readFloat64 (_detStart); fs.readFloat64 (_detInc); - fs.readFloat64 (_phmLen); + fs.readFloat64 (_viewDiameter); fs.readFloat64 (_focalLength); - fs.readFloat64 (_fieldOfView); + fs.readFloat64 (_fanBeamAngle); fs.readInt16 (_year); fs.readInt16 (_month); fs.readInt16 (_day); @@ -341,9 +341,9 @@ Projections::headerRead (fnetorderstream& fs) m_rotInc = _rotInc; m_detStart = _detStart; m_detInc = _detInc; - m_phmLen = _phmLen; - m_focalLength = _focalLength; - m_fieldOfView = _fieldOfView; + m_dFocalLength = _focalLength; + m_dViewDiameter = _viewDiameter; + m_dFanBeamAngle = _fanBeamAngle; m_year = _year; m_month = _month; m_day = _day; @@ -404,10 +404,16 @@ Projections::copyViewData (const std::string& filename, std::ostream& os, int st bool Projections::copyViewData (const char* const filename, std::ostream& os, int startView, int endView) { - frnetorderstream is (filename, ios::in | ios::binary); + frnetorderstream is (filename, std::ios::in | std::ios::binary); kuint16 sizeHeader, signature; kuint32 _nView, _nDet; + is.seekg (0); + if (is.fail()) { + sys_error (ERR_SEVERE, "Unable to read projection file %s", filename); + return false; + } + is.readInt16 (sizeHeader); is.readInt16 (signature); is.readInt32 (_nView); @@ -416,7 +422,7 @@ Projections::copyViewData (const char* const filename, std::ostream& os, int sta int nDet = _nDet; if (signature != m_signature) { - sys_error (ERR_FATAL, "Illegal signature in projection file %s", filename); + sys_error (ERR_SEVERE, "Illegal signature in projection file %s", filename); return false; } @@ -446,9 +452,9 @@ Projections::copyViewData (const char* const filename, std::ostream& os, int sta delete pViewData; if (is.fail()) - sys_error (ERR_FATAL, "Error reading projection file"); + sys_error (ERR_SEVERE, "Error reading projection file"); if (os.fail()) - sys_error (ERR_FATAL, "Error writing projection file"); + sys_error (ERR_SEVERE, "Error writing projection file"); return (! (is.fail() | os.fail())); } @@ -468,20 +474,20 @@ Projections::copyHeader (const char* const filename, std::ostream& os) is.readInt16 (signature); is.seekg (0); if (signature != m_signature) { - sys_error (ERR_FATAL, "Illegal signature in projection file %s", filename); + sys_error (ERR_SEVERE, "Illegal signature in projection file %s", filename); return false; } unsigned char* pHdrData = new unsigned char [sizeHeader]; is.read (reinterpret_cast(pHdrData), sizeHeader); if (is.fail()) { - sys_error (ERR_FATAL, "Error reading header"); + sys_error (ERR_SEVERE, "Error reading header"); return false; } os.write (reinterpret_cast(pHdrData), sizeHeader); if (os.fail()) { - sys_error (ERR_FATAL, "Error writing header"); + sys_error (ERR_SEVERE, "Error writing header"); return false; } @@ -625,12 +631,15 @@ Projections::printProjectionData (int startView, int endView) printf("Description: %s\n", m_remark.c_str()); printf("Geometry: %s\n", Scanner::convertGeometryIDToName (m_geometry)); printf("nView = %8d nDet = %8d\n", m_nView, m_nDet); - printf("focalLength = %8.4f fieldOfView = %8.4f\n", m_focalLength, m_fieldOfView); + printf("focalLength = %8.4f ViewDiameter = %8.4f\n", m_dFocalLength, m_dViewDiameter); + printf("fanBeamAngle= %8.4f\n", convertRadiansToDegrees(m_dFanBeamAngle)); printf("rotStart = %8.4f rotInc = %8.4f\n", m_rotStart, m_rotInc); printf("detStart = %8.4f detInc = %8.4f\n", m_detStart, m_detInc); if (m_projData != NULL) { if (startView < 0) startView = 0; + if (endView < 0) + endView = m_nView - 1; if (startView > m_nView - 1) startView = m_nView - 1; if (endView > m_nView - 1) @@ -649,16 +658,16 @@ void Projections::printScanInfo (std::ostringstream& os) const { os << "Number of detectors: " << m_nDet << "\n"; - os << " Number of views: " << m_nView<< "\n"; - os << " Remark: " << m_remark.c_str()<< "\n"; - os << " Geometry: " << Scanner::convertGeometryIDToName (m_geometry)<< "\n"; - os << " Focal Length: " << m_focalLength<< "\n"; - os << " Field Of View: " << m_fieldOfView<< "\n"; - os << " phmLen: " << m_phmLen<< "\n"; - os << " detStart: " << m_detStart<< "\n"; - os << " detInc: " << m_detInc<< "\n"; - os << " rotStart: " << m_rotStart<< "\n"; - os << " rotInc: " << m_rotInc<< "\n"; + os << "Number of views: " << m_nView<< "\n"; + os << "Description: " << m_remark.c_str()<< "\n"; + os << "Geometry: " << Scanner::convertGeometryIDToName (m_geometry)<< "\n"; + os << "Focal Length: " << m_dFocalLength<< "\n"; + os << "View Diameter: " << m_dViewDiameter<< "\n"; + os << "Fan Beam Angle: " << convertRadiansToDegrees(m_dFanBeamAngle) << "\n"; + os << "detStart: " << m_detStart<< "\n"; + os << "detInc: " << m_detInc<< "\n"; + os << "rotStart: " << m_rotStart<< "\n"; + os << "rotInc: " << m_rotInc<< "\n"; } @@ -678,10 +687,12 @@ Projections::convertPolar (ImageFile& rIF, int iInterpolationID) double** ppdView = adView.getArray(); double** ppdDet = adDet.getArray(); - calcArrayPolarCoordinates (nx, ny, ppdView, ppdDet); + if (! calcArrayPolarCoordinates (nx, ny, ppdView, ppdDet)) + return false; std::complex** ppcDetValue = new std::complex* [m_nView]; - for (unsigned int iView = 0; iView < m_nView; iView++) { + unsigned int iView; + for (iView = 0; iView < m_nView; iView++) { ppcDetValue[iView] = new std::complex [m_nDet]; for (unsigned int iDet = 0; iDet < m_nDet; iDet++) ppcDetValue[iView][iDet] = std::complex(getDetectorArray (iView).detValues()[iDet], 0); @@ -696,13 +707,71 @@ Projections::convertPolar (ImageFile& rIF, int iInterpolationID) return true; } -void + +bool +Projections::convertFFTPolar (ImageFile& rIF, int iInterpolationID, int iZeropad) +{ + unsigned int nx = rIF.nx(); + unsigned int ny = rIF.ny(); + ImageFileArray v = rIF.getArray(); + if (! rIF.isComplex()) + rIF.convertRealToComplex(); + ImageFileArray vImag = rIF.getImaginaryArray(); + + if (! v || nx == 0 || ny == 0) + return false; + +#ifndef HAVE_FFT + return false; +#else + Array2d adView (nx, ny); + Array2d adDet (nx, ny); + double** ppdView = adView.getArray(); + double** ppdDet = adDet.getArray(); + + std::complex** ppcDetValue = new std::complex* [m_nView]; + unsigned int iView; + double* pdDet = new double [m_nDet]; + fftw_complex* pcIn = new fftw_complex [m_nDet]; + fftw_plan plan = fftw_create_plan (m_nDet, FFTW_FORWARD, FFTW_IN_PLACE); + + for (iView = 0; iView < m_nView; iView++) { + unsigned int iDet; + for (iDet = 0; iDet < m_nDet; iDet++) { + pcIn[iDet].re = getDetectorArray(iView).detValues()[iDet]; + pcIn[iDet].im = 0; + } + fftw_one (plan, pcIn, NULL); + ppcDetValue[iView] = new std::complex [m_nDet]; + for (iDet = 0; iDet < m_nDet; iDet++) + ppcDetValue[iView][iDet] = std::complex (pcIn[iDet].re, pcIn[iDet].im); + Fourier::shuffleFourierToNaturalOrder (ppcDetValue[iView], m_nDet); + } + + fftw_destroy_plan (plan); + delete [] pcIn; + + bool bError = calcArrayPolarCoordinates (nx, ny, ppdView, ppdDet); + + if (! bError) + interpolatePolar (v, vImag, nx, ny, ppcDetValue, ppdView, ppdDet, m_nView, m_nDet, iInterpolationID); + + for (iView = 0; iView < m_nView; iView++) + delete [] ppcDetValue[iView]; + delete [] ppcDetValue; + + return bError; +#endif +} + + +bool Projections::calcArrayPolarCoordinates (unsigned int nx, unsigned int ny, double** ppdView, double** ppdDet) { - double xMin = -m_phmLen / 2; - double xMax = xMin + m_phmLen; - double yMin = -m_phmLen / 2; - double yMax = yMin + m_phmLen; + double xMin = -phmLen() / 2; + double xMax = xMin + phmLen(); + double yMin = -phmLen() / 2; + double yMax = yMin + phmLen(); double xInc = (xMax - xMin) / nx; // size of cells double yInc = (yMax - yMin) / ny; @@ -711,8 +780,10 @@ Projections::calcArrayPolarCoordinates (unsigned int nx, unsigned int ny, double if (m_geometry != Scanner::GEOMETRY_PARALLEL) { sys_error (ERR_WARNING, "convertPolar supports Parallel only"); + return false; } + // Calculates polar coordinates (view#, det#) for each point on phantom grid double x = xMin + xInc / 2; // Rectang coords of center of pixel for (unsigned int ix = 0; ix < nx; x += xInc, ix++) { double y = yMin + yInc / 2; @@ -731,6 +802,8 @@ Projections::calcArrayPolarCoordinates (unsigned int nx, unsigned int ny, double ppdDet[ix][iy] = (r / m_detInc) + iDetCenter; } } + + return true; } void @@ -743,8 +816,10 @@ Projections::interpolatePolar (ImageFileArray& v, ImageFileArray& vImag, if (iInterpolationID == POLAR_INTERP_NEAREST) { int iView = nearest (ppdView[ix][iy]); int iDet = nearest (ppdDet[ix][iy]); - if (iView == nView) - iView = nView - 1; + if (iView == nView) { + iView = 0; + // iDet = m_nDet - iDet; + } if (iDet >= 0 && iDet < nDet && iView >= 0 && iView < nView) { v[ix][iy] = ppcDetValue[iView][iDet].real(); if (vImag) @@ -755,9 +830,9 @@ Projections::interpolatePolar (ImageFileArray& v, ImageFileArray& vImag, v[ix][iy] = 0; } } else if (iInterpolationID == POLAR_INTERP_BILINEAR) { - int iFloorView = ppdView[ix][iy]; + int iFloorView = static_cast(ppdView[ix][iy]); double dFracView = ppdView[ix][iy] - iFloorView; - int iFloorDet = ppdDet[ix][iy]; + int iFloorDet = static_cast(ppdDet[ix][iy]); double dFracDet = ppdDet[ix][iy] - iFloorDet; if (iFloorDet >= 0 && iFloorView >= 0) { @@ -766,7 +841,7 @@ Projections::interpolatePolar (ImageFileArray& v, ImageFileArray& vImag, if (iFloorView < nView - 1) v2 = ppcDetValue[iFloorView + 1][iFloorDet]; else - v2 = v1; + v2 = ppcDetValue[0][iFloorDet]; if (iFloorDet < nDet - 1) v4 = ppcDetValue[iFloorView][iFloorDet+1]; else @@ -774,10 +849,12 @@ Projections::interpolatePolar (ImageFileArray& v, ImageFileArray& vImag, if (iFloorView < nView - 1 && iFloorDet < nDet - 1) v3 = ppcDetValue [iFloorView+1][iFloorDet+1]; else if (iFloorView < nView - 1) - v4 = v2; + v3 = v2; else - v4 = v1; - std::complex vInterp = v1 + dFracView * (v2 - v1) + dFracDet * (v4 - v1); + v3 = ppcDetValue[0][iFloorDet+1]; + std::complex vInterp = (1 - dFracView) * (1 - dFracDet) * v1 + + dFracView * (1 - dFracDet) * v2 + dFracView * dFracDet * v3 + + dFracDet * (1 - dFracView) * v4; v[ix][iy] = vInterp.real(); if (vImag) vImag[ix][iy] = vInterp.imag(); @@ -797,57 +874,4 @@ Projections::interpolatePolar (ImageFileArray& v, ImageFileArray& vImag, } } -bool -Projections::convertFFTPolar (ImageFile& rIF, int iInterpolationID, int iZeropad) -{ - unsigned int nx = rIF.nx(); - unsigned int ny = rIF.ny(); - ImageFileArray v = rIF.getArray(); - if (! rIF.isComplex()) - rIF.convertRealToComplex(); - ImageFileArray vImag = rIF.getImaginaryArray(); - - if (! v || nx == 0 || ny == 0) - return false; - - Array2d adView (nx, ny); - Array2d adDet (nx, ny); - double** ppdView = adView.getArray(); - double** ppdDet = adDet.getArray(); - - std::complex** ppcDetValue = new std::complex* [m_nView]; - unsigned int iView; - double* pdDet = new double [m_nDet]; - fftw_complex* pcIn = new fftw_complex [m_nDet]; - fftw_plan plan = fftw_create_plan (m_nDet, FFTW_FORWARD, FFTW_IN_PLACE); - - for (iView = 0; iView < m_nView; iView++) { - ppcDetValue[iView] = new std::complex [m_nDet]; - unsigned int iDet; - for (iDet = 0; iDet < m_nDet; iDet++) { - pcIn[iDet].re = getDetectorArray(iView).detValues()[iDet]; - pcIn[iDet].im = 0; - } - fftw_one (plan, pcIn, NULL); - for (iDet = 0; iDet < m_nDet; iDet++) - ppcDetValue[iView][iDet] = std::complex (pcIn[iDet].re, pcIn[iDet].im); - Fourier::shuffleFourierToNaturalOrder (ppcDetValue[iView], m_nDet); - } - - fftw_destroy_plan (plan); - - delete [] pcIn; - - calcArrayPolarCoordinates (nx, ny, ppdView, ppdDet); - - interpolatePolar (v, vImag, nx, ny, ppcDetValue, ppdView, ppdDet, m_nView, m_nDet, iInterpolationID); - - for (iView = 0; iView < m_nView; iView++) - delete [] ppcDetValue[iView]; - delete [] ppcDetValue; - - return true; -} - -