Skip to content

Commit 11c106b

Browse files
committed
[feat] drop image proc toolbox dep. self-implement volgrow,shrink,close,open
1 parent 2f2d2c6 commit 11c106b

8 files changed

+179
-55
lines changed

deislands2d.m

+2-2
Original file line numberDiff line numberDiff line change
@@ -26,11 +26,11 @@
2626

2727
cleanimg = zeros(size(img));
2828
if (sum(img(:)))
29-
img = imclose(img, strel('disk', 3));
29+
img = volclose(img, 1, true(5, 5));
3030
islands = bwislands(img);
3131
end
3232

33-
if (length(islands))
33+
if (~isempty(islands))
3434
% remove small islands of the foreground
3535
maxblock = -1;
3636
maxblockid = -1;

fillholes3d.m

+33-11
Original file line numberDiff line numberDiff line change
@@ -1,4 +1,4 @@
1-
function resimg = fillholes3d(img, maxgap)
1+
function resimg = fillholes3d(img, maxgap, varargin)
22
%
33
% resimg=fillholes3d(img,maxgap)
44
%
@@ -18,19 +18,41 @@
1818
% -- this function is part of iso2mesh toolbox (http://iso2mesh.sf.net)
1919
%
2020

21-
if (maxgap)
22-
resimg = imclose(img, strel(ones(maxgap, maxgap, maxgap)));
21+
if (nargin > 1 && maxgap)
22+
resimg = volclose(img, maxgap);
2323
else
2424
resimg = img;
2525
end
2626

27-
if (isoctavemesh)
28-
if (~exist('bwfill'))
29-
error('you need to install octave-image toolbox first');
30-
end
31-
for i = 1:size(resimg, 3)
32-
resimg(:, :, i) = bwfill(resimg(:, :, i), 'holes');
33-
end
27+
% if (exist('imfill', 'file'))
28+
% resimg = imfill(resimg, 'holes');
29+
% return;
30+
% end
31+
32+
newimg = ones(size(resimg) + 2);
33+
34+
oldimg = zeros(size(newimg));
35+
if (ndims(resimg) == 3)
36+
oldimg(2:end - 1, 2:end - 1, 2:end - 1) = resimg;
37+
newimg(2:end - 1, 2:end - 1, 2:end - 1) = 0;
3438
else
35-
resimg = imfill(resimg, 'holes');
39+
oldimg(2:end - 1, 2:end - 1) = resimg;
40+
newimg(2:end - 1, 2:end - 1) = 0;
41+
end
42+
43+
newsum = sum(newimg(:));
44+
oldsum = -1;
45+
46+
while (newsum ~= oldsum)
47+
newimg = (volgrow(newimg, 1, varargin{:}) & ~(oldimg > 0));
48+
oldsum = newsum;
49+
newsum = sum(newimg(:));
3650
end
51+
52+
if (ndims(resimg) == 3)
53+
resimg = newimg(2:end - 1, 2:end - 1, 2:end - 1);
54+
else
55+
resimg = newimg(2:end - 1, 2:end - 1);
56+
end
57+
58+
resimg = double(~resimg);

thickenbinvol.m

+2-14
Original file line numberDiff line numberDiff line change
@@ -1,4 +1,4 @@
1-
function vol = thickenbinvol(vol, layer)
1+
function vol = thickenbinvol(varargin)
22
%
33
% vol=thickenbinvol(vol,layer)
44
%
@@ -19,16 +19,4 @@
1919
% -- this function is part of iso2mesh toolbox (http://iso2mesh.sf.net)
2020
%
2121

22-
dim = size(vol);
23-
dxy = dim(1) * dim(2);
24-
fulllen = prod(dim);
25-
26-
offs = [1, -1, dim(1), -dim(1), dxy, -dxy];
27-
for i = 1:layer
28-
idx = find(vol);
29-
for j = 1:6
30-
idxnew = idx + offs(j);
31-
idxnew = idxnew(find(idxnew > 0 & idxnew < fulllen));
32-
vol(idxnew) = 1;
33-
end
34-
end
22+
vol = volgrow(varargin{:});

thinbinvol.m

+3-28
Original file line numberDiff line numberDiff line change
@@ -1,6 +1,6 @@
1-
function vol = thinbinvol(vol, layer, nobd)
1+
function vol = thinbinvol(varargin)
22
%
3-
% vol=thinbinvol(vol,layer,nobd)
3+
% vol=thinbinvol(vol,layer,mask)
44
%
55
% thinning a binary volume by a given pixel width
66
% this is similar to bwmorph(vol,'thin',n) except
@@ -21,29 +21,4 @@
2121
% -- this function is part of iso2mesh toolbox (http://iso2mesh.sf.net)
2222
%
2323

24-
dim = size(vol);
25-
dxy = dim(1) * dim(2);
26-
fulllen = prod(dim);
27-
28-
if (nargin < 3)
29-
nobd = 0;
30-
end
31-
32-
if (nobd == 1)
33-
bdmask = vol;
34-
if (ndims(vol) == 2)
35-
bdmask(2:end - 1, 2:end - 1) = 0;
36-
elseif (ndims(vol) == 3)
37-
bdmask(2:end - 1, 2:end - 1, 2:end - 1) = 0;
38-
end
39-
end
40-
41-
for i = 1:layer
42-
idx = find(~vol);
43-
idxnew = [idx + 1; idx - 1; idx + dim(1); idx - dim(1); idx + dxy; idx - dxy];
44-
idxnew = idxnew(find(idxnew > 0 & idxnew < fulllen));
45-
vol(idxnew) = 0;
46-
if (nobd)
47-
vol = vol | bdmask;
48-
end
49-
end
24+
vol = volshrink(varargin{:});

volclose.m

+26
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,26 @@
1+
function newvol = volclose(varargin)
2+
%
3+
% vol = volclose(vol, layer, mask)
4+
%
5+
% close a binary volume by applying layer-counter of growth operation
6+
% followed by layer-count of shrinkage operation
7+
%
8+
% author: Qianqian Fang, <q.fang at neu.edu>
9+
%
10+
% input:
11+
% vol: a volumetric binary image
12+
% layer: number of iterations for the thickenining
13+
% mask: (optional) a 3D neighborhood mask
14+
%
15+
% output:
16+
% vol: the volume image after closing
17+
%
18+
% -- this function is part of iso2mesh toolbox (http://iso2mesh.sf.net)
19+
%
20+
21+
if (isempty(varargin))
22+
error('must provide a volume');
23+
end
24+
25+
newvol = volgrow(varargin{:});
26+
newvol = volshrink(newvol, varargin{2:end});

volgrow.m

+43
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,43 @@
1+
function newvol = volgrow(vol, layer, mask)
2+
%
3+
% vol = volgrow(vol, layer, mask)
4+
%
5+
% thickening a binary image or volume by a given pixel width
6+
% this is similar to bwmorph(vol,'thicken',3) except
7+
% this does it in both 2d and 3d
8+
%
9+
% author: Qianqian Fang, <q.fang at neu.edu>
10+
%
11+
% input:
12+
% vol: a volumetric binary image
13+
% layer: number of iterations for the thickenining
14+
% mask: (optional) a 3D neighborhood mask
15+
%
16+
% output:
17+
% vol: the volume image after the thickening
18+
%
19+
% -- this function is part of iso2mesh toolbox (http://iso2mesh.sf.net)
20+
%
21+
22+
if (nargin < 2)
23+
layer = 1;
24+
end
25+
26+
if (nargin < 3)
27+
if (ndims(vol) == 3)
28+
mask = zeros(3, 3, 3);
29+
mask(2, 2, :) = 1;
30+
mask(:, 2, 2) = 1;
31+
mask(2, :, 2) = 1;
32+
else
33+
mask = [0 1 0; 1 1 1; 0 1 0];
34+
end
35+
end
36+
37+
newvol = vol;
38+
39+
for i = 1:layer
40+
newvol = (convn(newvol, mask, 'same') > 0);
41+
end
42+
43+
newvol = double(newvol);

volopen.m

+26
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,26 @@
1+
function newvol = volopen(varargin)
2+
%
3+
% vol = volopen(vol, layer, mask)
4+
%
5+
% open a binary volume by applying layer-counter of shrink operation
6+
% followed by layer-count of growth operation
7+
%
8+
% author: Qianqian Fang, <q.fang at neu.edu>
9+
%
10+
% input:
11+
% vol: a volumetric binary image
12+
% layer: number of iterations for the thickenining
13+
% mask: (optional) a 3D neighborhood mask
14+
%
15+
% output:
16+
% vol: the volume image after opening
17+
%
18+
% -- this function is part of iso2mesh toolbox (http://iso2mesh.sf.net)
19+
%
20+
21+
if (isempty(varargin))
22+
error('must provide a volume');
23+
end
24+
25+
newvol = volshrink(varargin{:});
26+
newvol = volgrow(newvol, varargin{2:end});

volshrink.m

+44
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,44 @@
1+
function newvol = volshrink(vol, layer, mask)
2+
%
3+
% vol = volshrink(vol, layer, mask)
4+
%
5+
% thinning a binary image or volume by a given pixel width
6+
% this is similar to bwmorph(vol,'thin',3) except
7+
% this does it in both 2d and 3d
8+
%
9+
% author: Qianqian Fang, <q.fang at neu.edu>
10+
%
11+
% input:
12+
% vol: a volumetric binary image
13+
% layer: number of iterations for the thickenining
14+
% mask: (optional) a 3D neighborhood mask
15+
%
16+
% output:
17+
% vol: the volume image after the thinning operations
18+
%
19+
% -- this function is part of iso2mesh toolbox (http://iso2mesh.sf.net)
20+
%
21+
22+
if (nargin < 2)
23+
layer = 1;
24+
end
25+
26+
if (nargin < 3)
27+
if (ndims(vol) == 3)
28+
mask = zeros(3, 3, 3);
29+
mask(2, 2, :) = 1;
30+
mask(:, 2, 2) = 1;
31+
mask(2, :, 2) = 1;
32+
else
33+
mask = [0 1 0; 1 1 1; 0 1 0];
34+
end
35+
end
36+
37+
totalmask = sum(mask(:));
38+
newvol = vol;
39+
40+
for i = 1:layer
41+
newvol = (convn(logical(newvol), logical(mask), 'same') == totalmask);
42+
end
43+
44+
newvol = double(newvol);

0 commit comments

Comments
 (0)