|
| 1 | +function newvol = laplacefill(vol, seeds, solver, tol, maxiter, varargin) |
| 2 | +% |
| 3 | +% newvol = laplacefill(vol, seeds, solver, tol, maxiter, ...) |
| 4 | +% |
| 5 | +% perform a hole-fill or flood-fill in a 2-D or 3-D binary image using |
| 6 | +% Laplace solver |
| 7 | +% |
| 8 | +% author: Qianqian Fang, <q.fang at neu.edu> |
| 9 | +% |
| 10 | +% input: |
| 11 | +% vol: a 2-D or 3-D binary image |
| 12 | +% seeds: seeds for flood-fill, 2 or 3-column matrix; if empty, perfrom |
| 13 | +% hole-fill |
| 14 | +% solver: (optional) linear solver, if not given, use bicgstab |
| 15 | +% tol: (optional) linear solver convergence tolerance, default 1e-5 |
| 16 | +% maxiter: (optional) linear solver max iteration, default 3000 |
| 17 | +% |
| 18 | +% additional parameters that are supported by any iterative solver, |
| 19 | +% such as qmr, pcg, gmres, tfqmr, mldivide, minres etc |
| 20 | +% |
| 21 | +% output: |
| 22 | +% newvol: the volume image after flood-fill or hole-fill |
| 23 | +% |
| 24 | +% example: |
| 25 | +% a = zeros(60, 80, 90); |
| 26 | +% a(20:40, 40:70, 30:80)=1; |
| 27 | +% a(25:35, 50:60, 50:60)=0; |
| 28 | +% newvol = laplacefill(a); |
| 29 | +% |
| 30 | +% -- this function is part of iso2mesh toolbox (http://iso2mesh.sf.net) |
| 31 | +% |
| 32 | + |
| 33 | +dims = size(vol) + 2; |
| 34 | + |
| 35 | +interiornode = prod(dims - 2); |
| 36 | + |
| 37 | +zeroidx = find(vol == 0); |
| 38 | + |
| 39 | +if (ndims(vol) == 3) |
| 40 | + [ix, iy, iz] = ndgrid(2:dims(1) - 1, 2:dims(2) - 1, 2:dims(3) - 1); |
| 41 | +else |
| 42 | + [ix, iy, iz] = ndgrid(2:dims(1) - 1, 2:dims(2) - 1, 2); |
| 43 | + dims(3) = 3; |
| 44 | +end |
| 45 | + |
| 46 | +ix = ix(zeroidx); |
| 47 | +iy = iy(zeroidx); |
| 48 | +iz = iz(zeroidx); |
| 49 | + |
| 50 | +idx = sub2ind(dims - 2, ix - 1, iy - 1, iz - 1); % index of all interior nodes in the subset of interior nodes |
| 51 | +idx = idx(:)'; |
| 52 | +nnzpos = {idx, idx(ix > 2), idx(ix < dims(1) - 1), idx(iy > 2), idx(iy < dims(2) - 1), idx(iz > 2), idx(iz < dims(3) - 1)}; |
| 53 | + |
| 54 | +nonzeroidx = find(vol > 0)'; |
| 55 | + |
| 56 | +seedidx = []; |
| 57 | + |
| 58 | +if (nargin > 1 && size(seeds, 2) == 3) |
| 59 | + seedidx = sub2ind(dims(1:3) - 2, seeds(:, 1), seeds(:, 2), seeds(:, 3)); |
| 60 | +end |
| 61 | + |
| 62 | +Amat = sparse([seedidx, nonzeroidx, nnzpos{1}, nnzpos{2} - 1, nnzpos{3} + 1, ... |
| 63 | + nnzpos{4} - dims(1) + 2, nnzpos{5} + dims(1) - 2, ... |
| 64 | + nnzpos{6} - (dims(1) - 2) * (dims(2) - 2), nnzpos{7} + (dims(1) - 2) * (dims(2) - 2)], ... |
| 65 | + [seedidx, nonzeroidx, nnzpos{:}], ... |
| 66 | + [ones(1, length(seedidx)) ones(1, length(nonzeroidx)) -ones(1, length(idx)) ... |
| 67 | + (1 / 6) * ones(1, sum(cellfun(@(x) length(x), nnzpos(2:end))))], ... |
| 68 | + interiornode, interiornode); |
| 69 | + |
| 70 | +if (~isempty(seedidx)) |
| 71 | + b = sparse(seedidx, ones(size(seedidx)), ones(size(seedidx)), interiornode, 1); |
| 72 | +else |
| 73 | + if (ndims(vol) == 3) |
| 74 | + boundnode = idx(ix == 2 | iy == 2 | iz == 2 | ix == dims(1) - 1 | iy == dims(2) - 1 | iz == dims(3) - 1); |
| 75 | + else |
| 76 | + boundnode = idx(ix == 2 | iy == 2 | ix == dims(1) - 1 | iy == dims(2) - 1); |
| 77 | + end |
| 78 | + b = sparse(boundnode, ones(size(boundnode)), -(1 / 6) * ones(size(boundnode)), interiornode, 1); |
| 79 | +end |
| 80 | + |
| 81 | +if (nargin < 4) |
| 82 | + tol = 1e-10; |
| 83 | +end |
| 84 | + |
| 85 | +if (nargin < 5) |
| 86 | + maxiter = 3000; |
| 87 | +end |
| 88 | + |
| 89 | +if (nargin > 2) |
| 90 | + mysolver = str2fun(solver); |
| 91 | + [newvol, flag] = mysolver(Amat, b, tol, maxiter, varargin{:}); |
| 92 | +else |
| 93 | + [newvol, flag] = bicgstab(Amat, b, tol, maxiter, varargin{:}); |
| 94 | +end |
| 95 | + |
| 96 | +if (flag) |
| 97 | + newvol = gmres(Amat, b, 100, tol, maxiter, varargin{:}); |
| 98 | +end |
| 99 | + |
| 100 | +newvol = reshape(full(newvol), size(vol)); |
| 101 | + |
| 102 | +if (isempty(seedidx)) |
| 103 | + newvol = ~newvol; |
| 104 | +end |
0 commit comments