-
Notifications
You must be signed in to change notification settings - Fork 48
/
smooth2.m
645 lines (500 loc) · 16.9 KB
/
smooth2.m
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
function [vert,conn,tria,tnum] = smooth2(varargin)
%SMOOTH2 "hill-climbing" mesh-smoothing for two-dimensional,
%2-simplex triangulations.
% [VERT,EDGE,TRIA,TNUM] = SMOOTH2(VERT,EDGE,TRIA,TNUM) re-
% turns a "smoothed" triangulation {VERT,TRIA}, incorpora-
% ting "optimised" vertex coordinates and mesh topology.
%
% VERT is a V-by-2 array of XY coordinates in the triangu-
% lation, EDGE is an array of constrained edges, TRIA is a
% T-by-3 array of triangles, and TNUM is a T-by-1 array of
% part indices. Each row of TRIA and EDGE define an eleme-
% nt. VERT(TRIA(II,1),:), VERT(TRIA(II,2),:) and VERT(TRIA
% (II,3),:) are the coordinates of the II-TH triangle. The
% edges in EDGE are defined in a similar manner. NUM is an
% array of part indexing, such that TNUM(II) is the index
% of the part in which the II-TH triangle resides.
%
% [VERT,EDGE,TRIA,TNUM] = SMOOTH2(... ,OPTS) passes an ad-
% ditional options structure OPTS, containing user-defined
% parameters, including:
%
% - OPTS.VTOL = {+1.0E-02} -- relative vertex movement tole-
% rance, smoothing is converged when (VNEW-VERT) <= VTOL *
% VLEN, where VLEN is a local effective length-scale.
%
% - OPTS.ITER = {+32} -- max. number of smoothing iterations
%
% - OPTS.DISP = {+ 4} -- smoothing verbosity. Set to INF for
% quiet execution.
%
% See also REFINE2, TRICOST, TRIDEMO
% This routine is loosely based on the DISTMESH algorithm,
% employing a "spring-based" analogy to redistribute mesh
% vertices. Such an approach is described in: P.O. Persson
% and Gilbert Strang. "A simple mesh generator in MATLAB."
% SIAM review 46(2) 2004, pp: 329--345. Details of the al-
% gorithm used here are somewhat different, with an alter-
% ative spring-based update employed, in addition to hill-
% climbing element quality guarantees, and vertex density
% controls.
%-----------------------------------------------------------
% Darren Engwirda : 2017 --
% Email : [email protected]
% Last updated : 13/02/2020
%-----------------------------------------------------------
vert = []; conn = []; tria = [] ;
tnum = [];
opts = []; hfun = []; harg = {} ;
%---------------------------------------------- extract args
if (nargin>=+1), vert = varargin{1}; end
if (nargin>=+2), conn = varargin{2}; end
if (nargin>=+3), tria = varargin{3}; end
if (nargin>=+4), tnum = varargin{4}; end
if (nargin>=+5), opts = varargin{5}; end
[opts] = makeopt(opts) ;
%---------------------------------------------- default CONN
if (isempty(conn))
[edge] = tricon2(tria);
ebnd = edge(:,4) < +1; %-- use bnd edge
conn = edge(ebnd,1:2);
end
%---------------------------------------------- default TNUM
if (isempty(tnum))
tnum = ones(size(tria, 1), 1) ;
end
%---------------------------------------------- basic checks
if ( ~isnumeric(vert) || ...
~isnumeric(conn) || ...
~isnumeric(tria) || ...
~isnumeric(tnum) || ...
~isstruct (opts) )
error('smooth2:incorrectInputClass' , ...
'Incorrect input class.') ;
end
%---------------------------------------------- basic checks
if (ndims(vert) ~= +2 || ...
ndims(conn) ~= +2 || ...
ndims(tria) ~= +2 || ...
ndims(tnum) ~= +2 )
error('smooth2:incorrectDimensions' , ...
'Incorrect input dimensions.');
end
if (size(vert,2)~= +2 || ...
size(conn,2)~= +2 || ...
size(tria,2)~= +3 || ...
size(tnum,2)~= +1 || ...
size(tria,1)~= size(tnum,1) )
error('smooth2:incorrectDimensions' , ...
'Incorrect input dimensions.');
end
nvrt = size(vert,1) ;
%---------------------------------------------- basic checks
if (min(min(conn(:,1:2))) < +1 || ...
max(max(conn(:,1:2))) > nvrt )
error('smooth2:invalidInputs', ...
'Invalid EDGE input array.') ;
end
if (min(min(tria(:,1:3))) < +1 || ...
max(max(tria(:,1:3))) > nvrt )
error('smooth2:invalidInputs', ...
'Invalid TRIA input array.') ;
end
%---------------------------------------------- output title
if (~isinf(opts.disp))
fprintf(1,'\n') ;
fprintf(1,' Smooth triangulation...\n') ;
fprintf(1,'\n') ;
fprintf(1,[...
' -------------------------------------------------------\n', ...
' |ITER.| |MOVE(X)| |DTRI(X)| \n', ...
' -------------------------------------------------------\n', ...
] ) ;
end
%---------------------------------------------- polygon bnds
node = vert; PSLG = conn; part = {};
pmax = max(tnum(:));
for ppos = +1 : pmax
tsel = tnum == ppos ;
tcur = tria(tsel,:) ;
[ecur,tcur] ...
= tricon2 (tcur) ;
ebnd = ecur(:,4)==0 ;
same = setset2( ...
PSLG,ecur(ebnd,1:2));
part{ppos} = find(same) ;
end
%---------------------------------------------- inflate bbox
vmin = min(vert,[],1);
vmax = max(vert,[],1);
vdel = vmax - 1.*vmin;
vmin = vmin - .5*vdel;
vmax = vmax + .5*vdel;
vbox = [
vmin(1), vmin(2)
vmax(1), vmin(2)
vmax(1), vmax(2)
vmin(1), vmax(2)
] ;
vert = [vert ; vbox] ;
%---------------------------------------------- DO MESH ITER
tnow = tic ;
tcpu = struct('full',0.,'dtri',0., ...
'tcon',0.,'iter',0.,'undo',0., ...
'keep',0.) ;
for iter = +1 : opts.iter
%------------------------------------------ inflate adj.
ttic = tic ;
[edge,tria] = tricon2(tria,conn) ;
tcpu.tcon = ...
tcpu.tcon + toc(ttic) ;
%------------------------------------------ compute scr.
oscr = triscr2(vert,tria) ;
%------------------------------------------ vert. iter's
ttic = tic ;
nvrt = size(vert,1);
nedg = size(edge,1);
IMAT = sparse( ...
edge(:,1),(1:nedg)',+1,nvrt,nedg) ;
JMAT = sparse( ...
edge(:,2),(1:nedg)',+1,nvrt,nedg) ;
EMAT = IMAT + JMAT ;
vdeg = sum(EMAT,2) ; %-- vertex |deg|
free = (vdeg == 0) ;
vold = vert ;
for isub = +1 : max(+2,min(+8,iter))
%-- compute HFUN at vert/midpoints
hvrt = evalhfn(vert, ...
edge,EMAT,hfun,harg) ;
hmid = hvrt(edge(:,1),:) ...
+ hvrt(edge(:,2),:) ;
hmid = hmid * +.5 ;
%-- calc. relative edge extensions
evec = vert(edge(:,2),:) ...
- vert(edge(:,1),:) ;
elen = ...
sqrt(sum(evec.^2,2)) ;
scal = +1.0-elen./hmid ;
scal = min (+1.0, scal);
scal = max (-1.0, scal);
%-- projected points from each end
ipos = vert(edge(:,1),:) ...
-.67*[scal,scal].*evec;
jpos = vert(edge(:,2),:) ...
+.67*[scal,scal].*evec;
%scal = ... %-- nlin. weight
% max(abs(scal).^.5,eps^.75);
scal = ...
max(abs(scal).^ 1,eps^.75);
%-- sum contributions edge-to-vert
vnew = ...
IMAT*([scal,scal] .* ipos) ...
+ JMAT*([scal,scal] .* jpos) ;
vsum = max(EMAT*scal,eps^.75);
vnew = vnew ./ [vsum,vsum] ;
%-- fixed points. edge projection?
vnew(conn(:),1:2) = ...
vert(conn(:),1:2) ;
vnew(vdeg==0,1:2) = ...
vert(vdeg==0,1:2) ;
%-- reset for the next local iter.
vert = vnew ;
end
tcpu.iter = ...
tcpu.iter + toc(ttic) ;
%------------------------------------------ hill-climber
ttic = tic ;
%-- unwind vert. upadte if score lower
nscr = ones(size(tria,1),1);
btri = true(size(tria,1),1);
umax = + 8 ;
for undo = +1 : umax
nscr(btri) = triscr2( ...
vert,tria(btri,:)) ;
%-- TRUE if tria needs "unwinding"
smin = +.70 ;
smax = +.90 ;
sdel = .025 ;
stol = smin+iter*sdel;
stol = min (smax,stol) ;
btri = nscr <= stol ...
& nscr < oscr ;
if (~any(btri)), break; end
%-- relax toward old vert. coord's
ivrt = ...
unique(tria(btri,1:3));
bvrt = ...
false(size(vert,1),1) ;
bvrt(ivrt) = true;
if (undo ~= umax)
bnew = +.75 ^ undo ;
bold = +1.0 - bnew ;
else
bnew = +0.0 ;
bold = +1.0 - bnew ;
end
vert(bvrt,:) = ...
bold * vold(bvrt,:) ...
+ bnew * vert(bvrt,:) ;
btri = any( ...
bvrt(tria(:,1:3)),2) ;
end
oscr = nscr ;
tcpu.undo = ...
tcpu.undo + toc(ttic) ;
%------------------------------------- test convergence!
ttic = tic ;
vdel = ...
sum((vert-vold).^2,2) ;
evec = vert(edge(:,2),:) ...
- vert(edge(:,1),:) ;
elen = ...
sqrt(sum(evec.^2,2)) ;
hvrt = evalhfn(vert, ...
edge,EMAT,hfun,harg) ;
hmid = hvrt(edge(:,1),:) ...
+ hvrt(edge(:,2),:) ;
hmid = hmid * 0.5 ;
scal = elen./hmid ;
emid = vert(edge(:,1),:) ...
+ vert(edge(:,2),:) ;
emid = emid * 0.5 ;
%------------------------------------- |deg|-based prune
keep = false(size(vert,1),1);
keep(vdeg>+4) = true ;
keep(conn(:)) = true ;
keep(free(:)) = true ;
%------------------------------------- 'density' control
lmax = +5. / +4. ;
lmin = +1. / lmax ;
less = scal<=lmin ;
more = scal>=lmax ;
vbnd = false(size(vert,1),1);
vbnd(conn(:,1)) = true ;
vbnd(conn(:,2)) = true ;
ebad = vbnd(edge(:,1)) ... %-- not at boundaries
| vbnd(edge(:,2)) ;
less(ebad(:)) = false;
more(ebad(:)) = false;
%------------------------------------- force as disjoint
lidx = find (less) ;
for lpos = 1 : length(lidx)
epos = lidx(lpos,1);
inod = edge(epos,1);
jnod = edge(epos,2);
%--------------------------------- if still disjoint
if (keep(inod) && ...
keep(jnod) )
keep(inod) = false ;
keep(jnod) = false ;
else
less(epos) = false ;
end
end
ebad = ...
keep(edge(less,1)) ...
& keep(edge(less,2)) ;
more(ebad(:)) = false ;
%------------------------------------- reindex vert/tria
redo = ...
zeros(size(vert,1),1);
itop = ...
length(find(keep));
iend = ...
length(find(less));
redo(keep) = (1:itop)';
redo(edge(less,1)) = ... %-- to new midpoints
(itop+1 : itop+iend)';
redo(edge(less,2)) = ...
(itop+1 : itop+iend)';
vnew =[vert(keep,:) ;
emid(less,:) ;
] ;
tnew = reshape( ...
redo(tria(:,+1:3)),[],3) ;
ttmp = sort(tnew,2) ; %-- filter collapsed
okay = all( ...
diff(ttmp,1,2)~=0,2) ;
okay = ...
okay & ttmp(:,1) > 0 ;
tnew = tnew(okay,:) ;
%------------------------------------- quality preserver
nscr = ...
triscr2 (vnew,tnew) ;
stol = +0.80 ;
tbad = nscr < stol ...
& nscr < oscr(okay) ;
vbad = ...
false(size(vnew,1),1);
vbad(tnew(tbad,:)) = true;
%------------------------------------- filter edge merge
lidx = find (less) ;
ebad = ...
vbad(redo(edge(lidx,1))) | ...
vbad(redo(edge(lidx,2))) ;
less(lidx(ebad)) = false ;
keep(edge(...
lidx(ebad),1:2)) = true ;
%------------------------------------- reindex vert/conn
redo = ...
zeros(size(vert,1),1);
itop = ...
length(find(keep));
iend = ...
length(find(less));
redo(keep) = (1:itop)';
redo(edge(less,1)) = ...
(itop+1 : itop+iend)';
redo(edge(less,2)) = ...
(itop+1 : itop+iend)';
vert =[vert(keep,:);
emid(less,:);
emid(more,:);
] ;
conn = reshape( ...
redo(conn(:,+1:2)),[],2) ;
tcpu.keep = ...
tcpu.keep + toc(ttic) ;
%------------------------------------- build current CDT
ttic = tic ;
[vert,conn,tria,tnum] = ...
deltri2 (vert, ...
conn,node,PSLG,part) ;
tcpu.dtri = ...
tcpu.dtri + toc(ttic) ;
%------------------------------------- dump-out progess!
vdel = vdel./(hvrt.*hvrt) ;
move = vdel > opts.vtol^2 ;
nmov = ...
length(find(move));
ntri = size(tria,1) ;
if (mod(iter,opts.disp)==+0)
fprintf(+1, ...
'%11i %18i %18i\n', ...
[iter,nmov,ntri]) ;
end
%------------------------------------- loop convergence!
if (nmov == +0), break; end
end
tria = tria( :,1:3);
%----------------------------------------- prune unused vert
keep = false(size(vert,1),1);
keep(tria(:)) = true ;
keep(conn(:)) = true ;
redo = zeros(size(vert,1),1);
redo(keep) = ...
(+1:length(find(keep)))';
conn = ...
reshape(redo(conn),[],2);
tria = ...
reshape(redo(tria),[],3);
vert = vert(keep,:);
tcpu.full = ...
tcpu.full + toc(tnow) ;
if (opts.dbug)
%----------------------------------------- print debug timer
fprintf(1,'\n') ;
fprintf(1,' Mesh smoothing timer...\n');
fprintf(1,'\n') ;
fprintf(1, ...
' FULL: %f \n', tcpu.full);
fprintf(1, ...
' DTRI: %f \n', tcpu.dtri);
fprintf(1, ...
' TCON: %f \n', tcpu.tcon);
fprintf(1, ...
' ITER: %f \n', tcpu.iter);
fprintf(1, ...
' UNDO: %f \n', tcpu.undo);
fprintf(1, ...
' KEEP: %f \n', tcpu.keep);
fprintf(1,'\n') ;
end
if (~isinf(opts.disp)), fprintf(1,'\n'); end
end
function [hvrt] = evalhfn(vert,edge,EMAT,hfun,harg)
%EVALHFN eval. the spacing-fun. at mesh vertices.
if (~isempty (hfun))
if (isnumeric(hfun))
hvrt = hfun * ...
ones(size(vert,1),1) ;
else
hvrt = feval( ...
hfun,vert,harg{:}) ;
end
else
%-- no HFUN - HVRT is mean edge-len. at vertices!
evec = vert(edge(:,2),:) - ...
vert(edge(:,1),:) ;
elen = sqrt(sum(evec.^2,2)) ;
hvrt = (EMAT*elen) ...
./ max(sum(EMAT,2),eps) ;
free = true(size(vert,1),1) ;
free(edge(:,1)) = false;
free(edge(:,2)) = false;
hvrt(free) = +inf;
end
end
function [opts] = makeopt(opts)
%MAKEOPT setup the options structure for SMOOTH2.
if (~isfield(opts,'iter'))
opts.iter = +32;
else
if (~isnumeric(opts.iter))
error('smooth2:incorrectInputClass', ...
'Incorrect input class.');
end
if (numel(opts.iter)~= +1)
error('smooth2:incorrectDimensions', ...
'Incorrect input dimensions.') ;
end
if (opts.iter <= +0)
error('smooth2:invalidOptionValues', ...
'Invalid OPT.ITER selection.') ;
end
end
if (~isfield(opts,'disp'))
opts.disp = + 4;
else
if (~isnumeric(opts.disp))
error('smooth2:incorrectInputClass', ...
'Incorrect input class.');
end
if (numel(opts.disp)~= +1)
error('smooth2:incorrectDimensions', ...
'Incorrect input dimensions.') ;
end
if (opts.disp <= +0)
error('smooth2:invalidOptionValues', ...
'Invalid OPT.DISP selection.') ;
end
end
if (~isfield(opts,'vtol'))
opts.vtol = +1.0E-02;
else
if (~isnumeric(opts.vtol))
error('smooth2:incorrectInputClass', ...
'Incorrect input class.');
end
if (numel(opts.vtol)~= +1)
error('smooth2:incorrectDimensions', ...
'Incorrect input dimensions.') ;
end
if (opts.vtol <= 0.)
error('smooth2:invalidOptionValues', ...
'Invalid OPT.VTOL selection.') ;
end
end
if (~isfield(opts,'dbug'))
opts.dbug = false;
else
if (~islogical(opts.dbug))
error('refine2:incorrectInputClass', ...
'Incorrect input class.');
end
if (numel(opts.dbug)~= +1)
error('refine2:incorrectDimensions', ...
'Incorrect input dimensions.') ;
end
end
end