diff options
author | Willem Jan Palenstijn <Willem.Jan.Palenstijn@cwi.nl> | 2018-07-05 16:13:37 +0200 |
---|---|---|
committer | Willem Jan Palenstijn <Willem.Jan.Palenstijn@cwi.nl> | 2018-07-17 11:29:22 +0200 |
commit | d58f6b51493c7931cbc02feb3ffeb0eca6ea3a4e (patch) | |
tree | 51a4b27c7825198ec119c4e420439f183584420a /src | |
parent | 0d06afc38d7a8443a079d25d72ec4b4b15353974 (diff) | |
download | astra-d58f6b51493c7931cbc02feb3ffeb0eca6ea3a4e.tar.gz astra-d58f6b51493c7931cbc02feb3ffeb0eca6ea3a4e.tar.bz2 astra-d58f6b51493c7931cbc02feb3ffeb0eca6ea3a4e.tar.xz astra-d58f6b51493c7931cbc02feb3ffeb0eca6ea3a4e.zip |
Refactor a few filter-related functions out of cuda code
Diffstat (limited to 'src')
-rw-r--r-- | src/CudaFilteredBackProjectionAlgorithm.cpp | 118 | ||||
-rw-r--r-- | src/Filters.cpp | 484 |
2 files changed, 488 insertions, 114 deletions
diff --git a/src/CudaFilteredBackProjectionAlgorithm.cpp b/src/CudaFilteredBackProjectionAlgorithm.cpp index 944f165..9ebbd60 100644 --- a/src/CudaFilteredBackProjectionAlgorithm.cpp +++ b/src/CudaFilteredBackProjectionAlgorithm.cpp @@ -27,10 +27,10 @@ along with the ASTRA Toolbox. If not, see <http://www.gnu.org/licenses/>. #include <astra/CudaFilteredBackProjectionAlgorithm.h> #include <astra/FanFlatProjectionGeometry2D.h> -#include <cstring> #include "astra/AstraObjectManager.h" #include "astra/CudaProjector2D.h" +#include "astra/Filters.h" #include "astra/cuda/2d/astra.h" #include "astra/cuda/2d/fbp.h" @@ -75,7 +75,7 @@ bool CCudaFilteredBackProjectionAlgorithm::initialize(const Config& _cfg) // filter type XMLNode node = _cfg.self.getSingleNode("FilterType"); if (node) - m_eFilter = _convertStringToFilter(node.getContent().c_str()); + m_eFilter = convertStringToFilter(node.getContent().c_str()); else m_eFilter = FILTER_RAMLAK; CC.markNodeParsed("FilterType"); @@ -215,6 +215,8 @@ bool CCudaFilteredBackProjectionAlgorithm::check() ASTRA_CONFIG_CHECK(m_pSinogram, "FBP_CUDA", "Invalid Projection Data Object."); ASTRA_CONFIG_CHECK(m_pReconstruction, "FBP_CUDA", "Invalid Reconstruction Data Object."); + ASTRA_CONFIG_CHECK(m_eFilter != FILTER_ERROR, "FBP_CUDA", "Invalid filter name."); + if((m_eFilter == FILTER_PROJECTION) || (m_eFilter == FILTER_SINOGRAM) || (m_eFilter == FILTER_RPROJECTION) || (m_eFilter == FILTER_RSINOGRAM)) { ASTRA_CONFIG_CHECK(m_pfFilter, "FBP_CUDA", "Invalid filter pointer."); @@ -235,116 +237,4 @@ bool CCudaFilteredBackProjectionAlgorithm::check() return true; } -static bool stringCompareLowerCase(const char * _stringA, const char * _stringB) -{ - int iCmpReturn = 0; - -#ifdef _MSC_VER - iCmpReturn = _stricmp(_stringA, _stringB); -#else - iCmpReturn = strcasecmp(_stringA, _stringB); -#endif - - return (iCmpReturn == 0); -} - -E_FBPFILTER CCudaFilteredBackProjectionAlgorithm::_convertStringToFilter(const char * _filterType) -{ - E_FBPFILTER output = FILTER_NONE; - - if(stringCompareLowerCase(_filterType, "ram-lak")) - { - output = FILTER_RAMLAK; - } - else if(stringCompareLowerCase(_filterType, "shepp-logan")) - { - output = FILTER_SHEPPLOGAN; - } - else if(stringCompareLowerCase(_filterType, "cosine")) - { - output = FILTER_COSINE; - } - else if(stringCompareLowerCase(_filterType, "hamming")) - { - output = FILTER_HAMMING; - } - else if(stringCompareLowerCase(_filterType, "hann")) - { - output = FILTER_HANN; - } - else if(stringCompareLowerCase(_filterType, "none")) - { - output = FILTER_NONE; - } - else if(stringCompareLowerCase(_filterType, "tukey")) - { - output = FILTER_TUKEY; - } - else if(stringCompareLowerCase(_filterType, "lanczos")) - { - output = FILTER_LANCZOS; - } - else if(stringCompareLowerCase(_filterType, "triangular")) - { - output = FILTER_TRIANGULAR; - } - else if(stringCompareLowerCase(_filterType, "gaussian")) - { - output = FILTER_GAUSSIAN; - } - else if(stringCompareLowerCase(_filterType, "barlett-hann")) - { - output = FILTER_BARTLETTHANN; - } - else if(stringCompareLowerCase(_filterType, "blackman")) - { - output = FILTER_BLACKMAN; - } - else if(stringCompareLowerCase(_filterType, "nuttall")) - { - output = FILTER_NUTTALL; - } - else if(stringCompareLowerCase(_filterType, "blackman-harris")) - { - output = FILTER_BLACKMANHARRIS; - } - else if(stringCompareLowerCase(_filterType, "blackman-nuttall")) - { - output = FILTER_BLACKMANNUTTALL; - } - else if(stringCompareLowerCase(_filterType, "flat-top")) - { - output = FILTER_FLATTOP; - } - else if(stringCompareLowerCase(_filterType, "kaiser")) - { - output = FILTER_KAISER; - } - else if(stringCompareLowerCase(_filterType, "parzen")) - { - output = FILTER_PARZEN; - } - else if(stringCompareLowerCase(_filterType, "projection")) - { - output = FILTER_PROJECTION; - } - else if(stringCompareLowerCase(_filterType, "sinogram")) - { - output = FILTER_SINOGRAM; - } - else if(stringCompareLowerCase(_filterType, "rprojection")) - { - output = FILTER_RPROJECTION; - } - else if(stringCompareLowerCase(_filterType, "rsinogram")) - { - output = FILTER_RSINOGRAM; - } - else - { - ASTRA_ERROR("Failed to convert \"%s\" into a filter.",_filterType); - } - - return output; -} diff --git a/src/Filters.cpp b/src/Filters.cpp new file mode 100644 index 0000000..30b2f24 --- /dev/null +++ b/src/Filters.cpp @@ -0,0 +1,484 @@ +/* +----------------------------------------------------------------------- +Copyright: 2010-2018, imec Vision Lab, University of Antwerp + 2014-2018, CWI, Amsterdam + +Contact: astra@astra-toolbox.com +Website: http://www.astra-toolbox.com/ + +This file is part of the ASTRA Toolbox. + + +The ASTRA Toolbox is free software: you can redistribute it and/or modify +it under the terms of the GNU General Public License as published by +the Free Software Foundation, either version 3 of the License, or +(at your option) any later version. + +The ASTRA Toolbox is distributed in the hope that it will be useful, +but WITHOUT ANY WARRANTY; without even the implied warranty of +MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the +GNU General Public License for more details. + +You should have received a copy of the GNU General Public License +along with the ASTRA Toolbox. If not, see <http://www.gnu.org/licenses/>. + +----------------------------------------------------------------------- +*/ + +#include "astra/Globals.h" +#include "astra/Logging.h" +#include "astra/Fourier.h" +#include "astra/Filters.h" + +#include <utility> +#include <cstring> + +namespace astra { + +float *genFilter(E_FBPFILTER _eFilter, float _fD, int _iProjectionCount, + int _iFFTRealDetectorCount, + int _iFFTFourierDetectorCount, float _fParameter /* = -1.0f */) +{ + float * pfFilt = new float[_iFFTFourierDetectorCount]; + float * pfW = new float[_iFFTFourierDetectorCount]; + + // We cache one Fourier transform for repeated FBP's of the same size + static float *pfData = 0; + static int iFilterCacheSize = 0; + + if (!pfData || iFilterCacheSize != _iFFTRealDetectorCount) { + // Compute filter in spatial domain + + delete[] pfData; + pfData = new float[2*_iFFTRealDetectorCount]; + int *ip = new int[int(2+sqrt(_iFFTRealDetectorCount)+1)]; + ip[0] = 0; + float32 *w = new float32[_iFFTRealDetectorCount/2]; + + for (int i = 0; i < _iFFTRealDetectorCount; ++i) { + pfData[2*i+1] = 0.0f; + + if (i & 1) { + int j = i; + if (2*j > _iFFTRealDetectorCount) + j = _iFFTRealDetectorCount - j; + float f = M_PI * j; + pfData[2*i] = -1 / (f*f); + } else { + pfData[2*i] = 0.0f; + } + } + + pfData[0] = 0.25f; + + cdft(2*_iFFTRealDetectorCount, -1, pfData, ip, w); + delete[] ip; + delete[] w; + + iFilterCacheSize = _iFFTRealDetectorCount; + } + + for(int iDetectorIndex = 0; iDetectorIndex < _iFFTFourierDetectorCount; iDetectorIndex++) + { + float fRelIndex = (float)iDetectorIndex / (float)_iFFTRealDetectorCount; + + pfFilt[iDetectorIndex] = 2.0f * pfData[2*iDetectorIndex]; + pfW[iDetectorIndex] = M_PI * 2.0f * fRelIndex; + } + + switch(_eFilter) + { + case FILTER_RAMLAK: + { + // do nothing + break; + } + case FILTER_SHEPPLOGAN: + { + // filt(2:end) = filt(2:end) .* (sin(w(2:end)/(2*d))./(w(2:end)/(2*d))) + for(int iDetectorIndex = 1; iDetectorIndex < _iFFTFourierDetectorCount; iDetectorIndex++) + { + pfFilt[iDetectorIndex] = pfFilt[iDetectorIndex] * (sinf(pfW[iDetectorIndex] / 2.0f / _fD) / (pfW[iDetectorIndex] / 2.0f / _fD)); + } + break; + } + case FILTER_COSINE: + { + // filt(2:end) = filt(2:end) .* cos(w(2:end)/(2*d)) + for(int iDetectorIndex = 1; iDetectorIndex < _iFFTFourierDetectorCount; iDetectorIndex++) + { + pfFilt[iDetectorIndex] = pfFilt[iDetectorIndex] * cosf(pfW[iDetectorIndex] / 2.0f / _fD); + } + break; + } + case FILTER_HAMMING: + { + // filt(2:end) = filt(2:end) .* (.54 + .46 * cos(w(2:end)/d)) + for(int iDetectorIndex = 1; iDetectorIndex < _iFFTFourierDetectorCount; iDetectorIndex++) + { + pfFilt[iDetectorIndex] = pfFilt[iDetectorIndex] * ( 0.54f + 0.46f * cosf(pfW[iDetectorIndex] / _fD)); + } + break; + } + case FILTER_HANN: + { + // filt(2:end) = filt(2:end) .*(1+cos(w(2:end)./d)) / 2 + for(int iDetectorIndex = 1; iDetectorIndex < _iFFTFourierDetectorCount; iDetectorIndex++) + { + pfFilt[iDetectorIndex] = pfFilt[iDetectorIndex] * (1.0f + cosf(pfW[iDetectorIndex] / _fD)) / 2.0f; + } + break; + } + case FILTER_TUKEY: + { + float fAlpha = _fParameter; + if(_fParameter < 0.0f) fAlpha = 0.5f; + float fN = (float)_iFFTFourierDetectorCount; + float fHalfN = fN / 2.0f; + float fEnumTerm = fAlpha * fHalfN; + float fDenom = (1.0f - fAlpha) * fHalfN; + float fBlockStart = fHalfN - fEnumTerm; + float fBlockEnd = fHalfN + fEnumTerm; + + for(int iDetectorIndex = 1; iDetectorIndex < _iFFTFourierDetectorCount; iDetectorIndex++) + { + float fAbsSmallN = fabs((float)iDetectorIndex); + float fStoredValue = 0.0f; + + if((fBlockStart <= fAbsSmallN) && (fAbsSmallN <= fBlockEnd)) + { + fStoredValue = 1.0f; + } + else + { + float fEnum = fAbsSmallN - fEnumTerm; + float fCosInput = M_PI * fEnum / fDenom; + fStoredValue = 0.5f * (1.0f + cosf(fCosInput)); + } + + pfFilt[iDetectorIndex] *= fStoredValue; + } + + break; + } + case FILTER_LANCZOS: + { + float fDenum = (float)(_iFFTFourierDetectorCount - 1); + + for(int iDetectorIndex = 1; iDetectorIndex < _iFFTFourierDetectorCount; iDetectorIndex++) + { + float fSmallN = (float)iDetectorIndex; + float fX = 2.0f * fSmallN / fDenum - 1.0f; + float fSinInput = M_PI * fX; + float fStoredValue = 0.0f; + + if(fabsf(fSinInput) > 0.001f) + { + fStoredValue = sin(fSinInput)/fSinInput; + } + else + { + fStoredValue = 1.0f; + } + + pfFilt[iDetectorIndex] *= fStoredValue; + } + + break; + } + case FILTER_TRIANGULAR: + { + float fNMinusOne = (float)(_iFFTFourierDetectorCount - 1); + + for(int iDetectorIndex = 1; iDetectorIndex < _iFFTFourierDetectorCount; iDetectorIndex++) + { + float fSmallN = (float)iDetectorIndex; + float fAbsInput = fSmallN - fNMinusOne / 2.0f; + float fParenInput = fNMinusOne / 2.0f - fabsf(fAbsInput); + float fStoredValue = 2.0f / fNMinusOne * fParenInput; + + pfFilt[iDetectorIndex] *= fStoredValue; + } + + break; + } + case FILTER_GAUSSIAN: + { + float fSigma = _fParameter; + if(_fParameter < 0.0f) fSigma = 0.4f; + float fN = (float)_iFFTFourierDetectorCount; + float fQuotient = (fN - 1.0f) / 2.0f; + + for(int iDetectorIndex = 1; iDetectorIndex < _iFFTFourierDetectorCount; iDetectorIndex++) + { + float fSmallN = (float)iDetectorIndex; + float fEnum = fSmallN - fQuotient; + float fDenom = fSigma * fQuotient; + float fPower = -0.5f * (fEnum / fDenom) * (fEnum / fDenom); + float fStoredValue = expf(fPower); + + pfFilt[iDetectorIndex] *= fStoredValue; + } + + break; + } + case FILTER_BARTLETTHANN: + { + const float fA0 = 0.62f; + const float fA1 = 0.48f; + const float fA2 = 0.38f; + float fNMinusOne = (float)(_iFFTFourierDetectorCount) - 1.0f; + + for(int iDetectorIndex = 1; iDetectorIndex < _iFFTFourierDetectorCount; iDetectorIndex++) + { + float fSmallN = (float)iDetectorIndex; + float fAbsInput = fSmallN / fNMinusOne - 0.5f; + float fFirstTerm = fA1 * fabsf(fAbsInput); + float fCosInput = 2.0f * M_PI * fSmallN / fNMinusOne; + float fSecondTerm = fA2 * cosf(fCosInput); + float fStoredValue = fA0 - fFirstTerm - fSecondTerm; + + pfFilt[iDetectorIndex] *= fStoredValue; + } + + break; + } + case FILTER_BLACKMAN: + { + float fAlpha = _fParameter; + if(_fParameter < 0.0f) fAlpha = 0.16f; + float fA0 = (1.0f - fAlpha) / 2.0f; + float fA1 = 0.5f; + float fA2 = fAlpha / 2.0f; + float fNMinusOne = (float)(_iFFTFourierDetectorCount - 1); + + for(int iDetectorIndex = 1; iDetectorIndex < _iFFTFourierDetectorCount; iDetectorIndex++) + { + float fSmallN = (float)iDetectorIndex; + float fCosInput1 = 2.0f * M_PI * 0.5f * fSmallN / fNMinusOne; + float fCosInput2 = 4.0f * M_PI * 0.5f * fSmallN / fNMinusOne; + float fStoredValue = fA0 - fA1 * cosf(fCosInput1) + fA2 * cosf(fCosInput2); + + pfFilt[iDetectorIndex] *= fStoredValue; + } + + break; + } + case FILTER_NUTTALL: + { + const float fA0 = 0.355768f; + const float fA1 = 0.487396f; + const float fA2 = 0.144232f; + const float fA3 = 0.012604f; + float fNMinusOne = (float)(_iFFTFourierDetectorCount) - 1.0f; + + for(int iDetectorIndex = 1; iDetectorIndex < _iFFTFourierDetectorCount; iDetectorIndex++) + { + float fSmallN = (float)iDetectorIndex; + float fBaseCosInput = M_PI * fSmallN / fNMinusOne; + float fFirstTerm = fA1 * cosf(2.0f * fBaseCosInput); + float fSecondTerm = fA2 * cosf(4.0f * fBaseCosInput); + float fThirdTerm = fA3 * cosf(6.0f * fBaseCosInput); + float fStoredValue = fA0 - fFirstTerm + fSecondTerm - fThirdTerm; + + pfFilt[iDetectorIndex] *= fStoredValue; + } + + break; + } + case FILTER_BLACKMANHARRIS: + { + const float fA0 = 0.35875f; + const float fA1 = 0.48829f; + const float fA2 = 0.14128f; + const float fA3 = 0.01168f; + float fNMinusOne = (float)(_iFFTFourierDetectorCount) - 1.0f; + + for(int iDetectorIndex = 1; iDetectorIndex < _iFFTFourierDetectorCount; iDetectorIndex++) + { + float fSmallN = (float)iDetectorIndex; + float fBaseCosInput = M_PI * fSmallN / fNMinusOne; + float fFirstTerm = fA1 * cosf(2.0f * fBaseCosInput); + float fSecondTerm = fA2 * cosf(4.0f * fBaseCosInput); + float fThirdTerm = fA3 * cosf(6.0f * fBaseCosInput); + float fStoredValue = fA0 - fFirstTerm + fSecondTerm - fThirdTerm; + + pfFilt[iDetectorIndex] *= fStoredValue; + } + + break; + } + case FILTER_BLACKMANNUTTALL: + { + const float fA0 = 0.3635819f; + const float fA1 = 0.4891775f; + const float fA2 = 0.1365995f; + const float fA3 = 0.0106411f; + float fNMinusOne = (float)(_iFFTFourierDetectorCount) - 1.0f; + + for(int iDetectorIndex = 1; iDetectorIndex < _iFFTFourierDetectorCount; iDetectorIndex++) + { + float fSmallN = (float)iDetectorIndex; + float fBaseCosInput = M_PI * fSmallN / fNMinusOne; + float fFirstTerm = fA1 * cosf(2.0f * fBaseCosInput); + float fSecondTerm = fA2 * cosf(4.0f * fBaseCosInput); + float fThirdTerm = fA3 * cosf(6.0f * fBaseCosInput); + float fStoredValue = fA0 - fFirstTerm + fSecondTerm - fThirdTerm; + + pfFilt[iDetectorIndex] *= fStoredValue; + } + + break; + } + case FILTER_FLATTOP: + { + const float fA0 = 1.0f; + const float fA1 = 1.93f; + const float fA2 = 1.29f; + const float fA3 = 0.388f; + const float fA4 = 0.032f; + float fNMinusOne = (float)(_iFFTFourierDetectorCount) - 1.0f; + + for(int iDetectorIndex = 1; iDetectorIndex < _iFFTFourierDetectorCount; iDetectorIndex++) + { + float fSmallN = (float)iDetectorIndex; + float fBaseCosInput = M_PI * fSmallN / fNMinusOne; + float fFirstTerm = fA1 * cosf(2.0f * fBaseCosInput); + float fSecondTerm = fA2 * cosf(4.0f * fBaseCosInput); + float fThirdTerm = fA3 * cosf(6.0f * fBaseCosInput); + float fFourthTerm = fA4 * cosf(8.0f * fBaseCosInput); + float fStoredValue = fA0 - fFirstTerm + fSecondTerm - fThirdTerm + fFourthTerm; + + pfFilt[iDetectorIndex] *= fStoredValue; + } + + break; + } + case FILTER_KAISER: + { + float fAlpha = _fParameter; + if(_fParameter < 0.0f) fAlpha = 3.0f; + float fPiTimesAlpha = M_PI * fAlpha; + float fNMinusOne = (float)(_iFFTFourierDetectorCount - 1); + float fDenom = (float)j0((double)fPiTimesAlpha); + + for(int iDetectorIndex = 1; iDetectorIndex < _iFFTFourierDetectorCount; iDetectorIndex++) + { + float fSmallN = (float)iDetectorIndex; + float fSquareInput = 2.0f * fSmallN / fNMinusOne - 1; + float fSqrtInput = 1.0f - fSquareInput * fSquareInput; + float fBesselInput = fPiTimesAlpha * sqrt(fSqrtInput); + float fEnum = (float)j0((double)fBesselInput); + float fStoredValue = fEnum / fDenom; + + pfFilt[iDetectorIndex] *= fStoredValue; + } + + break; + } + case FILTER_PARZEN: + { + for(int iDetectorIndex = 1; iDetectorIndex < _iFFTFourierDetectorCount; iDetectorIndex++) + { + float fSmallN = (float)iDetectorIndex; + float fQ = fSmallN / (float)(_iFFTFourierDetectorCount - 1); + float fStoredValue = 0.0f; + + if(fQ <= 0.5f) + { + fStoredValue = 1.0f - 6.0f * fQ * fQ * (1.0f - fQ); + } + else + { + float fCubedValue = 1.0f - fQ; + fStoredValue = 2.0f * fCubedValue * fCubedValue * fCubedValue; + } + + pfFilt[iDetectorIndex] *= fStoredValue; + } + + break; + } + default: + { + ASTRA_ERROR("Cannot serve requested filter"); + } + } + + // filt(w>pi*d) = 0; + float fPiTimesD = M_PI * _fD; + for(int iDetectorIndex = 0; iDetectorIndex < _iFFTFourierDetectorCount; iDetectorIndex++) + { + float fWValue = pfW[iDetectorIndex]; + + if(fWValue > fPiTimesD) + { + pfFilt[iDetectorIndex] = 0.0f; + } + } + + delete[] pfW; + + return pfFilt; +} + +static bool stringCompareLowerCase(const char * _stringA, const char * _stringB) +{ + int iCmpReturn = 0; + +#ifdef _MSC_VER + iCmpReturn = _stricmp(_stringA, _stringB); +#else + iCmpReturn = strcasecmp(_stringA, _stringB); +#endif + + return (iCmpReturn == 0); +} + +struct FilterNameMapEntry { + const char *m_name; + E_FBPFILTER m_type; +}; + +E_FBPFILTER convertStringToFilter(const char * _filterType) +{ + + static const FilterNameMapEntry map[] = { + { "ram-lak", FILTER_RAMLAK }, + { "shepp-logan", FILTER_SHEPPLOGAN }, + { "cosine", FILTER_COSINE }, + { "hamming", FILTER_HAMMING}, + { "hann", FILTER_HANN}, + { "tukey", FILTER_TUKEY }, + { "lanczos", FILTER_LANCZOS}, + { "triangular", FILTER_TRIANGULAR}, + { "gaussian", FILTER_GAUSSIAN}, + { "barlett-hann", FILTER_BARTLETTHANN }, + { "blackman", FILTER_BLACKMAN}, + { "nuttall", FILTER_NUTTALL }, + { "blackman-harris", FILTER_BLACKMANHARRIS }, + { "blackman-nuttall", FILTER_BLACKMANNUTTALL }, + { "flat-top", FILTER_FLATTOP }, + { "kaiser", FILTER_KAISER }, + { "parzen", FILTER_PARZEN }, + { "projection", FILTER_PROJECTION }, + { "sinogram", FILTER_SINOGRAM }, + { "rprojection", FILTER_RPROJECTION }, + { "rsinogram", FILTER_RSINOGRAM }, + { "none", FILTER_NONE }, + { 0, FILTER_ERROR } }; + + const FilterNameMapEntry *i; + + for (i = &map[0]; i->m_name; ++i) + if (stringCompareLowerCase(_filterType, i->m_name)) + return i->m_type; + + ASTRA_ERROR("Failed to convert \"%s\" into a filter.",_filterType); + + return FILTER_ERROR; +} + + + +} |