# HG changeset patch # User bmailhe # Date 1346752836 -3600 # Node ID 5c8bcdadb380bbf09004f07face4b516b3f485dc # Parent 1fbd28dfb99efbcf84159e7b9c634309c55ac1bd added UNLocBox interface, example and a 128x128 Lena diff -r 1fbd28dfb99e -r 5c8bcdadb380 config/SMALL_solve_config.m --- a/config/SMALL_solve_config.m Wed Aug 29 10:39:14 2012 +0100 +++ b/config/SMALL_solve_config.m Tue Sep 04 11:00:36 2012 +0100 @@ -74,6 +74,14 @@ end [y, cost] = wrapper_mm_solver(b, A, solver.param); + +elseif (strcmpi(solver.toolbox, 'UNLocBox')) + if ~isa(Problem.A,'float') + % MMbox does not accept implicit dictionary definition + A = opToMatrix(Problem.A, 1); + end + + y = unloc_solver(b, A, solver.param,solver.name); %% % Please do not make any changes to the 'SMALL_solve_config.m' file diff -r 1fbd28dfb99e -r 5c8bcdadb380 data/images/lena_little.png Binary file data/images/lena_little.png has changed diff -r 1fbd28dfb99e -r 5c8bcdadb380 examples/Image Denoising/SMALL_ImgDenoise_UNLocBox.m --- /dev/null Thu Jan 01 00:00:00 1970 +0000 +++ b/examples/Image Denoising/SMALL_ImgDenoise_UNLocBox.m Tue Sep 04 11:00:36 2012 +0100 @@ -0,0 +1,85 @@ +%% Demostration of the use of the UNLocBox solver to solve an image denoising problem of the SMALLBox form +% +% argmin_x ||x||_1 such that ||Ax-b||_2 < sima_t +% +% The dictionary A is first learn using SMALLbox +% +% Queen Mary University +% 29 August 2012 +% Nathanael Perraudin +% nathanael.perraudin@epfl.ch + + +% % Initialisation +% clear all; +% close all; +% clc; + + +%% Load an image and create the problem + +% Defining Image Denoising Problem as Dictionary Learning +% Problem. As an input we set the number of training patches. + +SMALL.Problem = generateImageDenoiseProblem('', 40000); + + +%% +% Use KSVD Dictionary Learning Algorithm to Learn overcomplete dictionary + +% Initialising Dictionary structure +% Setting Dictionary structure fields (toolbox, name, param, D and time) +% to zero values + +SMALL.DL=SMALL_init_DL(); + +% Defining the parameters needed for dictionary learning + +SMALL.DL.toolbox = 'KSVD'; +SMALL.DL.name = 'ksvd'; + +% Defining the parameters for KSVD +% In this example we are learning 256 atoms in 20 iterations, so that +% every patch in the training set can be represented with target error in +% L2-norm (EData) +% Type help ksvd in MATLAB prompt for more options. + +Edata=sqrt(prod(SMALL.Problem.blocksize)) * SMALL.Problem.sigma * SMALL.Problem.gain; +maxatoms = floor(prod(SMALL.Problem.blocksize)/2); + +SMALL.DL.param=struct(... + 'Edata', Edata,... + 'initdict', SMALL.Problem.initdict,... + 'dictsize', SMALL.Problem.p,... + 'iternum', 20,... + 'memusage', 'high'); + +% Learn the dictionary + +SMALL.DL = SMALL_learn(SMALL.Problem, SMALL.DL(1)); + +% Set SMALL.Problem.A dictionary +% (backward compatiblity with SPARCO: solver structure communicate +% only with Problem structure, ie no direct communication between DL and +% solver structures) + +SMALL.Problem.A = SMALL.DL.D; +SMALL.Problem.reconstruct = @(x) ImageDenoise_reconstruct(x, SMALL.Problem); + +%% Solving the problem with UNLocBox +% This might not be the better way to solve the problem... + + +% Set the different parameter + +SMALL.solver=SMALL_init_solver; % Initialisation +SMALL.solver.toolbox='UNLocBox'; % select the UNLocBox solver +SMALL.solver.name='Douglas_Rachford'; % 'Forward_Backard' 'ADMM' Warning forward backward still need some review... +SMALL.solver.param.sigma=1.15*sqrt(SMALL.Problem.m*SMALL.Problem.n)*SMALL.Problem.sigma; % set the radius of the ball +SMALL.solver.param.max_iter=100; + + +SMALL.solver=SMALL_solve(SMALL.Problem, SMALL.solver); + + +SMALL_ImgDeNoiseResult(SMALL); \ No newline at end of file diff -r 1fbd28dfb99e -r 5c8bcdadb380 util/small_to_unloc.m --- /dev/null Thu Jan 01 00:00:00 1970 +0000 +++ b/util/small_to_unloc.m Tue Sep 04 11:00:36 2012 +0100 @@ -0,0 +1,43 @@ +function [b2,f1,f2,param] = small_to_unloc(b,A,param) + + + if nargin<3, param=struct; end + % set parameters + if ~isfield(param, 'verbose'), param.verbose = 0; end + if ~isfield(param, 'T'), param.T = 1; end + if ~isfield(param, 'gamma'), param.gamma = 1/(1.1*norm(A)^2); end + if ~isfield(param, 'sigma'), param.sigma = 1; end + + + %set function f2 + param_f2.verbose=param.verbose; + param_f2.A=@(x) A*x; + param_f2.At=@(x) A'*x; + param_f2.y=b; + param_f2.tight=0; + param_f2.nu=norm(A,2)^2; + param_f2.epsilon=param.sigma; + + f2.prox= @(x,l) fast_proj_B2(x,l,param_f2); % douglas rachford, admm + f2.grad= @(x) 2*A'*(A*x-b); %forward backward + f2.x0=A'*b; + f2.norm=@(x) 0; + + + + + %set function f1 + param_f1.verbose=param.verbose; + + f1.x0=A'*b; + f1.norm= @(x) sum(sum(abs(x))); + f1.prox= @(x,l) prox_L1(x,param.T*l,param_f1); + + + + %set initial point + b2=A'*b; + + +end + diff -r 1fbd28dfb99e -r 5c8bcdadb380 util/unloc_solver.m --- /dev/null Thu Jan 01 00:00:00 1970 +0000 +++ b/util/unloc_solver.m Tue Sep 04 11:00:36 2012 +0100 @@ -0,0 +1,31 @@ +function [ sol ] = unloc_solver( b, A, param,name ) +%UNLOC_SOLVER Solve the minimisation problem using UNLOCBOX toolbox +% +% argmin_x ||x||_1 such that ||Ax-b||_2 < sima_t + +if nargin<3, param=struct; end + +% set parameters +if ~isfield(param, 'verbose'), param.verbose = 1; end +if ~isfield(param, 'T'), param.T = 256; end +if ~isfield(param, 'sigma'), param.sigma = 1; end +if ~isfield(param, 'max_iter'), param.max_iter = 100; end +if ~isfield(param, 'epsilon'), param.epsilon = 1e-3; end + +[b2,f1,f2,param] = small_to_unloc(b,A,param); + +if strcmpi(name, 'Douglas_Rachford') + sol=douglas_rachford(b2,f2,f1,param); +elseif strcmpi(name, 'Forward_Backard') + sol=forward_backward(b2,f1,f2,param); +elseif strcmpi(name, 'ADMM') + opL= @(x) x; + sol=admm(b2,f1,f2,opL,param); +else + error('UNLocBox tells you: Unknown solver name!') +end; + + + +end +