summaryrefslogtreecommitdiffstats
path: root/Wrappers/Matlab
diff options
context:
space:
mode:
Diffstat (limited to 'Wrappers/Matlab')
-rw-r--r--Wrappers/Matlab/demos/demoMatlab_3Ddenoise.m45
-rw-r--r--Wrappers/Matlab/demos/demoMatlab_denoise.m37
-rw-r--r--Wrappers/Matlab/mex_compile/compileCPU_mex.m5
-rw-r--r--Wrappers/Matlab/mex_compile/compileGPU_mex.m6
-rw-r--r--Wrappers/Matlab/mex_compile/regularisers_CPU/FGP_TV.c2
-rw-r--r--Wrappers/Matlab/mex_compile/regularisers_CPU/FGP_TV.c~91
-rw-r--r--Wrappers/Matlab/mex_compile/regularisers_CPU/FGP_dTV.c113
-rw-r--r--Wrappers/Matlab/mex_compile/regularisers_GPU/FGP_dTV_GPU.cpp111
8 files changed, 305 insertions, 105 deletions
diff --git a/Wrappers/Matlab/demos/demoMatlab_3Ddenoise.m b/Wrappers/Matlab/demos/demoMatlab_3Ddenoise.m
index 71082e7..dc49d9c 100644
--- a/Wrappers/Matlab/demos/demoMatlab_3Ddenoise.m
+++ b/Wrappers/Matlab/demos/demoMatlab_3Ddenoise.m
@@ -1,5 +1,6 @@
% Volume (3D) denoising demo using CCPi-RGL
-
+clear
+close all
addpath('../mex_compile/installed');
addpath('../../../data/');
@@ -14,31 +15,65 @@ vol3D(vol3D < 0) = 0;
figure; imshow(vol3D(:,:,15), [0 1]); title('Noisy image');
%%
-fprintf('Denoise using ROF-TV model (CPU) \n');
+fprintf('Denoise a volume using the ROF-TV model (CPU) \n');
lambda_rof = 0.03; % regularisation parameter
tau_rof = 0.0025; % time-marching constant
iter_rof = 300; % number of ROF iterations
tic; u_rof = ROF_TV(single(vol3D), lambda_rof, iter_rof, tau_rof); toc;
figure; imshow(u_rof(:,:,15), [0 1]); title('ROF-TV denoised volume (CPU)');
%%
-% fprintf('Denoise using ROF-TV model (GPU) \n');
+% fprintf('Denoise a volume using the ROF-TV model (GPU) \n');
% lambda_rof = 0.03; % regularisation parameter
% tau_rof = 0.0025; % time-marching constant
% iter_rof = 300; % number of ROF iterations
% tic; u_rofG = ROF_TV_GPU(single(vol3D), lambda_rof, iter_rof, tau_rof); toc;
% figure; imshow(u_rofG(:,:,15), [0 1]); title('ROF-TV denoised volume (GPU)');
%%
-fprintf('Denoise using FGP-TV model (CPU) \n');
+fprintf('Denoise a volume using the FGP-TV model (CPU) \n');
lambda_fgp = 0.03; % regularisation parameter
iter_fgp = 300; % number of FGP iterations
epsil_tol = 1.0e-05; % tolerance
tic; u_fgp = FGP_TV(single(vol3D), lambda_fgp, iter_fgp, epsil_tol); toc;
figure; imshow(u_fgp(:,:,15), [0 1]); title('FGP-TV denoised volume (CPU)');
%%
-% fprintf('Denoise using FGP-TV model (GPU) \n');
+% fprintf('Denoise a volume using the FGP-TV model (GPU) \n');
% lambda_fgp = 0.03; % regularisation parameter
% iter_fgp = 300; % number of FGP iterations
% epsil_tol = 1.0e-05; % tolerance
% tic; u_fgpG = FGP_TV_GPU(single(vol3D), lambda_fgp, iter_fgp, epsil_tol); toc;
% figure; imshow(u_fgpG(:,:,15), [0 1]); title('FGP-TV denoised volume (GPU)');
%%
+fprintf('Denoise a volume using the FGP-dTV model (CPU) \n');
+
+% create another volume (reference) with slightly less amount of noise
+vol3D_ref = zeros(N,N,slices, 'single');
+for i = 1:slices
+vol3D_ref(:,:,i) = Im + .01*randn(size(Im));
+end
+vol3D_ref(vol3D_ref < 0) = 0;
+% vol3D_ref = zeros(size(Im),'single'); % pass zero reference (dTV -> TV)
+
+lambda_fgp = 0.03; % regularisation parameter
+iter_fgp = 300; % number of FGP iterations
+epsil_tol = 1.0e-05; % tolerance
+eta = 0.2; % Reference image gradient smoothing constant
+tic; u_fgp_dtv = FGP_dTV(single(vol3D), single(vol3D_ref), lambda_fgp, iter_fgp, epsil_tol, eta); toc;
+figure; imshow(u_fgp_dtv(:,:,15), [0 1]); title('FGP-dTV denoised volume (CPU)');
+%%
+fprintf('Denoise a volume using the FGP-dTV model (GPU) \n');
+
+% create another volume (reference) with slightly less amount of noise
+vol3D_ref = zeros(N,N,slices, 'single');
+for i = 1:slices
+vol3D_ref(:,:,i) = Im + .01*randn(size(Im));
+end
+vol3D_ref(vol3D_ref < 0) = 0;
+% vol3D_ref = zeros(size(Im),'single'); % pass zero reference (dTV -> TV)
+
+lambda_fgp = 0.03; % regularisation parameter
+iter_fgp = 300; % number of FGP iterations
+epsil_tol = 1.0e-05; % tolerance
+eta = 0.2; % Reference image gradient smoothing constant
+tic; u_fgp_dtv_g = FGP_dTV_GPU(single(vol3D), single(vol3D_ref), lambda_fgp, iter_fgp, epsil_tol, eta); toc;
+figure; imshow(u_fgp_dtv_g(:,:,15), [0 1]); title('FGP-dTV denoised volume (GPU)');
+%% \ No newline at end of file
diff --git a/Wrappers/Matlab/demos/demoMatlab_denoise.m b/Wrappers/Matlab/demos/demoMatlab_denoise.m
index 7f87fbb..145f2ff 100644
--- a/Wrappers/Matlab/demos/demoMatlab_denoise.m
+++ b/Wrappers/Matlab/demos/demoMatlab_denoise.m
@@ -1,5 +1,6 @@
% Image (2D) denoising demo using CCPi-RGL
-
+clear
+close all
addpath('../mex_compile/installed');
addpath('../../../data/');
@@ -8,31 +9,55 @@ u0 = Im + .05*randn(size(Im)); u0(u0 < 0) = 0;
figure; imshow(u0, [0 1]); title('Noisy image');
%%
-fprintf('Denoise using ROF-TV model (CPU) \n');
+fprintf('Denoise using the ROF-TV model (CPU) \n');
lambda_rof = 0.03; % regularisation parameter
tau_rof = 0.0025; % time-marching constant
iter_rof = 2000; % number of ROF iterations
tic; u_rof = ROF_TV(single(u0), lambda_rof, iter_rof, tau_rof); toc;
figure; imshow(u_rof, [0 1]); title('ROF-TV denoised image (CPU)');
%%
-% fprintf('Denoise using ROF-TV model (GPU) \n');
+% fprintf('Denoise using the ROF-TV model (GPU) \n');
% lambda_rof = 0.03; % regularisation parameter
% tau_rof = 0.0025; % time-marching constant
% iter_rof = 2000; % number of ROF iterations
% tic; u_rofG = ROF_TV_GPU(single(u0), lambda_rof, iter_rof, tau_rof); toc;
% figure; imshow(u_rofG, [0 1]); title('ROF-TV denoised image (GPU)');
%%
-fprintf('Denoise using FGP-TV model (CPU) \n');
+fprintf('Denoise using the FGP-TV model (CPU) \n');
lambda_fgp = 0.03; % regularisation parameter
iter_fgp = 1000; % number of FGP iterations
-epsil_tol = 1.0e-05; % tolerance
+epsil_tol = 1.0e-06; % tolerance
tic; u_fgp = FGP_TV(single(u0), lambda_fgp, iter_fgp, epsil_tol); toc;
figure; imshow(u_fgp, [0 1]); title('FGP-TV denoised image (CPU)');
%%
-% fprintf('Denoise using FGP-TV model (GPU) \n');
+% fprintf('Denoise using the FGP-TV model (GPU) \n');
% lambda_fgp = 0.03; % regularisation parameter
% iter_fgp = 1000; % number of FGP iterations
% epsil_tol = 1.0e-05; % tolerance
% tic; u_fgpG = FGP_TV_GPU(single(u0), lambda_fgp, iter_fgp, epsil_tol); toc;
% figure; imshow(u_fgpG, [0 1]); title('FGP-TV denoised image (GPU)');
%%
+fprintf('Denoise using the FGP-dTV model (CPU) \n');
+% create another image (reference) with slightly less amount of noise
+u_ref = Im + .01*randn(size(Im)); u_ref(u_ref < 0) = 0;
+% u_ref = zeros(size(Im),'single'); % pass zero reference (dTV -> TV)
+
+lambda_fgp = 0.03; % regularisation parameter
+iter_fgp = 1000; % number of FGP iterations
+epsil_tol = 1.0e-06; % tolerance
+eta = 0.2; % Reference image gradient smoothing constant
+tic; u_fgp_dtv = FGP_dTV(single(u0), single(u_ref), lambda_fgp, iter_fgp, epsil_tol, eta); toc;
+figure; imshow(u_fgp_dtv, [0 1]); title('FGP-dTV denoised image (CPU)');
+%%
+% fprintf('Denoise using the FGP-dTV model (GPU) \n');
+% % create another image (reference) with slightly less amount of noise
+% u_ref = Im + .01*randn(size(Im)); u_ref(u_ref < 0) = 0;
+% % u_ref = zeros(size(Im),'single'); % pass zero reference (dTV -> TV)
+%
+% lambda_fgp = 0.03; % regularisation parameter
+% iter_fgp = 1000; % number of FGP iterations
+% epsil_tol = 1.0e-06; % tolerance
+% eta = 0.2; % Reference image gradient smoothing constant
+% tic; u_fgp_dtvG = FGP_dTV_GPU(single(u0), single(u_ref), lambda_fgp, iter_fgp, epsil_tol, eta); toc;
+% figure; imshow(u_fgp_dtvG, [0 1]); title('FGP-dTV denoised image (GPU)');
+%%
diff --git a/Wrappers/Matlab/mex_compile/compileCPU_mex.m b/Wrappers/Matlab/mex_compile/compileCPU_mex.m
index 8da81ad..71f345a 100644
--- a/Wrappers/Matlab/mex_compile/compileCPU_mex.m
+++ b/Wrappers/Matlab/mex_compile/compileCPU_mex.m
@@ -11,7 +11,10 @@ movefile ROF_TV.mex* ../installed/
mex FGP_TV.c FGP_TV_core.c utils.c CFLAGS="\$CFLAGS -fopenmp -Wall -std=c99" LDFLAGS="\$LDFLAGS -fopenmp"
movefile FGP_TV.mex* ../installed/
-delete ROF_TV_core* FGP_TV_core* utils.c utils.h CCPiDefines.h
+mex FGP_dTV.c FGP_dTV_core.c utils.c CFLAGS="\$CFLAGS -fopenmp -Wall -std=c99" LDFLAGS="\$LDFLAGS -fopenmp"
+movefile FGP_dTV.mex* ../installed/
+
+delete ROF_TV_core* FGP_TV_core* FGP_dTV_core* utils* CCPiDefines.h
fprintf('%s \n', 'All successfully compiled!');
diff --git a/Wrappers/Matlab/mex_compile/compileGPU_mex.m b/Wrappers/Matlab/mex_compile/compileGPU_mex.m
index 45236fa..f58e9bc 100644
--- a/Wrappers/Matlab/mex_compile/compileGPU_mex.m
+++ b/Wrappers/Matlab/mex_compile/compileGPU_mex.m
@@ -23,7 +23,11 @@ movefile ROF_TV_GPU.mex* ../installed/
mex -g -I/usr/local/cuda-7.5/include -L/usr/local/cuda-7.5/lib64 -lcudart -lcufft -lmwgpu FGP_TV_GPU.cpp TV_FGP_GPU_core.o
movefile FGP_TV_GPU.mex* ../installed/
-delete TV_ROF_GPU_core* TV_FGP_GPU_core* CCPiDefines.h
+!/usr/local/cuda/bin/nvcc -O0 -c dTV_FGP_GPU_core.cu -Xcompiler -fPIC -I~/SOFT/MATLAB9/extern/include/
+mex -g -I/usr/local/cuda-7.5/include -L/usr/local/cuda-7.5/lib64 -lcudart -lcufft -lmwgpu FGP_dTV_GPU.cpp dTV_FGP_GPU_core.o
+movefile FGP_dTV_GPU.mex* ../installed/
+
+delete TV_ROF_GPU_core* TV_FGP_GPU_core* dTV_FGP_GPU_core* CCPiDefines.h
fprintf('%s \n', 'All successfully compiled!');
cd ../../
diff --git a/Wrappers/Matlab/mex_compile/regularisers_CPU/FGP_TV.c b/Wrappers/Matlab/mex_compile/regularisers_CPU/FGP_TV.c
index ba06cc7..aae1cb7 100644
--- a/Wrappers/Matlab/mex_compile/regularisers_CPU/FGP_TV.c
+++ b/Wrappers/Matlab/mex_compile/regularisers_CPU/FGP_TV.c
@@ -52,7 +52,7 @@ void mexFunction(
dim_array = mxGetDimensions(prhs[0]);
/*Handling Matlab input data*/
- if ((nrhs < 2) || (nrhs > 7)) mexErrMsgTxt("At least 2 parameters is required, all parameters are: Image(2D/3D), Regularization parameter. The full list of parameters: Image(2D/3D), Regularization parameter, iterations number, tolerance, penalty type ('iso' or 'l1'), nonnegativity switch, print switch");
+ if ((nrhs < 2) || (nrhs > 7)) mexErrMsgTxt("At least 2 parameters is required, all parameters are: Image(2D/3D), Regularization parameter, Regularization parameter, iterations number, tolerance, penalty type ('iso' or 'l1'), nonnegativity switch, print switch");
Input = (float *) mxGetData(prhs[0]); /*noisy image (2D/3D) */
lambda = (float) mxGetScalar(prhs[1]); /* regularization parameter */
diff --git a/Wrappers/Matlab/mex_compile/regularisers_CPU/FGP_TV.c~ b/Wrappers/Matlab/mex_compile/regularisers_CPU/FGP_TV.c~
deleted file mode 100644
index 30d61cd..0000000
--- a/Wrappers/Matlab/mex_compile/regularisers_CPU/FGP_TV.c~
+++ /dev/null
@@ -1,91 +0,0 @@
-/*
- * This work is part of the Core Imaging Library developed by
- * Visual Analytics and Imaging System Group of the Science Technology
- * Facilities Council, STFC
- *
- * Copyright 2017 Daniil Kazantsev
- * Copyright 2017 Srikanth Nagella, Edoardo Pasca
- *
- * Licensed under the Apache License, Version 2.0 (the "License");
- * you may not use this file except in compliance with the License.
- * You may obtain a copy of the License at
- * http://www.apache.org/licenses/LICENSE-2.0
- * Unless required by applicable law or agreed to in writing, software
- * distributed under the License is distributed on an "AS IS" BASIS,
- * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
- * See the License for the specific language governing permissions and
- * limitations under the License.
- */
-#include "matrix.h"
-#include "mex.h"
-#include "FGP_TV_core.h"
-
-/* C-OMP implementation of FGP-TV [1] denoising/regularization model (2D/3D case)
- *
- * Input Parameters:
- * 1. Noisy image/volume
- * 2. lambdaPar - regularization parameter
- * 3. Number of iterations
- * 4. eplsilon: tolerance constant
- * 5. TV-type: methodTV - 'iso' (0) or 'l1' (1)
- * 6. nonneg: 'nonnegativity (0 is OFF by default)
- * 7. print information: 0 (off) or 1 (on)
- *
- * Output:
- * [1] Filtered/regularized image
- *
- * This function is based on the Matlab's code and paper by
- * [1] Amir Beck and Marc Teboulle, "Fast Gradient-Based Algorithms for Constrained Total Variation Image Denoising and Deblurring Problems"
- */
-
-
-void mexFunction(
- int nlhs, mxArray *plhs[],
- int nrhs, const mxArray *prhs[])
-
-{
- int number_of_dims, iter, dimX, dimY, dimZ, methTV, printswitch;
- const int *dim_array;
- float *Input, *Output, lambda, epsil;
-
- number_of_dims = mxGetNumberOfDimensions(prhs[0]);
- dim_array = mxGetDimensions(prhs[0]);
-
- /*Handling Matlab input data*/
- if ((nrhs < 2) || (nrhs > 6)) mexErrMsgTxt("At least 2 parameters is required: Image(2D/3D), Regularization parameter. The full list of parameters: Image(2D/3D), Regularization parameter, iterations number, tolerance, penalty type ('iso' or 'l1'), print switch");
-
- Input = (float *) mxGetData(prhs[0]); /*noisy image (2D/3D) */
- lambda = (float) mxGetScalar(prhs[1]); /* regularization parameter */
- iter = 300; /* default iterations number */
- epsil = 0.0001; /* default tolerance constant */
- methTV = 0; /* default isotropic TV penalty */
- printswitch = 0; /*default print is switched off - 0 */
-
- if (mxGetClassID(prhs[0]) != mxSINGLE_CLASS) {mexErrMsgTxt("The input image must be in a single precision"); }
-
- if ((nrhs == 3) || (nrhs == 4) || (nrhs == 5) || (nrhs == 6)) iter = (int) mxGetScalar(prhs[2]); /* iterations number */
- if ((nrhs == 4) || (nrhs == 5) || (nrhs == 6)) epsil = (float) mxGetScalar(prhs[3]); /* tolerance constant */
- if ((nrhs == 5) || (nrhs == 6)) {
- char *penalty_type;
- penalty_type = mxArrayToString(prhs[4]); /* choosing TV penalty: 'iso' or 'l1', 'iso' is the default */
- if ((strcmp(penalty_type, "l1") != 0) && (strcmp(penalty_type, "iso") != 0)) mexErrMsgTxt("Choose TV type: 'iso' or 'l1',");
- if (strcmp(penalty_type, "l1") == 0) methTV = 1; /* enable 'l1' penalty */
- mxFree(penalty_type);
- }
- if (nrhs == 6) {
- printswitch = (int) mxGetScalar(prhs[5]);
- if ((printswitch != 0) || (printswitch != 1)) {mexErrMsgTxt("Print can be enabled by choosing 1 or off - 0"); }
- }
-
- /*Handling Matlab output data*/
- dimX = dim_array[0]; dimY = dim_array[1]; dimZ = dim_array[2];
-
- if (number_of_dims == 2) {
- dimZ = 1; /*2D case*/
- Output = (float*)mxGetPr(plhs[0] = mxCreateNumericArray(2, dim_array, mxSINGLE_CLASS, mxREAL));
- }
- if (number_of_dims == 3) Output = (float*)mxGetPr(plhs[0] = mxCreateNumericArray(3, dim_array, mxSINGLE_CLASS, mxREAL));
-
-
- TV_FGP_CPU_main(Input, Output, lambda, iter, epsil, methTV, nonneg, printswitch, dimX, dimY, dimZ)
-}
diff --git a/Wrappers/Matlab/mex_compile/regularisers_CPU/FGP_dTV.c b/Wrappers/Matlab/mex_compile/regularisers_CPU/FGP_dTV.c
new file mode 100644
index 0000000..bb868c7
--- /dev/null
+++ b/Wrappers/Matlab/mex_compile/regularisers_CPU/FGP_dTV.c
@@ -0,0 +1,113 @@
+/*
+ * This work is part of the Core Imaging Library developed by
+ * Visual Analytics and Imaging System Group of the Science Technology
+ * Facilities Council, STFC
+ *
+ * Copyright 2017 Daniil Kazantsev
+ * Copyright 2017 Srikanth Nagella, Edoardo Pasca
+ *
+ * Licensed under the Apache License, Version 2.0 (the "License");
+ * you may not use this file except in compliance with the License.
+ * You may obtain a copy of the License at
+ * http://www.apache.org/licenses/LICENSE-2.0
+ * Unless required by applicable law or agreed to in writing, software
+ * distributed under the License is distributed on an "AS IS" BASIS,
+ * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
+ * See the License for the specific language governing permissions and
+ * limitations under the License.
+ */
+#include "matrix.h"
+#include "mex.h"
+#include "FGP_dTV_core.h"
+
+/* C-OMP implementation of FGP-dTV [1,2] denoising/regularization model (2D/3D case)
+ * which employs structural similarity of the level sets of two images/volumes, see [1,2]
+ * The current implementation updates image 1 while image 2 is being fixed.
+ *
+ * Input Parameters:
+ * 1. Noisy image/volume [REQUIRED]
+ * 2. Additional reference image/volume of the same dimensions as (1) [REQUIRED]
+ * 3. lambdaPar - regularization parameter [REQUIRED]
+ * 4. Number of iterations [OPTIONAL]
+ * 5. eplsilon: tolerance constant [OPTIONAL]
+ * 6. eta: smoothing constant to calculate gradient of the reference [OPTIONAL] *
+ * 7. TV-type: methodTV - 'iso' (0) or 'l1' (1) [OPTIONAL]
+ * 8. nonneg: 'nonnegativity (0 is OFF by default) [OPTIONAL]
+ * 9. print information: 0 (off) or 1 (on) [OPTIONAL]
+ *
+ * Output:
+ * [1] Filtered/regularized image/volume
+ *
+ * This function is based on the Matlab's codes and papers by
+ * [1] Amir Beck and Marc Teboulle, "Fast Gradient-Based Algorithms for Constrained Total Variation Image Denoising and Deblurring Problems"
+ * [2] M. J. Ehrhardt and M. M. Betcke, Multi-Contrast MRI Reconstruction with Structure-Guided Total Variation, SIAM Journal on Imaging Sciences 9(3), pp. 1084–1106
+ */
+
+
+void mexFunction(
+ int nlhs, mxArray *plhs[],
+ int nrhs, const mxArray *prhs[])
+
+{
+ int number_of_dims, iter, dimX, dimY, dimZ, methTV, printswitch, nonneg;
+ const int *dim_array;
+ const int *dim_array2;
+ float *Input, *InputRef, *Output=NULL, lambda, epsil, eta;
+
+ number_of_dims = mxGetNumberOfDimensions(prhs[0]);
+ dim_array = mxGetDimensions(prhs[0]);
+ dim_array2 = mxGetDimensions(prhs[1]);
+
+ /*Handling Matlab input data*/
+ if ((nrhs < 3) || (nrhs > 9)) mexErrMsgTxt("At least 3 parameters is required, all parameters are: Image(2D/3D), Reference(2D/3D), Regularization parameter, iterations number, tolerance, smoothing constant, penalty type ('iso' or 'l1'), nonnegativity switch, print switch");
+
+ Input = (float *) mxGetData(prhs[0]); /*noisy image (2D/3D) */
+ InputRef = (float *) mxGetData(prhs[1]); /* reference image (2D/3D) */
+ lambda = (float) mxGetScalar(prhs[2]); /* regularization parameter */
+ iter = 300; /* default iterations number */
+ epsil = 0.0001; /* default tolerance constant */
+ eta = 0.01; /* default smoothing constant */
+ methTV = 0; /* default isotropic TV penalty */
+ nonneg = 0; /* default nonnegativity switch, off - 0 */
+ printswitch = 0; /*default print is switched, off - 0 */
+
+
+ if (mxGetClassID(prhs[0]) != mxSINGLE_CLASS) {mexErrMsgTxt("The input image must be in a single precision"); }
+ if (mxGetClassID(prhs[1]) != mxSINGLE_CLASS) {mexErrMsgTxt("The input image must be in a single precision"); }
+
+ /*Handling Matlab output data*/
+ dimX = dim_array[0]; dimY = dim_array[1]; dimZ = dim_array[2];
+ if (number_of_dims == 2) { if ((dimX != dim_array2[0]) || (dimY != dim_array2[1])) mexErrMsgTxt("The input images have different dimensionalities");}
+ if (number_of_dims == 3) { if ((dimX != dim_array2[0]) || (dimY != dim_array2[1]) || (dimZ != dim_array2[2])) mexErrMsgTxt("The input volumes have different dimensionalities");}
+
+
+ if ((nrhs == 4) || (nrhs == 5) || (nrhs == 6) || (nrhs == 7) || (nrhs == 8) || (nrhs == 9)) iter = (int) mxGetScalar(prhs[3]); /* iterations number */
+ if ((nrhs == 5) || (nrhs == 6) || (nrhs == 7) || (nrhs == 8) || (nrhs == 9)) epsil = (float) mxGetScalar(prhs[4]); /* tolerance constant */
+ if ((nrhs == 6) || (nrhs == 7) || (nrhs == 8) || (nrhs == 9)) {
+ eta = (float) mxGetScalar(prhs[5]); /* smoothing constant for the gradient of InputRef */
+ }
+ if ((nrhs == 7) || (nrhs == 8) || (nrhs == 9)) {
+ char *penalty_type;
+ penalty_type = mxArrayToString(prhs[6]); /* choosing TV penalty: 'iso' or 'l1', 'iso' is the default */
+ if ((strcmp(penalty_type, "l1") != 0) && (strcmp(penalty_type, "iso") != 0)) mexErrMsgTxt("Choose TV type: 'iso' or 'l1',");
+ if (strcmp(penalty_type, "l1") == 0) methTV = 1; /* enable 'l1' penalty */
+ mxFree(penalty_type);
+ }
+ if ((nrhs == 8) || (nrhs == 9)) {
+ nonneg = (int) mxGetScalar(prhs[7]);
+ if ((nonneg != 0) && (nonneg != 1)) mexErrMsgTxt("Nonnegativity constraint can be enabled by choosing 1 or off - 0");
+ }
+ if (nrhs == 9) {
+ printswitch = (int) mxGetScalar(prhs[8]);
+ if ((printswitch != 0) && (printswitch != 1)) mexErrMsgTxt("Print can be enabled by choosing 1 or off - 0");
+ }
+
+ if (number_of_dims == 2) {
+ dimZ = 1; /*2D case*/
+ Output = (float*)mxGetPr(plhs[0] = mxCreateNumericArray(2, dim_array, mxSINGLE_CLASS, mxREAL));
+ }
+ if (number_of_dims == 3) Output = (float*)mxGetPr(plhs[0] = mxCreateNumericArray(3, dim_array, mxSINGLE_CLASS, mxREAL));
+
+ /* running the function */
+ dTV_FGP_CPU_main(Input, InputRef, Output, lambda, iter, epsil, eta, methTV, nonneg, printswitch, dimX, dimY, dimZ);
+} \ No newline at end of file
diff --git a/Wrappers/Matlab/mex_compile/regularisers_GPU/FGP_dTV_GPU.cpp b/Wrappers/Matlab/mex_compile/regularisers_GPU/FGP_dTV_GPU.cpp
new file mode 100644
index 0000000..5b80616
--- /dev/null
+++ b/Wrappers/Matlab/mex_compile/regularisers_GPU/FGP_dTV_GPU.cpp
@@ -0,0 +1,111 @@
+/*
+ * This work is part of the Core Imaging Library developed by
+ * Visual Analytics and Imaging System Group of the Science Technology
+ * Facilities Council, STFC
+ *
+ * Copyright 2017 Daniil Kazantsev
+ * Copyright 2017 Srikanth Nagella, Edoardo Pasca
+ *
+ * Licensed under the Apache License, Version 2.0 (the "License");
+ * you may not use this file except in compliance with the License.
+ * You may obtain a copy of the License at
+ * http://www.apache.org/licenses/LICENSE-2.0
+ * Unless required by applicable law or agreed to in writing, software
+ * distributed under the License is distributed on an "AS IS" BASIS,
+ * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
+ * See the License for the specific language governing permissions and
+ * limitations under the License.
+ */
+#include "matrix.h"
+#include "mex.h"
+#include "dTV_FGP_GPU_core.h"
+
+/* CUDA implementation of FGP-dTV [1,2] denoising/regularization model (2D/3D case)
+ * which employs structural similarity of the level sets of two images/volumes, see [1,2]
+ * The current implementation updates image 1 while image 2 is being fixed.
+ *
+ * Input Parameters:
+ * 1. Noisy image/volume [REQUIRED]
+ * 2. Additional reference image/volume of the same dimensions as (1) [REQUIRED]
+ * 3. lambdaPar - regularization parameter [REQUIRED]
+ * 4. Number of iterations [OPTIONAL]
+ * 5. eplsilon: tolerance constant [OPTIONAL]
+ * 6. eta: smoothing constant to calculate gradient of the reference [OPTIONAL] *
+ * 7. TV-type: methodTV - 'iso' (0) or 'l1' (1) [OPTIONAL]
+ * 8. nonneg: 'nonnegativity (0 is OFF by default) [OPTIONAL]
+ * 9. print information: 0 (off) or 1 (on) [OPTIONAL]
+ *
+ * Output:
+ * [1] Filtered/regularized image/volume
+ *
+ * This function is based on the Matlab's codes and papers by
+ * [1] Amir Beck and Marc Teboulle, "Fast Gradient-Based Algorithms for Constrained Total Variation Image Denoising and Deblurring Problems"
+ * [2] M. J. Ehrhardt and M. M. Betcke, Multi-Contrast MRI Reconstruction with Structure-Guided Total Variation, SIAM Journal on Imaging Sciences 9(3), pp. 1084–1106
+ */
+void mexFunction(
+ int nlhs, mxArray *plhs[],
+ int nrhs, const mxArray *prhs[])
+
+{
+ int number_of_dims, iter, dimX, dimY, dimZ, methTV, printswitch, nonneg;
+ const int *dim_array;
+ const int *dim_array2;
+ float *Input, *InputRef, *Output=NULL, lambda, epsil, eta;
+
+ number_of_dims = mxGetNumberOfDimensions(prhs[0]);
+ dim_array = mxGetDimensions(prhs[0]);
+ dim_array2 = mxGetDimensions(prhs[1]);
+
+ /*Handling Matlab input data*/
+ if ((nrhs < 3) || (nrhs > 9)) mexErrMsgTxt("At least 3 parameters is required, all parameters are: Image(2D/3D), Reference(2D/3D), Regularization parameter, iterations number, tolerance, smoothing constant, penalty type ('iso' or 'l1'), nonnegativity switch, print switch");
+
+ Input = (float *) mxGetData(prhs[0]); /*noisy image (2D/3D) */
+ InputRef = (float *) mxGetData(prhs[1]); /* reference image (2D/3D) */
+ lambda = (float) mxGetScalar(prhs[2]); /* regularization parameter */
+ iter = 300; /* default iterations number */
+ epsil = 0.0001; /* default tolerance constant */
+ eta = 0.01; /* default smoothing constant */
+ methTV = 0; /* default isotropic TV penalty */
+ nonneg = 0; /* default nonnegativity switch, off - 0 */
+ printswitch = 0; /*default print is switched, off - 0 */
+
+
+ if (mxGetClassID(prhs[0]) != mxSINGLE_CLASS) {mexErrMsgTxt("The input image must be in a single precision"); }
+ if (mxGetClassID(prhs[1]) != mxSINGLE_CLASS) {mexErrMsgTxt("The input image must be in a single precision"); }
+
+ /*Handling Matlab output data*/
+ dimX = dim_array[0]; dimY = dim_array[1]; dimZ = dim_array[2];
+ if (number_of_dims == 2) { if ((dimX != dim_array2[0]) || (dimY != dim_array2[1])) mexErrMsgTxt("The input images have different dimensionalities");}
+ if (number_of_dims == 3) { if ((dimX != dim_array2[0]) || (dimY != dim_array2[1]) || (dimZ != dim_array2[2])) mexErrMsgTxt("The input volumes have different dimensionalities");}
+
+
+ if ((nrhs == 4) || (nrhs == 5) || (nrhs == 6) || (nrhs == 7) || (nrhs == 8) || (nrhs == 9)) iter = (int) mxGetScalar(prhs[3]); /* iterations number */
+ if ((nrhs == 5) || (nrhs == 6) || (nrhs == 7) || (nrhs == 8) || (nrhs == 9)) epsil = (float) mxGetScalar(prhs[4]); /* tolerance constant */
+ if ((nrhs == 6) || (nrhs == 7) || (nrhs == 8) || (nrhs == 9)) {
+ eta = (float) mxGetScalar(prhs[5]); /* smoothing constant for the gradient of InputRef */
+ }
+ if ((nrhs == 7) || (nrhs == 8) || (nrhs == 9)) {
+ char *penalty_type;
+ penalty_type = mxArrayToString(prhs[6]); /* choosing TV penalty: 'iso' or 'l1', 'iso' is the default */
+ if ((strcmp(penalty_type, "l1") != 0) && (strcmp(penalty_type, "iso") != 0)) mexErrMsgTxt("Choose TV type: 'iso' or 'l1',");
+ if (strcmp(penalty_type, "l1") == 0) methTV = 1; /* enable 'l1' penalty */
+ mxFree(penalty_type);
+ }
+ if ((nrhs == 8) || (nrhs == 9)) {
+ nonneg = (int) mxGetScalar(prhs[7]);
+ if ((nonneg != 0) && (nonneg != 1)) mexErrMsgTxt("Nonnegativity constraint can be enabled by choosing 1 or off - 0");
+ }
+ if (nrhs == 9) {
+ printswitch = (int) mxGetScalar(prhs[8]);
+ if ((printswitch != 0) && (printswitch != 1)) mexErrMsgTxt("Print can be enabled by choosing 1 or off - 0");
+ }
+
+ if (number_of_dims == 2) {
+ dimZ = 1; /*2D case*/
+ Output = (float*)mxGetPr(plhs[0] = mxCreateNumericArray(2, dim_array, mxSINGLE_CLASS, mxREAL));
+ }
+ if (number_of_dims == 3) Output = (float*)mxGetPr(plhs[0] = mxCreateNumericArray(3, dim_array, mxSINGLE_CLASS, mxREAL));
+
+ /* running the function */
+ dTV_FGP_GPU_main(Input, InputRef, Output, lambda, iter, epsil, eta, methTV, nonneg, printswitch, dimX, dimY, dimZ);
+} \ No newline at end of file