diff options
| author | Willem Jan Palenstijn <Willem.Jan.Palenstijn@cwi.nl> | 2014-10-01 16:52:21 +0200 | 
|---|---|---|
| committer | Willem Jan Palenstijn <Willem.Jan.Palenstijn@cwi.nl> | 2014-10-01 16:52:21 +0200 | 
| commit | 9831569ef1730b1003f8ebe4378ce31973fcdf9f (patch) | |
| tree | 6b6d7a26cd74564631b0268b9ad58cbb6c0023e6 | |
| parent | 78a179737366188aea505173e3a935b08c3cef55 (diff) | |
| download | astra-9831569ef1730b1003f8ebe4378ce31973fcdf9f.tar.gz astra-9831569ef1730b1003f8ebe4378ce31973fcdf9f.tar.bz2 astra-9831569ef1730b1003f8ebe4378ce31973fcdf9f.tar.xz astra-9831569ef1730b1003f8ebe4378ce31973fcdf9f.zip | |
Support flexible 2D volume geometry
| -rw-r--r-- | cuda/2d/astra.cu | 266 | ||||
| -rw-r--r-- | cuda/2d/astra.h | 36 | ||||
| -rw-r--r-- | src/CudaForwardProjectionAlgorithm.cpp | 61 | 
3 files changed, 231 insertions, 132 deletions
| diff --git a/cuda/2d/astra.cu b/cuda/2d/astra.cu index f4d4717..087905d 100644 --- a/cuda/2d/astra.cu +++ b/cuda/2d/astra.cu @@ -43,6 +43,10 @@ $Id$  #include <cuda.h>  #include "../../include/astra/Logger.h" +#include "../../include/astra/VolumeGeometry2D.h" +#include "../../include/astra/ParallelProjectionGeometry2D.h" +#include "../../include/astra/FanFlatProjectionGeometry2D.h" +#include "../../include/astra/FanFlatVecProjectionGeometry2D.h"  // For fan beam FBP weighting  #include "../3d/fdk.h" @@ -628,7 +632,7 @@ bool astraCudaFP(const float* pfVolume, float* pfSinogram,                   unsigned int iProjAngles, unsigned int iProjDets,                   const float *pfAngles, const float *pfOffsets,                   float fDetSize, unsigned int iDetSuperSampling, -                 int iGPUIndex) +                 float fOutputScale, int iGPUIndex)  {  	SDimensions dims; @@ -687,7 +691,7 @@ bool astraCudaFP(const float* pfVolume, float* pfSinogram,  	}  	zeroProjectionData(D_sinoData, sinoPitch, dims); -	ok = FP(D_volumeData, volumePitch, D_sinoData, sinoPitch, dims, pfAngles, pfOffsets, 1.0f); +	ok = FP(D_volumeData, volumePitch, D_sinoData, sinoPitch, dims, pfAngles, pfOffsets, fOutputScale);  	if (!ok) {  		cudaFree(D_volumeData);  		cudaFree(D_sinoData); @@ -711,19 +715,18 @@ bool astraCudaFP(const float* pfVolume, float* pfSinogram,  bool astraCudaFanFP(const float* pfVolume, float* pfSinogram,                      unsigned int iVolWidth, unsigned int iVolHeight,                      unsigned int iProjAngles, unsigned int iProjDets, -                    const float *pfAngles, float fOriginSourceDistance, -                    float fOriginDetectorDistance, float fPixelSize, -                    float fDetSize, -                    unsigned int iDetSuperSampling, +                    const SFanProjection *pAngles, +                    unsigned int iDetSuperSampling, float fOutputScale,                      int iGPUIndex)  {  	SDimensions dims; -	if (iProjAngles == 0 || iProjDets == 0 || pfAngles == 0) +	if (iProjAngles == 0 || iProjDets == 0 || pAngles == 0)  		return false;  	dims.iProjAngles = iProjAngles;  	dims.iProjDets = iProjDets; +	dims.fDetScale = 1.0f; // TODO?  	if (iDetSuperSampling == 0)  		return false; @@ -774,27 +777,7 @@ bool astraCudaFanFP(const float* pfVolume, float* pfSinogram,  	zeroProjectionData(D_sinoData, sinoPitch, dims); -	// TODO: Turn this geometry conversion into a util function -	SFanProjection* projs = new SFanProjection[dims.iProjAngles]; - -	float fSrcX0 = 0.0f; -	float fSrcY0 = -fOriginSourceDistance / fPixelSize; -	float fDetUX0 = fDetSize / fPixelSize; -	float fDetUY0 = 0.0f; -	float fDetSX0 = dims.iProjDets * fDetUX0 / -2.0f; -	float fDetSY0 = fOriginDetectorDistance / fPixelSize; - -#define ROTATE0(name,i,alpha) do { projs[i].f##name##X = f##name##X0 * cos(alpha) - f##name##Y0 * sin(alpha); projs[i].f##name##Y = f##name##X0 * sin(alpha) + f##name##Y0 * cos(alpha); } while(0) -	for (int i = 0; i < dims.iProjAngles; ++i) { -		ROTATE0(Src, i, pfAngles[i]); -		ROTATE0(DetS, i, pfAngles[i]); -		ROTATE0(DetU, i, pfAngles[i]); -	} - -#undef ROTATE0 - -	ok = FanFP(D_volumeData, volumePitch, D_sinoData, sinoPitch, dims, projs, 1.0f); -	delete[] projs; +	ok = FanFP(D_volumeData, volumePitch, D_sinoData, sinoPitch, dims, pAngles, fOutputScale);  	if (!ok) {  		cudaFree(D_volumeData); @@ -819,94 +802,205 @@ bool astraCudaFanFP(const float* pfVolume, float* pfSinogram,  } -bool astraCudaFanFP(const float* pfVolume, float* pfSinogram, -                    unsigned int iVolWidth, unsigned int iVolHeight, -                    unsigned int iProjAngles, unsigned int iProjDets, -                    const SFanProjection *pAngles, -                    unsigned int iDetSuperSampling, -                    int iGPUIndex) +bool convertAstraGeometry(const CVolumeGeometry2D* pVolGeom, +                          const CParallelProjectionGeometry2D* pProjGeom, +                          float*& detectorOffsets, float*& projectionAngles, +                          float& detSize, float& outputScale)  { -	SDimensions dims; +	assert(pVolGeom); +	assert(pProjGeom); +	assert(pProjGeom->getProjectionAngles()); -	if (iProjAngles == 0 || iProjDets == 0 || pAngles == 0) -		return false; +	const float EPS = 0.00001f; -	dims.iProjAngles = iProjAngles; -	dims.iProjDets = iProjDets; -	dims.fDetScale = 1.0f; // TODO? +	int nth = pProjGeom->getProjectionAngleCount(); -	if (iDetSuperSampling == 0) +	// Check if pixels are square +	if (abs(pVolGeom->getPixelLengthX() - pVolGeom->getPixelLengthY()) > EPS)  		return false; -	dims.iRaysPerDet = iDetSuperSampling; -	if (iVolWidth <= 0 || iVolHeight <= 0) -		return false; +	// Scale volume pixels to 1x1 +	detSize = pProjGeom->getDetectorWidth() / pVolGeom->getPixelLengthX(); -	dims.iVolWidth = iVolWidth; -	dims.iVolHeight = iVolHeight; - -	if (iGPUIndex != -1) { -		cudaSetDevice(iGPUIndex); -		cudaError_t err = cudaGetLastError(); +	// Copy angles +	float *angles = new float[nth]; +	for (int i = 0; i < nth; ++i) +		angles[i] = pProjGeom->getProjectionAngles()[i]; +	projectionAngles = angles; -		// Ignore errors caused by calling cudaSetDevice multiple times -		if (err != cudaSuccess && err != cudaErrorSetOnActiveProcess) -			return false; +	// Check if we need to translate +	bool offCenter = false; +	if (abs(pVolGeom->getWindowMinX() + pVolGeom->getWindowMaxX()) > EPS || +	    abs(pVolGeom->getWindowMinY() + pVolGeom->getWindowMaxY()) > EPS) +	{ +		offCenter = true;  	} -	bool ok; +	// If there are existing detector offsets, or if we need to translate, +	// we need to return offsets +	if (pProjGeom->getExtraDetectorOffset() || offCenter) +	{ +		float* offset = new float[nth]; + +		if (pProjGeom->getExtraDetectorOffset()) { +			for (int i = 0; i < nth; ++i) +				offset[i] = pProjGeom->getExtraDetectorOffset()[i]; +		} else { +			for (int i = 0; i < nth; ++i) +				offset[i] = 0.0f; +		} -	float* D_volumeData; -	unsigned int volumePitch; +		if (offCenter) { +			float dx = (pVolGeom->getWindowMinX() + pVolGeom->getWindowMaxX()) / 2; +			float dy = (pVolGeom->getWindowMinY() + pVolGeom->getWindowMaxY()) / 2; -	ok = allocateVolumeData(D_volumeData, volumePitch, dims); -	if (!ok) -		return false; +			// CHECKME: Is d in pixels or in units? -	float* D_sinoData; -	unsigned int sinoPitch; +			for (int i = 0; i < nth; ++i) { +				float d = dx * cos(angles[i]) + dy * sin(angles[i]); +				offset[i] += d; +			} +		} -	ok = allocateProjectionData(D_sinoData, sinoPitch, dims); -	if (!ok) { -		cudaFree(D_volumeData); -		return false; +		// CHECKME: Order of scaling and translation + +		// Scale volume pixels to 1x1 +		for (int i = 0; i < nth; ++i) { +			//offset[i] /= pVolGeom->getPixelLengthX(); +			//offset[i] *= detSize; +		} + + +		detectorOffsets = offset; +	} else { +		detectorOffsets = 0;  	} -	ok = copyVolumeToDevice(pfVolume, dims.iVolWidth, -	                        dims, -	                        D_volumeData, volumePitch); -	if (!ok) { -		cudaFree(D_volumeData); -		cudaFree(D_sinoData); -		return false; +	outputScale = pVolGeom->getPixelLengthX(); +	outputScale *= outputScale; + +	return true; +} + +static void convertAstraGeometry_internal(const CVolumeGeometry2D* pVolGeom, +                          unsigned int iProjectionAngleCount, +                          astraCUDA::SFanProjection*& pProjs, +                          float& outputScale) +{ +	// Translate +	float dx = (pVolGeom->getWindowMinX() + pVolGeom->getWindowMaxX()) / 2; +	float dy = (pVolGeom->getWindowMinY() + pVolGeom->getWindowMaxY()) / 2; + +	for (int i = 0; i < iProjectionAngleCount; ++i) { +		pProjs[i].fSrcX -= dx; +		pProjs[i].fSrcY -= dy; +		pProjs[i].fDetSX -= dx; +		pProjs[i].fDetSY -= dy;  	} -	zeroProjectionData(D_sinoData, sinoPitch, dims); +	// CHECKME: Order of scaling and translation -	ok = FanFP(D_volumeData, volumePitch, D_sinoData, sinoPitch, dims, pAngles, 1.0f); +	// Scale +	float factor = 1.0f / pVolGeom->getPixelLengthX(); +	for (int i = 0; i < iProjectionAngleCount; ++i) { +		pProjs[i].fSrcX *= factor; +		pProjs[i].fSrcY *= factor; +		pProjs[i].fDetSX *= factor; +		pProjs[i].fDetSY *= factor; +		pProjs[i].fDetUX *= factor; +		pProjs[i].fDetUY *= factor; -	if (!ok) { -		cudaFree(D_volumeData); -		cudaFree(D_sinoData); -		return false;  	} -	ok = copySinogramFromDevice(pfSinogram, dims.iProjDets, -	                            dims, -	                            D_sinoData, sinoPitch); -	if (!ok) { -		cudaFree(D_volumeData); -		cudaFree(D_sinoData); +	// CHECKME: Check factor +	outputScale = pVolGeom->getPixelLengthX(); +//	outputScale *= outputScale; +} + + +bool convertAstraGeometry(const CVolumeGeometry2D* pVolGeom, +                          const CFanFlatProjectionGeometry2D* pProjGeom, +                          astraCUDA::SFanProjection*& pProjs, +                          float& outputScale) +{ +	assert(pVolGeom); +	assert(pProjGeom); +	assert(pProjGeom->getProjectionAngles()); + +	const float EPS = 0.00001f; + +	int nth = pProjGeom->getProjectionAngleCount(); + +	// Check if pixels are square +	if (abs(pVolGeom->getPixelLengthX() - pVolGeom->getPixelLengthY()) > EPS)  		return false; + +	// TODO: Deprecate this. +//	if (pProjGeom->getExtraDetectorOffset()) +//		return false; + + +	float fOriginSourceDistance = pProjGeom->getOriginSourceDistance(); +	float fOriginDetectorDistance = pProjGeom->getOriginDetectorDistance(); +	float fDetSize = pProjGeom->getDetectorWidth(); +	const float *pfAngles = pProjGeom->getProjectionAngles(); + +	pProjs = new SFanProjection[nth]; + +	float fSrcX0 = 0.0f; +	float fSrcY0 = -fOriginSourceDistance; +	float fDetUX0 = fDetSize; +	float fDetUY0 = 0.0f; +	float fDetSX0 = pProjGeom->getDetectorCount() * fDetUX0 / -2.0f; +	float fDetSY0 = fOriginDetectorDistance; + +#define ROTATE0(name,i,alpha) do { pProjs[i].f##name##X = f##name##X0 * cos(alpha) - f##name##Y0 * sin(alpha); pProjs[i].f##name##Y = f##name##X0 * sin(alpha) + f##name##Y0 * cos(alpha); } while(0) +	for (int i = 0; i < nth; ++i) { +		ROTATE0(Src, i, pfAngles[i]); +		ROTATE0(DetS, i, pfAngles[i]); +		ROTATE0(DetU, i, pfAngles[i]);  	} -	cudaFree(D_volumeData); -	cudaFree(D_sinoData); +#undef ROTATE0 + +	convertAstraGeometry_internal(pVolGeom, nth, pProjs, outputScale);  	return true;  } +bool convertAstraGeometry(const CVolumeGeometry2D* pVolGeom, +                          const CFanFlatVecProjectionGeometry2D* pProjGeom, +                          astraCUDA::SFanProjection*& pProjs, +                          float& outputScale) +{ +	assert(pVolGeom); +	assert(pProjGeom); +	assert(pProjGeom->getProjectionVectors()); + +	const float EPS = 0.00001f; + +	int nx = pVolGeom->getGridColCount(); +	int ny = pVolGeom->getGridRowCount(); +	int nth = pProjGeom->getProjectionAngleCount(); + +	// Check if pixels are square +	if (abs(pVolGeom->getPixelLengthX() - pVolGeom->getPixelLengthY()) > EPS) +		return false; + +	pProjs = new SFanProjection[nth]; + +	// Copy vectors +	for (int i = 0; i < nth; ++i) +		pProjs[i] = pProjGeom->getProjectionVectors()[i]; + +	convertAstraGeometry_internal(pVolGeom, nth, pProjs, outputScale); + +	return true; +} + + +  } diff --git a/cuda/2d/astra.h b/cuda/2d/astra.h index 42b15e8..d2d7059 100644 --- a/cuda/2d/astra.h +++ b/cuda/2d/astra.h @@ -42,6 +42,11 @@ enum Cuda2DProjectionKernel {  	ker2d_default = 0  }; +class CParallelProjectionGeometry2D; +class CFanFlatProjectionGeometry2D; +class CFanFlatVecProjectionGeometry2D; +class CVolumeGeometry2D; +  class AstraFBP_internal;  class _AstraExport AstraFBP { @@ -195,24 +200,31 @@ _AstraExport bool astraCudaFP(const float* pfVolume, float* pfSinogram,                   unsigned int iProjAngles, unsigned int iProjDets,                   const float *pfAngles, const float *pfOffsets,                   float fDetSize = 1.0f, unsigned int iDetSuperSampling = 1, -                 int iGPUIndex = 0); - -// Do a single forward projection, fan beam -_AstraExport bool astraCudaFanFP(const float* pfVolume, float* pfSinogram, -                    unsigned int iVolWidth, unsigned int iVolHeight, -                    unsigned int iProjAngles, unsigned int iProjDets, -                    const float *pfAngles, float fOriginSourceDistance, -                    float fOriginDetectorDistance, float fPixelSize = 1.0f, -                    float fDetSize = 1.0f, -                    unsigned int iDetSuperSampling = 1, -                    int iGPUIndex = 0); +                 float fOutputScale = 1.0f, int iGPUIndex = 0);  _AstraExport bool astraCudaFanFP(const float* pfVolume, float* pfSinogram,                      unsigned int iVolWidth, unsigned int iVolHeight,                      unsigned int iProjAngles, unsigned int iProjDets,                      const SFanProjection *pAngles,                      unsigned int iDetSuperSampling = 1, -                    int iGPUIndex = 0); +                    float fOutputScale = 1.0f, int iGPUIndex = 0); + + +_AstraExport bool convertAstraGeometry(const CVolumeGeometry2D* pVolGeom, +                    const CParallelProjectionGeometry2D* pProjGeom, +                    float*& pfDetectorOffsets, float*& pfProjectionAngles, +                    float& fDetSize, float& fOutputScale); + +_AstraExport bool convertAstraGeometry(const CVolumeGeometry2D* pVolGeom, +                    const CFanFlatProjectionGeometry2D* pProjGeom, +                    astraCUDA::SFanProjection*& pProjs, +                    float& outputScale); + +_AstraExport bool convertAstraGeometry(const CVolumeGeometry2D* pVolGeom, +                    const CFanFlatVecProjectionGeometry2D* pProjGeom, +                    astraCUDA::SFanProjection*& pProjs, +                    float& outputScale); +  }  #endif diff --git a/src/CudaForwardProjectionAlgorithm.cpp b/src/CudaForwardProjectionAlgorithm.cpp index bd9cffb..b182251 100644 --- a/src/CudaForwardProjectionAlgorithm.cpp +++ b/src/CudaForwardProjectionAlgorithm.cpp @@ -214,52 +214,45 @@ void CCudaForwardProjectionAlgorithm::run(int)  	bool ok = false;  	if (parProjGeom) { + +		float *offsets, *angles, detSize, outputScale; +		ok = convertAstraGeometry(pVolGeom, parProjGeom, offsets, angles, detSize, outputScale); + +		ASTRA_ASSERT(ok); // FIXME + +		// FIXME: Output scaling +  		ok = astraCudaFP(m_pVolume->getDataConst(), m_pSinogram->getData(),  		                 pVolGeom->getGridColCount(), pVolGeom->getGridRowCount(),  		                 parProjGeom->getProjectionAngleCount(),  		                 parProjGeom->getDetectorCount(), -		                 parProjGeom->getProjectionAngles(), -		                 parProjGeom->getExtraDetectorOffset(), parProjGeom->getDetectorWidth() / pVolGeom->getPixelLengthX(), -		                 m_iDetectorSuperSampling, m_iGPUIndex); +		                 angles, offsets, detSize, +		                 m_iDetectorSuperSampling, 1.0f * outputScale, m_iGPUIndex); -	} else if (fanProjGeom) { +		delete[] offsets; +		delete[] angles; -		ok = astraCudaFanFP(m_pVolume->getDataConst(), m_pSinogram->getData(), -		                    pVolGeom->getGridColCount(), pVolGeom->getGridRowCount(), -		                    fanProjGeom->getProjectionAngleCount(), -		                    fanProjGeom->getDetectorCount(), -		                    fanProjGeom->getProjectionAngles(), -		                    fanProjGeom->getOriginSourceDistance(), -		                    fanProjGeom->getOriginDetectorDistance(), -		                    pVolGeom->getPixelLengthX(), -		                    fanProjGeom->getDetectorWidth(), -		                    m_iDetectorSuperSampling, m_iGPUIndex); - -	} else if (fanVecProjGeom) { - -		// Rescale projs to fPixelSize == 1 -		float fPixelSize = pVolGeom->getPixelLengthX(); -		const astraCUDA::SFanProjection* projs; -		projs = fanVecProjGeom->getProjectionVectors(); - -		astraCUDA::SFanProjection* scaledProjs = new astraCUDA::SFanProjection[fanVecProjGeom->getProjectionAngleCount()]; -#define SCALE(name,i,alpha) do { scaledProjs[i].f##name##X = projs[i].f##name##X * alpha; scaledProjs[i].f##name##Y = projs[i].f##name##Y * alpha; } while (0) -		for (unsigned int i = 0; i < fanVecProjGeom->getProjectionAngleCount(); ++i) { -			SCALE(Src,i,1.0f/fPixelSize); -			SCALE(DetS,i,1.0f/fPixelSize); -			SCALE(DetU,i,1.0f/fPixelSize); +	} else if (fanProjGeom || fanVecProjGeom) { + +		astraCUDA::SFanProjection* projs; +		float outputScale; + +		if (fanProjGeom) { +			ok = convertAstraGeometry(pVolGeom, fanProjGeom, projs, outputScale); +		} else { +			ok = convertAstraGeometry(pVolGeom, fanVecProjGeom, projs, outputScale);  		} +		ASTRA_ASSERT(ok);  		ok = astraCudaFanFP(m_pVolume->getDataConst(), m_pSinogram->getData(),  		                    pVolGeom->getGridColCount(), pVolGeom->getGridRowCount(), -		                    fanVecProjGeom->getProjectionAngleCount(), -		                    fanVecProjGeom->getDetectorCount(), -		                    scaledProjs, -		                    /* 1.0f / pVolGeom->getPixelLengthX(), */ -		                    m_iDetectorSuperSampling, m_iGPUIndex); +		                    m_pSinogram->getGeometry()->getProjectionAngleCount(), +		                    m_pSinogram->getGeometry()->getDetectorCount(), +		                    projs, +		                    m_iDetectorSuperSampling, outputScale, m_iGPUIndex); -		delete[] scaledProjs; +		delete[] projs;  	} else { | 
