summaryrefslogtreecommitdiffstats
path: root/Wrappers/Matlab
diff options
context:
space:
mode:
authorDaniil Kazantsev <dkazanc3@googlemail.com>2019-02-19 22:14:48 +0000
committerGitHub <noreply@github.com>2019-02-19 22:14:48 +0000
commit8e71dd67abef0caddb24caa365321d3874529254 (patch)
tree108c7e25e3d741ca04ef45aa29eb61a9732075f4 /Wrappers/Matlab
parent8f2e86726669b9dadb3c788e0ea681d397a2eeb7 (diff)
parent53d5508915709245d0573e0335de83fc24313b5a (diff)
downloadregularization-8e71dd67abef0caddb24caa365321d3874529254.tar.gz
regularization-8e71dd67abef0caddb24caa365321d3874529254.tar.bz2
regularization-8e71dd67abef0caddb24caa365321d3874529254.tar.xz
regularization-8e71dd67abef0caddb24caa365321d3874529254.zip
Merge pull request #99 from vais-ral/TGV_bugfix
TGV bugfix
Diffstat (limited to 'Wrappers/Matlab')
-rw-r--r--Wrappers/Matlab/demos/demoMatlab_3Ddenoise.m16
-rw-r--r--Wrappers/Matlab/demos/demoMatlab_denoise.m16
-rw-r--r--Wrappers/Matlab/mex_compile/regularisers_GPU/TGV_GPU.cpp32
3 files changed, 38 insertions, 26 deletions
diff --git a/Wrappers/Matlab/demos/demoMatlab_3Ddenoise.m b/Wrappers/Matlab/demos/demoMatlab_3Ddenoise.m
index 0c331a4..ac8e1ba 100644
--- a/Wrappers/Matlab/demos/demoMatlab_3Ddenoise.m
+++ b/Wrappers/Matlab/demos/demoMatlab_3Ddenoise.m
@@ -8,7 +8,7 @@ addpath(Path2);
addpath(Path3);
N = 512;
-slices = 7;
+slices = 15;
vol3D = zeros(N,N,slices, 'single');
Ideal3D = zeros(N,N,slices, 'single');
Im = double(imread('lena_gray_512.tif'))/255; % loading image
@@ -17,9 +17,7 @@ vol3D(:,:,i) = Im + .05*randn(size(Im));
Ideal3D(:,:,i) = Im;
end
vol3D(vol3D < 0) = 0;
-figure; imshow(vol3D(:,:,15), [0 1]); title('Noisy image');
-
-
+figure; imshow(vol3D(:,:,7), [0 1]); title('Noisy image');
lambda_reg = 0.03; % regularsation parameter for all methods
%%
fprintf('Denoise a volume using the ROF-TV model (CPU) \n');
@@ -143,6 +141,16 @@ rmseTGV = RMSE(Ideal3D(:),u_tgv(:));
fprintf('%s %f \n', 'RMSE error for TGV is:', rmseTGV);
figure; imshow(u_tgv(:,:,3), [0 1]); title('TGV denoised volume (CPU)');
%%
+% fprintf('Denoise using the TGV model (GPU) \n');
+% lambda_TGV = 0.03; % regularisation parameter
+% alpha1 = 1.0; % parameter to control the first-order term
+% alpha0 = 2.0; % parameter to control the second-order term
+% iter_TGV = 500; % number of Primal-Dual iterations for TGV
+% tic; u_tgv_gpu = TGV_GPU(single(vol3D), lambda_TGV, alpha1, alpha0, iter_TGV); toc;
+% rmseTGV = RMSE(Ideal3D(:),u_tgv_gpu(:));
+% fprintf('%s %f \n', 'RMSE error for TGV is:', rmseTGV);
+% figure; imshow(u_tgv_gpu(:,:,3), [0 1]); title('TGV denoised volume (GPU)');
+%%
%>>>>>>>>>>>>>> MULTI-CHANNEL priors <<<<<<<<<<<<<<< %
fprintf('Denoise a volume using the FGP-dTV model (CPU) \n');
diff --git a/Wrappers/Matlab/demos/demoMatlab_denoise.m b/Wrappers/Matlab/demos/demoMatlab_denoise.m
index 14d3096..62e5834 100644
--- a/Wrappers/Matlab/demos/demoMatlab_denoise.m
+++ b/Wrappers/Matlab/demos/demoMatlab_denoise.m
@@ -5,7 +5,9 @@ fsep = '/';
Path1 = sprintf(['..' fsep 'mex_compile' fsep 'installed'], 1i);
Path2 = sprintf(['..' fsep '..' fsep '..' fsep 'data' fsep], 1i);
Path3 = sprintf(['..' fsep 'supp'], 1i);
-addpath(Path1); addpath(Path2); addpath(Path3);
+addpath(Path1);
+addpath(Path2);
+addpath(Path3);
Im = double(imread('lena_gray_512.tif'))/255; % loading image
u0 = Im + .05*randn(size(Im)); u0(u0 < 0) = 0;
@@ -29,7 +31,7 @@ figure; imshow(u_rof, [0 1]); title('ROF-TV denoised image (CPU)');
% figure; imshow(u_rofG, [0 1]); title('ROF-TV denoised image (GPU)');
%%
fprintf('Denoise using the FGP-TV model (CPU) \n');
-iter_fgp = 1000; % number of FGP iterations
+iter_fgp = 1300; % number of FGP iterations
epsil_tol = 1.0e-06; % tolerance
tic; u_fgp = FGP_TV(single(u0), lambda_reg, iter_fgp, epsil_tol); toc;
energyfunc_val_fgp = TV_energy(single(u_fgp),single(u0),lambda_reg, 1); % get energy function value
@@ -39,8 +41,8 @@ figure; imshow(u_fgp, [0 1]); title('FGP-TV denoised image (CPU)');
%%
% fprintf('Denoise using the FGP-TV model (GPU) \n');
-% iter_fgp = 1000; % number of FGP iterations
-% epsil_tol = 1.0e-05; % tolerance
+% iter_fgp = 1300; % number of FGP iterations
+% epsil_tol = 1.0e-06; % tolerance
% tic; u_fgpG = FGP_TV_GPU(single(u0), lambda_reg, iter_fgp, epsil_tol); toc;
% figure; imshow(u_fgpG, [0 1]); title('FGP-TV denoised image (GPU)');
%%
@@ -63,17 +65,17 @@ fprintf('Denoise using the TGV model (CPU) \n');
lambda_TGV = 0.045; % regularisation parameter
alpha1 = 1.0; % parameter to control the first-order term
alpha0 = 2.0; % parameter to control the second-order term
-iter_TGV = 2000; % number of Primal-Dual iterations for TGV
+iter_TGV = 1500; % number of Primal-Dual iterations for TGV
tic; u_tgv = TGV(single(u0), lambda_TGV, alpha1, alpha0, iter_TGV); toc;
rmseTGV = (RMSE(u_tgv(:),Im(:)));
fprintf('%s %f \n', 'RMSE error for TGV is:', rmseTGV);
figure; imshow(u_tgv, [0 1]); title('TGV denoised image (CPU)');
-%%
+
% fprintf('Denoise using the TGV model (GPU) \n');
% lambda_TGV = 0.045; % regularisation parameter
% alpha1 = 1.0; % parameter to control the first-order term
% alpha0 = 2.0; % parameter to control the second-order term
-% iter_TGV = 2000; % number of Primal-Dual iterations for TGV
+% iter_TGV = 1500; % number of Primal-Dual iterations for TGV
% tic; u_tgv_gpu = TGV_GPU(single(u0), lambda_TGV, alpha1, alpha0, iter_TGV); toc;
% rmseTGV_gpu = (RMSE(u_tgv_gpu(:),Im(:)));
% fprintf('%s %f \n', 'RMSE error for TGV is:', rmseTGV_gpu);
diff --git a/Wrappers/Matlab/mex_compile/regularisers_GPU/TGV_GPU.cpp b/Wrappers/Matlab/mex_compile/regularisers_GPU/TGV_GPU.cpp
index edb551d..1173282 100644
--- a/Wrappers/Matlab/mex_compile/regularisers_GPU/TGV_GPU.cpp
+++ b/Wrappers/Matlab/mex_compile/regularisers_GPU/TGV_GPU.cpp
@@ -21,18 +21,18 @@ limitations under the License.
#include "TGV_GPU_core.h"
/* CUDA implementation of Primal-Dual denoising method for
- * Total Generilized Variation (TGV)-L2 model [1] (2D case only)
+ * Total Generilized Variation (TGV)-L2 model [1] (2D/3D)
*
* Input Parameters:
- * 1. Noisy image (2D) (required)
- * 2. lambda - regularisation parameter (required)
- * 3. parameter to control the first-order term (alpha1) (default - 1)
- * 4. parameter to control the second-order term (alpha0) (default - 0.5)
- * 5. Number of Chambolle-Pock (Primal-Dual) iterations (default is 300)
+ * 1. Noisy image/volume (2D/3D)
+ * 2. lambda - regularisation parameter
+ * 3. parameter to control the first-order term (alpha1)
+ * 4. parameter to control the second-order term (alpha0)
+ * 5. Number of Chambolle-Pock (Primal-Dual) iterations
* 6. Lipshitz constant (default is 12)
*
* Output:
- * Filtered/regulariaed image
+ * Filtered/regularised image
*
* References:
* [1] K. Bredies "Total Generalized Variation"
@@ -44,7 +44,7 @@ void mexFunction(
{
int number_of_dims, iter;
- mwSize dimX, dimY;
+ mwSize dimX, dimY, dimZ;
const mwSize *dim_array;
float *Input, *Output=NULL, lambda, alpha0, alpha1, L2;
@@ -57,8 +57,8 @@ void mexFunction(
Input = (float *) mxGetData(prhs[0]); /*noisy image (2D) */
lambda = (float) mxGetScalar(prhs[1]); /* regularisation parameter */
alpha1 = 1.0f; /* parameter to control the first-order term */
- alpha0 = 0.5f; /* parameter to control the second-order term */
- iter = 300; /* Iterations number */
+ alpha0 = 2.0f; /* parameter to control the second-order term */
+ iter = 500; /* Iterations number */
L2 = 12.0f; /* Lipshitz constant */
if (mxGetClassID(prhs[0]) != mxSINGLE_CLASS) {mexErrMsgTxt("The input image must be in a single precision"); }
@@ -68,12 +68,14 @@ void mexFunction(
if (nrhs == 6) L2 = (float) mxGetScalar(prhs[5]); /* Lipshitz constant */
/*Handling Matlab output data*/
- dimX = dim_array[0]; dimY = dim_array[1];
+ dimX = dim_array[0]; dimY = dim_array[1]; dimZ = dim_array[2];
if (number_of_dims == 2) {
- Output = (float*)mxGetPr(plhs[0] = mxCreateNumericArray(2, dim_array, mxSINGLE_CLASS, mxREAL));
- /* running the function */
- TGV_GPU_main(Input, Output, lambda, alpha1, alpha0, iter, L2, dimX, dimY);
+ dimZ = 1; /*2D case*/
+ Output = (float*)mxGetPr(plhs[0] = mxCreateNumericArray(2, dim_array, mxSINGLE_CLASS, mxREAL));
}
- if (number_of_dims == 3) {mexErrMsgTxt("Only 2D images accepted");}
+ if (number_of_dims == 3) Output = (float*)mxGetPr(plhs[0] = mxCreateNumericArray(3, dim_array, mxSINGLE_CLASS, mxREAL));
+
+ /* running the function */
+ TGV_GPU_main(Input, Output, lambda, alpha1, alpha0, iter, L2, dimX, dimY, dimZ);
}