-
Notifications
You must be signed in to change notification settings - Fork 0
Expand file tree
/
Copy pathvem_p2_elasticity.py
More file actions
1728 lines (1435 loc) · 63.3 KB
/
vem_p2_elasticity.py
File metadata and controls
1728 lines (1435 loc) · 63.3 KB
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
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
"""
Second-order (P2) VEM for 2D Linear Elasticity on Polygonal Meshes.
Extension of vem_elasticity.py (P1, lowest-order) to P2 (second-order).
Key change: vertex-only DOFs -> vertex + edge midpoint DOFs.
P1 (vem_elasticity.py):
- DOFs: 2 per vertex node (u_x, u_y)
- Polynomial basis: 6 functions (3 rigid body + 3 strain modes)
- Accuracy: O(h) for H1, O(h^2) for L2
P2 (this file):
- DOFs: 2 per vertex + 2 per edge midpoint = 4*n_v per element
- Polynomial basis: 12 functions (dim P2^2 in 2D = 2 x 6)
- Accuracy: O(h^2) for H1, O(h^3) for L2
P2 scalar monomials (in scaled coordinates xhat, yhat):
{1, xhat, yhat, xhat^2, xhat*yhat, yhat^2} (6 total)
Vector P2 basis (12 total):
(m_i, 0) for i=1..6 and (0, m_i) for i=1..6
Split: 3 rigid body + 9 strain modes.
Boundary integrals use Simpson's rule (vertex + midpoint) instead of trapezoidal.
References:
- Artioli, Beirao da Veiga, Lovadina, Sacco (2017) "Arbitrary order 2D virtual
elements for polygonal meshes: Part I and II" (arXiv:1701.04306)
- Beirao da Veiga et al. (2013) "Basic principles of VEM"
- Gain, Talischi, Paulino (2014) "VEM for 3D elasticity"
- Sutton (2017) "The VEM in 50 lines of MATLAB"
"""
import numpy as np
import scipy.sparse as sp
import matplotlib.pyplot as plt
from matplotlib.patches import Polygon as MplPolygon
from matplotlib.collections import PatchCollection
from scipy.spatial import Voronoi
# ── Mesh Utilities ────────────────────────────────────────────────────────
def generate_voronoi_mesh(n_cells, domain=(0, 1, 0, 1), seed=42):
"""
Generate a clipped Voronoi mesh on a rectangular domain.
Parameters
----------
n_cells : int — approximate number of cells
domain : (xmin, xmax, ymin, ymax)
seed : int — random seed
Returns
-------
vertices : (N, 2) array
elements : list of int arrays (0-based connectivity, vertex-only)
boundary : array of boundary node indices
"""
rng = np.random.RandomState(seed)
xmin, xmax, ymin, ymax = domain
Lx, Ly = xmax - xmin, ymax - ymin
# Generate seed points inside domain
pts = rng.rand(n_cells, 2) * [Lx, Ly] + [xmin, ymin]
# Mirror for bounded Voronoi
pts_mirror = np.vstack([
pts,
np.column_stack([2 * xmin - pts[:, 0], pts[:, 1]]),
np.column_stack([2 * xmax - pts[:, 0], pts[:, 1]]),
np.column_stack([pts[:, 0], 2 * ymin - pts[:, 1]]),
np.column_stack([pts[:, 0], 2 * ymax - pts[:, 1]]),
])
vor = Voronoi(pts_mirror)
# Clip to domain
tol = 1e-12
vertices_list = []
vert_map = {}
elements = []
def add_vertex(v):
key = (round(v[0], 10), round(v[1], 10))
if key not in vert_map:
vert_map[key] = len(vertices_list)
vertices_list.append(v.copy())
return vert_map[key]
for i in range(n_cells):
region_idx = vor.point_region[i]
region = vor.regions[region_idx]
if -1 in region or len(region) < 3:
continue
poly = vor.vertices[region]
# Clip to domain
clipped = _clip_polygon_to_box(poly, xmin, xmax, ymin, ymax)
if clipped is None or len(clipped) < 3:
continue
# Order CCW
cx = clipped[:, 0].mean()
cy = clipped[:, 1].mean()
angles = np.arctan2(clipped[:, 1] - cy, clipped[:, 0] - cx)
order = np.argsort(angles)
clipped = clipped[order]
el_ids = []
for v in clipped:
el_ids.append(add_vertex(v))
elements.append(np.array(el_ids, dtype=int))
vertices = np.array(vertices_list)
# Boundary nodes
boundary = []
for i, v in enumerate(vertices):
if (abs(v[0] - xmin) < tol or abs(v[0] - xmax) < tol or
abs(v[1] - ymin) < tol or abs(v[1] - ymax) < tol):
boundary.append(i)
boundary = np.array(boundary, dtype=int)
return vertices, elements, boundary
def _clip_polygon_to_box(poly, xmin, xmax, ymin, ymax):
"""Sutherland-Hodgman polygon clipping to axis-aligned box."""
output = list(poly)
edges = [
(lambda p: p[0] - xmin, lambda p, q: _intersect_edge(p, q, 0, xmin)),
(lambda p: xmax - p[0], lambda p, q: _intersect_edge(p, q, 0, xmax)),
(lambda p: p[1] - ymin, lambda p, q: _intersect_edge(p, q, 1, ymin)),
(lambda p: ymax - p[1], lambda p, q: _intersect_edge(p, q, 1, ymax)),
]
for inside_fn, intersect_fn in edges:
if len(output) == 0:
return None
inp = output
output = []
for i in range(len(inp)):
curr = inp[i]
prev = inp[i - 1]
curr_in = inside_fn(curr) >= -1e-14
prev_in = inside_fn(prev) >= -1e-14
if curr_in:
if not prev_in:
output.append(intersect_fn(prev, curr))
output.append(curr)
elif prev_in:
output.append(intersect_fn(prev, curr))
if len(output) < 3:
return None
return np.array(output)
def _intersect_edge(p, q, axis, val):
"""Intersection of line segment p->q with axis=val."""
t = (val - p[axis]) / (q[axis] - p[axis] + 1e-30)
result = p + t * (q - p)
result[axis] = val
return result
# ── Edge Midpoint Augmentation ────────────────────────────────────────────
def add_edge_midpoints(vertices, elements):
"""
Augment mesh with edge midpoint nodes for P2 VEM.
For each polygon element, compute midpoint of each edge.
Shared edges between adjacent elements get a single midpoint node.
Parameters
----------
vertices : (N, 2) array — original vertex coordinates
elements : list of int arrays — P1 connectivity (vertex-only)
Returns
-------
new_vertices : (N + N_mid, 2) array — vertices + midpoints
new_elements : list of int arrays — each element: [v0..v_{n-1}, m0..m_{n-1}]
where m_i is the midpoint of edge (v_i, v_{i+1 mod n})
edge_midpoint_map : dict — (min_v, max_v) -> midpoint node index
"""
new_vertices = list(vertices)
edge_midpoint_map = {}
new_elements = []
for el in elements:
el_int = el.astype(int)
n_v = len(el_int)
midpoint_ids = []
for i in range(n_v):
v1 = el_int[i]
v2 = el_int[(i + 1) % n_v]
edge_key = (min(v1, v2), max(v1, v2))
if edge_key not in edge_midpoint_map:
mid = 0.5 * (vertices[v1] + vertices[v2])
mid_idx = len(new_vertices)
new_vertices.append(mid)
edge_midpoint_map[edge_key] = mid_idx
midpoint_ids.append(edge_midpoint_map[edge_key])
# Element DOF ordering: [v0, v1, ..., v_{n-1}, m0, m1, ..., m_{n-1}]
new_el = np.concatenate([el_int, np.array(midpoint_ids, dtype=int)])
new_elements.append(new_el)
new_vertices = np.array(new_vertices)
return new_vertices, new_elements, edge_midpoint_map
# ── P2 Polynomial Basis (P1-compatible ordering) ────────────────────────
#
# P1 basis (first 6, identical to vem_elasticity.py):
# 0: (1, 0) — translation x
# 1: (0, 1) — translation y
# 2: (-yhat, xhat) — rigid rotation (zero strain)
# 3: (xhat, 0) — eps_xx = 1/h
# 4: (0, yhat) — eps_yy = 1/h
# 5: (yhat, xhat) — 2*eps_xy = 2/h (symmetric shear)
#
# P2 extension (6 quadratic modes with non-constant strain):
# 6: (xhat^2, 0) — eps_xx = 2*xhat/h
# 7: (0, yhat^2) — eps_yy = 2*yhat/h
# 8: (xhat*yhat, 0) — eps_xx = yhat/h, 2*eps_xy = xhat/h
# 9: (0, xhat*yhat) — eps_yy = xhat/h, 2*eps_xy = yhat/h
# 10: (yhat^2, 0) — 2*eps_xy = 2*yhat/h
# 11: (0, xhat^2) — 2*eps_xy = 2*xhat/h
#
# Key insight: rotation mode MUST be at index 2 so that the vorticity
# B-row (row 2) matches the D-column (col 2). The old (m_i,0)/(0,m_i)
# interleaved ordering placed (xhat,0) at index 2, which has ZERO
# vorticity, causing G[2,2]=0 → singular G matrix.
#
# For quadratic modes (6-11), div(sigma(p_alpha)) != 0, so the B matrix
# needs a volume correction term in addition to the boundary integral.
def p2_polynomial_basis_2d():
"""
Information about the 12 P2 basis functions for 2D vector elasticity.
P1-compatible ordering: first 6 = P1 basis, next 6 = quadratic extension.
"""
return {
'n_polys': 12,
'n_scalar_monomials': 6,
'n_rigid': 3,
'n_strain': 9,
'basis': [
'(1, 0)', '(0, 1)', '(-yhat, xhat)',
'(xhat, 0)', '(0, yhat)', '(yhat, xhat)',
'(xhat^2, 0)', '(0, yhat^2)', '(xhat*yhat, 0)',
'(0, xhat*yhat)', '(yhat^2, 0)', '(0, xhat^2)',
],
}
def _eval_p2_vector_basis(x, y, xc, yc, h):
"""
Evaluate all 12 P2 vector basis functions at a point.
P1-compatible ordering (see module-level comment).
Returns: (12, 2) array where row i = (p_i^x, p_i^y).
"""
xh = (x - xc) / h
yh = (y - yc) / h
return np.array([
[1.0, 0.0], # 0: (1, 0)
[0.0, 1.0], # 1: (0, 1)
[-yh, xh], # 2: (-yhat, xhat) rotation
[xh, 0.0], # 3: (xhat, 0)
[0.0, yh], # 4: (0, yhat)
[yh, xh], # 5: (yhat, xhat) shear
[xh**2, 0.0], # 6: (xhat^2, 0)
[0.0, yh**2], # 7: (0, yhat^2)
[xh*yh, 0.0], # 8: (xhat*yhat, 0)
[0.0, xh*yh], # 9: (0, xhat*yhat)
[yh**2, 0.0], # 10: (yhat^2, 0)
[0.0, xh**2], # 11: (0, xhat^2)
])
def _eval_p2_strain(x, y, xc, yc, h):
"""
Evaluate strain (Voigt: eps_xx, eps_yy, 2*eps_xy) of each of the 12 P2
vector basis functions at a point. P1-compatible ordering.
Returns: (12, 3) array where row i = strain of p_i in Voigt notation.
"""
xh = (x - xc) / h
yh = (y - yc) / h
ih = 1.0 / h
return np.array([
[0.0, 0.0, 0.0], # 0: (1, 0)
[0.0, 0.0, 0.0], # 1: (0, 1)
[0.0, 0.0, 0.0], # 2: (-yhat, xhat) — rotation, zero strain
[ih, 0.0, 0.0], # 3: (xhat, 0) — eps_xx=1/h
[0.0, ih, 0.0], # 4: (0, yhat) — eps_yy=1/h
[0.0, 0.0, 2*ih], # 5: (yhat, xhat) — 2eps_xy=1/h+1/h=2/h
[2*xh*ih, 0.0, 0.0], # 6: (xhat^2, 0)
[0.0, 2*yh*ih, 0.0], # 7: (0, yhat^2)
[yh*ih, 0.0, xh*ih], # 8: (xhat*yhat, 0)
[0.0, xh*ih, yh*ih], # 9: (0, xhat*yhat)
[0.0, 0.0, 2*yh*ih], # 10: (yhat^2, 0)
[0.0, 0.0, 2*xh*ih], # 11: (0, xhat^2)
])
def _compute_div_sigma(C, h):
"""
Compute div(sigma(p_alpha)) for each of the 12 P2 basis functions.
For linear modes (0-5), div=0. For quadratic modes (6-11), div is constant.
Uses: div(sigma)_x = d(sigma_xx)/dx + d(sigma_xy)/dy
div(sigma)_y = d(sigma_xy)/dx + d(sigma_yy)/dy
Returns: (12, 2) array where row i = [div_x, div_y] of sigma(p_i).
"""
h2 = h * h
C00, C01, C11, C22 = C[0, 0], C[0, 1], C[1, 1], C[2, 2]
div_sigma = np.zeros((12, 2))
# Modes 0-5: constant or zero strain → div(sigma) = 0
# Mode 6: (xhat^2, 0) — sigma_xx = C00*2xhat/h, sigma_yy = C01*2xhat/h, sigma_xy = 0
# div_x = d(C00*2xhat/h)/dx = 2*C00/h^2
div_sigma[6] = [2 * C00 / h2, 0.0]
# Mode 7: (0, yhat^2) — sigma_yy = C11*2yhat/h, sigma_xx = C01*2yhat/h, sigma_xy = 0
# div_y = d(C11*2yhat/h)/dy = 2*C11/h^2
div_sigma[7] = [0.0, 2 * C11 / h2]
# Mode 8: (xhat*yhat, 0) — eps=(yhat/h, 0, xhat/h)
# sigma_xx = C00*yhat/h, sigma_yy = C01*yhat/h, sigma_xy = C22*xhat/h
# div_x = d(sigma_xx)/dx + d(sigma_xy)/dy = 0 + 0 = 0
# div_y = d(sigma_xy)/dx + d(sigma_yy)/dy = C22/h^2 + C01/h^2
div_sigma[8] = [0.0, (C22 + C01) / h2]
# Mode 9: (0, xhat*yhat) — eps=(0, xhat/h, yhat/h)
# sigma_xx = C01*xhat/h, sigma_yy = C11*xhat/h, sigma_xy = C22*yhat/h
# div_x = d(sigma_xx)/dx + d(sigma_xy)/dy = C01/h^2 + C22/h^2 (wait...)
# Actually: sigma_xx = C01*(xhat/h), d(sigma_xx)/dx = C01/h^2
# sigma_xy = C22*(yhat/h), d(sigma_xy)/dy = C22/h^2
# div_x = C01/h^2 + C22/h^2 ... but wait, sigma_xy dep on y, dsigma_xy/dy = C22/h^2? No!
# sigma_xy = C22*yhat/h. d/dy of yhat/h = 1/h^2. So d(sigma_xy)/dy = C22/h^2.
# Actually: d(sigma_xy)/dy is part of div_x? No!
# div_x = d(sigma_xx)/dx + d(sigma_xy)/dy = C01/h^2 + 0 (sigma_xy = C22*yhat/h, d/dy = C22/h^2)
# Wait: sigma_xy depends on yhat. d(sigma_xy)/dy = C22*(1/h)*(1/h) = C22/h^2.
# So div_x = C01/h^2 + C22/h^2? Hmm let me recheck.
# sigma_xx = C[0,1]*eps_yy + ... = C01*(xhat/h). d(sigma_xx)/dx = C01/h^2. ✓
# sigma_xy = C22*(yhat/h). d(sigma_xy)/dy = C22/h^2. ✓
# div_x = C01/h^2 + C22/h^2. Hmm but that seems wrong...
# Let me recompute: eps = (0, xhat/h, yhat/h).
# sigma = C @ eps: sigma_xx = C[0,0]*0 + C[0,1]*xhat/h + C[0,2]*yhat/h = C01*xhat/h
# sigma_yy = C[1,1]*xhat/h, sigma_xy = C[2,2]*yhat/h.
# d(sigma_xx)/dx = C01*(1/h)*(d xhat/dx) = C01/h^2. ✓
# d(sigma_xy)/dy = C22*(1/h)*(d yhat/dy) = C22/h^2. ✓
# div_x = C01/h^2 + C22/h^2 ... but actually this is wrong because
# d(sigma_xy)/dy contributes to div_x only: div_x = d(sigma_xx)/dx + d(sigma_xy)/dy. ✓
# Wait no: div_x = d(sigma_xx)/dx + d(sigma_xy)/dy. For 2D stress:
# equilibrium: d(sigma_xx)/dx + d(sigma_xy)/dy = -f_x (first eq)
# d(sigma_xy)/dx + d(sigma_yy)/dy = -f_y (second eq)
# For mode 9: sigma_xy = C22*yhat/h. d(sigma_xy)/dy = C22/h^2
# So div_x = d(sigma_xx)/dx + d(sigma_xy)/dy = C01/h^2 + C22/h^2. Hmm.
# But sigma_xy depends on y, NOT x. So d(sigma_xy)/dx = 0.
# div_y = d(sigma_xy)/dx + d(sigma_yy)/dy = 0 + 0 = 0
# (sigma_yy = C11*xhat/h depends on x only, d/dy = 0.)
div_sigma[9] = [(C01 + C22) / h2, 0.0]
# Mode 10: (yhat^2, 0) — eps=(0, 0, 2*yhat/h)
# sigma_xy = C22*2*yhat/h, all others 0
# div_x = d(sigma_xx)/dx + d(sigma_xy)/dy = 0 + 2*C22/h^2
div_sigma[10] = [2 * C22 / h2, 0.0]
# Mode 11: (0, xhat^2) — eps=(0, 0, 2*xhat/h)
# sigma_xy = C22*2*xhat/h, all others 0
# div_y = d(sigma_xy)/dx + d(sigma_yy)/dy = 2*C22/h^2 + 0
div_sigma[11] = [0.0, 2 * C22 / h2]
return div_sigma
def _compute_strain_energy_matrix(verts, xc, yc, h, C):
"""
Compute the 12x12 strain energy matrix a_K[alpha, beta] analytically.
a_K[α, β] = ∫_K σ(p_α) : ε(p_β) dK = ∫_K ε(p_α)^T C ε(p_β) dK
Uses sub-triangulation from centroid with 3-point Gauss quadrature
on each triangle (exact for degree 2 polynomials → sufficient for
the product of two linear strains).
Returns: (12, 12) symmetric PSD matrix.
"""
n_polys = 12
n_v = len(verts)
a_K = np.zeros((n_polys, n_polys))
# 3-point Gauss quadrature on reference triangle (exact for degree 2)
# Points: (1/6, 1/6), (2/3, 1/6), (1/6, 2/3), weights: 1/6 each
gauss_pts = np.array([[1.0/6, 1.0/6], [2.0/3, 1.0/6], [1.0/6, 2.0/3]])
gauss_wts = np.array([1.0/6, 1.0/6, 1.0/6])
for i in range(n_v):
j = (i + 1) % n_v
# Triangle: centroid, vertex i, vertex j
t0 = np.array([xc, yc])
t1 = verts[i]
t2 = verts[j]
# Triangle area
tri_area = 0.5 * abs((t1[0]-t0[0])*(t2[1]-t0[1]) - (t2[0]-t0[0])*(t1[1]-t0[1]))
if tri_area < 1e-30:
continue
for gp, gw in zip(gauss_pts, gauss_wts):
# Map from reference triangle to physical
x_g = t0[0]*(1-gp[0]-gp[1]) + t1[0]*gp[0] + t2[0]*gp[1]
y_g = t0[1]*(1-gp[0]-gp[1]) + t1[1]*gp[0] + t2[1]*gp[1]
strain_all = _eval_p2_strain(x_g, y_g, xc, yc, h) # (12, 3)
# a_K[alpha, beta] += w * tri_area * eps_alpha^T @ C @ eps_beta
for alpha in range(n_polys):
eps_a = strain_all[alpha]
if np.all(np.abs(eps_a) < 1e-30):
continue
sig_a = C @ eps_a
for beta in range(alpha, n_polys):
eps_b = strain_all[beta]
val = 2.0 * tri_area * gw * np.dot(sig_a, eps_b)
a_K[alpha, beta] += val
if beta != alpha:
a_K[beta, alpha] += val
return a_K
# ── P2 VEM Core ───────────────────────────────────────────────────────────
def vem_p2_elasticity(vertices, elements, E_field, nu, bc_fixed_dofs, bc_vals,
load_dofs=None, load_vals=None, stabilization_alpha=1.0):
"""
Second-order (P2) VEM for 2D plane-stress linear elasticity.
Parameters
----------
vertices : (N, 2) array — node coordinates (vertices + midpoints)
elements : list of int arrays — P2 connectivity:
[v0, v1, ..., v_{n-1}, m0, m1, ..., m_{n-1}]
First n_v entries are vertices, next n_v are edge midpoints.
E_field : float or (N_el,) array — Young's modulus per element
nu : float — Poisson's ratio
bc_fixed_dofs : array of int — constrained DOF indices (global)
bc_vals : array — prescribed values for fixed DOFs
load_dofs : array of int — DOFs with applied point loads
load_vals : array — load values
stabilization_alpha : float — stabilization parameter (default 1.0)
Returns
-------
u : (2*N,) displacement vector
"""
n_nodes = vertices.shape[0]
n_dofs = 2 * n_nodes
n_polys = 12 # dim P_2^2 in 2D
# Sparse assembly via COO triplets
row_idx = []
col_idx = []
val_data = []
F_global = np.zeros(n_dofs)
for el_id in range(len(elements)):
el_nodes = elements[el_id].astype(int)
n_total = len(el_nodes)
n_v = n_total // 2 # half vertex, half midpoint
n_el_dofs = 2 * n_total # = 4 * n_v
# Vertex and midpoint coordinates
vert_ids = el_nodes[:n_v]
mid_ids = el_nodes[n_v:]
verts = vertices[vert_ids]
mids = vertices[mid_ids]
all_coords = vertices[el_nodes] # (2*n_v, 2)
# ── Element E ──
E_el = E_field[el_id] if hasattr(E_field, '__len__') else E_field
# Plane stress constitutive matrix (Voigt: [sigma_xx, sigma_yy, sigma_xy])
C = (E_el / (1.0 - nu**2)) * np.array([
[1.0, nu, 0.0],
[nu, 1.0, 0.0],
[0.0, 0.0, (1.0 - nu) / 2.0]
])
# ── Geometry (based on vertices only) ──
area_comp = (verts[:, 0] * np.roll(verts[:, 1], -1)
- np.roll(verts[:, 0], -1) * verts[:, 1])
area = 0.5 * abs(np.sum(area_comp))
if area < 1e-30:
continue
centroid = np.sum(
(np.roll(verts, -1, axis=0) + verts) * area_comp[:, None],
axis=0) / (6.0 * area)
# Diameter
h = max(np.linalg.norm(verts[i] - verts[j])
for i in range(n_v) for j in range(i + 1, n_v))
if h < 1e-30:
continue
xc, yc = centroid
# ── D matrix (n_el_dofs x 12) ──
D = np.zeros((n_el_dofs, n_polys))
for i in range(n_total):
x_i, y_i = all_coords[i]
basis_vals = _eval_p2_vector_basis(x_i, y_i, xc, yc, h)
D[2 * i, :] = basis_vals[:, 0]
D[2 * i + 1, :] = basis_vals[:, 1]
# ── B matrix (12 x n_el_dofs) ──
B = np.zeros((n_polys, n_el_dofs))
# --- Rows 0-1: translations (average displacement) ---
for i in range(n_total):
B[0, 2 * i] = 1.0 / n_total
B[1, 2 * i + 1] = 1.0 / n_total
# --- Row 2: rigid rotation via vorticity boundary integral ---
# 2*area*omega_avg = oint (u_y*n_x - u_x*n_y) ds
for i in range(n_v):
j = (i + 1) % n_v
nx = verts[j][1] - verts[i][1]
ny = verts[i][0] - verts[j][0]
edge_len = np.sqrt(nx**2 + ny**2)
if edge_len < 1e-30:
continue
nx_u, ny_u = nx / edge_len, ny / edge_len
w_v = edge_len / 6.0
w_m = 4.0 * edge_len / 6.0
B[2, 2 * i] += -ny_u * w_v
B[2, 2 * i + 1] += nx_u * w_v
B[2, 2 * (n_v + i)] += -ny_u * w_m
B[2, 2 * (n_v + i) + 1] += nx_u * w_m
B[2, 2 * j] += -ny_u * w_v
B[2, 2 * j + 1] += nx_u * w_v
B[2, :] /= (2.0 * area)
# --- Rows 3-11: strain modes via boundary integrals ---
# B[alpha, :] . u = oint sigma(p_alpha) . n . u_h ds
# For quadratic modes (6-11), we also need volume correction:
# - oint ... = int_K eps(u_h) : C : eps(p_alpha) dK + int_K div(sigma(p_alpha)) . u_h dK
# So: int_K eps:C:eps dK = oint ... - int_K div(sigma) . u_h dK
# The B row should give the strain energy, so:
# B[alpha, :] . u = boundary_integral - volume_correction
# Boundary integral part (Simpson's rule on each edge)
for i in range(n_v):
j = (i + 1) % n_v
nx = verts[j][1] - verts[i][1]
ny = verts[i][0] - verts[j][0]
edge_len = np.sqrt(nx**2 + ny**2)
if edge_len < 1e-30:
continue
nx_u, ny_u = nx / edge_len, ny / edge_len
pts = [verts[i], mids[i], verts[j]]
dof_indices = [i, n_v + i, j]
sw = [edge_len / 6.0, 4.0 * edge_len / 6.0, edge_len / 6.0]
for pt_idx in range(3):
px, py = pts[pt_idx]
node_idx = dof_indices[pt_idx]
w = sw[pt_idx]
strain_all = _eval_p2_strain(px, py, xc, yc, h)
for alpha in range(9):
poly_idx = 3 + alpha
eps_voigt = strain_all[poly_idx]
sigma = C @ eps_voigt
tx = sigma[0] * nx_u + sigma[2] * ny_u
ty = sigma[2] * nx_u + sigma[1] * ny_u
B[poly_idx, 2 * node_idx] += w * tx
B[poly_idx, 2 * node_idx + 1] += w * ty
# Volume correction for quadratic modes (6-11):
# B[alpha, :] . u should = boundary - int_K div(sigma) . u_h dK
# Approximate int_K div(sigma) . u_h dK ≈ div_sigma . (area/n_total) * sum(u_i)
div_sigma = _compute_div_sigma(C, h)
for alpha_idx in range(6, 12):
dx, dy = div_sigma[alpha_idx]
if abs(dx) + abs(dy) > 1e-30:
for i in range(n_total):
B[alpha_idx, 2 * i] -= dx * area / n_total
B[alpha_idx, 2 * i + 1] -= dy * area / n_total
# ── Projector ──
G = B @ D # 12 x 12
cond = np.linalg.cond(G)
if cond > 1e12:
G += 1e-10 * np.eye(n_polys)
projector = np.linalg.solve(G, B)
# Consistency: use analytically computed strain energy (guaranteed PSD)
a_K_matrix = _compute_strain_energy_matrix(verts, xc, yc, h, C)
K_cons = projector.T @ a_K_matrix @ projector
# Stabilization (Wriggers projection)
I_minus_PiD = np.eye(n_el_dofs) - D @ projector
trace_k = np.trace(K_cons)
stab_param = stabilization_alpha * trace_k / n_el_dofs if trace_k > 0 else E_el
K_stab = stab_param * (I_minus_PiD.T @ I_minus_PiD)
K_local = K_cons + K_stab
# ── Assemble (sparse triplet) ──
gdofs = np.zeros(n_el_dofs, dtype=int)
for i in range(n_total):
gdofs[2 * i] = 2 * el_nodes[i]
gdofs[2 * i + 1] = 2 * el_nodes[i] + 1
ii, jj = np.meshgrid(gdofs, gdofs, indexing='ij')
row_idx.append(ii.ravel())
col_idx.append(jj.ravel())
val_data.append(K_local.ravel())
# Build sparse global stiffness matrix
row_idx = np.concatenate(row_idx)
col_idx = np.concatenate(col_idx)
val_data = np.concatenate(val_data)
K_global = sp.csr_matrix((val_data, (row_idx, col_idx)), shape=(n_dofs, n_dofs))
# ── Point loads ──
if load_dofs is not None and load_vals is not None:
F_global[load_dofs] += load_vals
# ── Solve with BCs ──
u = np.zeros(n_dofs)
bc_set = set(bc_fixed_dofs)
internal = np.array([i for i in range(n_dofs) if i not in bc_set])
u[bc_fixed_dofs] = bc_vals
F_global -= K_global[:, bc_fixed_dofs].toarray() @ bc_vals
K_ii = K_global[np.ix_(internal, internal)]
u[internal] = sp.linalg.spsolve(K_ii, F_global[internal])
return u
# ── Visualization ─────────────────────────────────────────────────────────
def plot_p2_elasticity(vertices, elements, u, field='magnitude',
deform_scale=1.0, title=None, ax=None, save=None):
"""Plot deformed P2 mesh colored by displacement field (vertex cells only)."""
ux = u[0::2]
uy = u[1::2]
deformed = vertices + deform_scale * np.column_stack([ux, uy])
if field == 'magnitude':
vals = np.sqrt(ux**2 + uy**2)
cbar_label = '|u|'
elif field == 'ux':
vals = ux
cbar_label = '$u_x$'
else:
vals = uy
cbar_label = '$u_y$'
own_fig = ax is None
if own_fig:
fig, ax = plt.subplots(1, 1, figsize=(7, 6))
else:
fig = ax.figure
patches = []
patch_colors = []
for el in elements:
el_int = el.astype(int)
n_v = len(el_int) // 2
vert_ids = el_int[:n_v] # Only use vertices for polygon shape
poly = MplPolygon(deformed[vert_ids], closed=True)
patches.append(poly)
# Color by average displacement of all nodes
patch_colors.append(np.mean(vals[el_int]))
pc = PatchCollection(patches, cmap='viridis', edgecolor='k', linewidth=0.3)
pc.set_array(np.array(patch_colors))
ax.add_collection(pc)
ax.set_xlim(deformed[:, 0].min() - 0.05, deformed[:, 0].max() + 0.05)
ax.set_ylim(deformed[:, 1].min() - 0.05, deformed[:, 1].max() + 0.05)
ax.set_aspect('equal')
fig.colorbar(pc, ax=ax, label=cbar_label, shrink=0.8)
if title:
ax.set_title(title)
if save and own_fig:
plt.tight_layout()
plt.savefig(save, dpi=150, bbox_inches='tight')
print(f" Saved: {save}")
plt.close()
return ax
# ── Demo: P2 vs P1 Comparison ────────────────────────────────────────────
def demo_p2_vs_p1():
"""
Compare P1 and P2 VEM on a Voronoi mesh with DI-dependent E field.
- 20-cell Voronoi mesh, domain (0,1)x(0,1)
- E(DI) with spatial gradient
- Bottom fixed, top loaded
- 2x2 comparison figure
"""
print("=" * 60)
print("Demo: P2 vs P1 VEM Comparison")
print("=" * 60)
import vem_elasticity as p1_module
# ── Generate Voronoi mesh ──
vertices, elements, boundary = generate_voronoi_mesh(20, seed=42)
n_nodes_p1 = vertices.shape[0]
n_els = len(elements)
print(f" Mesh: {n_els} elements, {n_nodes_p1} P1 nodes")
# ── E(DI) field ──
E_max = 1000.0
E_min = 30.0
n_hill = 2
nu = 0.3
center = np.array([0.5, 0.5])
max_dist = 0.5 * np.sqrt(2)
E_per_element = np.zeros(n_els)
for i, el in enumerate(elements):
el_int = el.astype(int)
el_centroid = vertices[el_int].mean(axis=0)
dist = np.linalg.norm(el_centroid - center)
DI = 0.9 - 0.8 * (dist / max_dist)
DI = np.clip(DI, 0.05, 0.95)
E_per_element[i] = E_min + (E_max - E_min) * (1.0 - DI) ** n_hill
print(f" E range: [{E_per_element.min():.0f}, {E_per_element.max():.0f}] Pa")
# ── Boundary conditions ──
tol = 1e-6
# --- P1 solve ---
print(" Solving P1...")
bottom_p1 = np.where(vertices[:, 1] < tol)[0]
bc_dofs_p1 = np.concatenate([2 * bottom_p1, 2 * bottom_p1 + 1])
bc_vals_p1 = np.zeros(len(bc_dofs_p1))
top_p1 = np.where(vertices[:, 1] > 1.0 - tol)[0]
load_per_node_p1 = -0.5 / max(len(top_p1), 1)
load_dofs_p1 = 2 * top_p1 + 1
load_vals_p1 = np.full(len(top_p1), load_per_node_p1)
u_p1 = p1_module.vem_elasticity(
vertices, elements, E_per_element, nu,
bc_dofs_p1, bc_vals_p1, load_dofs_p1, load_vals_p1)
ux_p1 = u_p1[0::2]
uy_p1 = u_p1[1::2]
mag_p1 = np.sqrt(ux_p1**2 + uy_p1**2)
print(f" P1 max |u|: {mag_p1.max():.6f}")
# --- P2 solve ---
print(" Solving P2...")
p2_vertices, p2_elements, edge_map = add_edge_midpoints(vertices, elements)
n_nodes_p2 = p2_vertices.shape[0]
print(f" P2 nodes: {n_nodes_p2} ({n_nodes_p1} vertices + {n_nodes_p2 - n_nodes_p1} midpoints)")
# BCs for P2: fix all nodes on bottom edge (vertices + midpoints)
bottom_p2 = np.where(p2_vertices[:, 1] < tol)[0]
bc_dofs_p2 = np.concatenate([2 * bottom_p2, 2 * bottom_p2 + 1])
bc_vals_p2 = np.zeros(len(bc_dofs_p2))
top_p2 = np.where(p2_vertices[:, 1] > 1.0 - tol)[0]
load_per_node_p2 = -0.5 / max(len(top_p2), 1)
load_dofs_p2 = 2 * top_p2 + 1
load_vals_p2 = np.full(len(top_p2), load_per_node_p2)
u_p2 = vem_p2_elasticity(
p2_vertices, p2_elements, E_per_element, nu,
bc_dofs_p2, bc_vals_p2, load_dofs_p2, load_vals_p2)
ux_p2 = u_p2[0::2]
uy_p2 = u_p2[1::2]
mag_p2 = np.sqrt(ux_p2**2 + uy_p2**2)
print(f" P2 max |u|: {mag_p2.max():.6f}")
# ── 2x2 comparison plot ──
fig, axes = plt.subplots(2, 2, figsize=(14, 12))
# (a) P1 deformed mesh + |u| colormap
ax = axes[0, 0]
deformed_p1 = vertices + 200 * np.column_stack([ux_p1, uy_p1])
patches_p1 = []
colors_p1 = []
for el in elements:
el_int = el.astype(int)
poly = MplPolygon(deformed_p1[el_int], closed=True)
patches_p1.append(poly)
colors_p1.append(np.mean(mag_p1[el_int]))
pc1 = PatchCollection(patches_p1, cmap='viridis', edgecolor='k', linewidth=0.3)
pc1.set_array(np.array(colors_p1))
ax.add_collection(pc1)
ax.set_xlim(deformed_p1[:, 0].min() - 0.05, deformed_p1[:, 0].max() + 0.05)
ax.set_ylim(deformed_p1[:, 1].min() - 0.05, deformed_p1[:, 1].max() + 0.05)
ax.set_aspect('equal')
ax.set_title('(a) P1 VEM: deformed mesh + |u|')
fig.colorbar(pc1, ax=ax, label='|u|', shrink=0.8)
# (b) P2 deformed mesh + |u| colormap (show vertex polygons only)
ax = axes[0, 1]
deformed_p2 = p2_vertices + 200 * np.column_stack([ux_p2, uy_p2])
patches_p2 = []
colors_p2 = []
for el in p2_elements:
el_int = el.astype(int)
n_v = len(el_int) // 2
vert_ids = el_int[:n_v]
poly = MplPolygon(deformed_p2[vert_ids], closed=True)
patches_p2.append(poly)
colors_p2.append(np.mean(mag_p2[el_int]))
pc2 = PatchCollection(patches_p2, cmap='viridis', edgecolor='k', linewidth=0.3)
pc2.set_array(np.array(colors_p2))
ax.add_collection(pc2)
ax.set_xlim(deformed_p2[:, 0].min() - 0.05, deformed_p2[:, 0].max() + 0.05)
ax.set_ylim(deformed_p2[:, 1].min() - 0.05, deformed_p2[:, 1].max() + 0.05)
ax.set_aspect('equal')
ax.set_title('(b) P2 VEM: deformed mesh + |u|')
fig.colorbar(pc2, ax=ax, label='|u|', shrink=0.8)
# (c) Difference field |u_P2 - u_P1_interp|
# Interpolate P1 to P2 nodes: P1 vertices map directly, midpoints = avg of endpoints
ax = axes[1, 0]
u_p1_interp = np.zeros(2 * n_nodes_p2)
# Vertex nodes: direct copy
u_p1_interp[:2 * n_nodes_p1] = u_p1.copy()
# Midpoint nodes: average of endpoint values
for (v1, v2), mid_idx in edge_map.items():
u_p1_interp[2 * mid_idx] = 0.5 * (u_p1[2 * v1] + u_p1[2 * v2])
u_p1_interp[2 * mid_idx + 1] = 0.5 * (u_p1[2 * v1 + 1] + u_p1[2 * v2 + 1])
diff_ux = ux_p2 - u_p1_interp[0::2]
diff_uy = uy_p2 - u_p1_interp[1::2]
diff_mag = np.sqrt(diff_ux**2 + diff_uy**2)
patches_diff = []
colors_diff = []
for el in p2_elements:
el_int = el.astype(int)
n_v = len(el_int) // 2
vert_ids = el_int[:n_v]
poly = MplPolygon(vertices[vert_ids], closed=True)
patches_diff.append(poly)
colors_diff.append(np.mean(diff_mag[el_int]))
pc3 = PatchCollection(patches_diff, cmap='hot_r', edgecolor='k', linewidth=0.3)
pc3.set_array(np.array(colors_diff))
ax.add_collection(pc3)
ax.set_xlim(-0.05, 1.05)
ax.set_ylim(-0.05, 1.05)
ax.set_aspect('equal')
ax.set_title('(c) Difference |u_P2 - u_P1_interp|')
fig.colorbar(pc3, ax=ax, label='|u_P2 - u_P1|', shrink=0.8)
# (d) Stress comparison: sigma_xx along y=0.5
ax = axes[1, 1]
# Sample stress at element centroids near y=0.5
y_target = 0.5
y_tol = 0.15
x_vals_p1, sxx_p1 = _sample_stress_along_line(
vertices, elements, u_p1, E_per_element, nu, y_target, y_tol, order=1)
x_vals_p2, sxx_p2 = _sample_stress_along_line(
p2_vertices, p2_elements, u_p2, E_per_element, nu, y_target, y_tol, order=2)
if len(x_vals_p1) > 0:
sort_p1 = np.argsort(x_vals_p1)
ax.plot(x_vals_p1[sort_p1], sxx_p1[sort_p1], 'o-', label='P1', markersize=4)
if len(x_vals_p2) > 0:
sort_p2 = np.argsort(x_vals_p2)
ax.plot(x_vals_p2[sort_p2], sxx_p2[sort_p2], 's--', label='P2', markersize=4)
ax.set_xlabel('x')
ax.set_ylabel(r'$\sigma_{xx}$')
ax.set_title(r'(d) $\sigma_{xx}$ along y=0.5')
ax.legend()
ax.grid(True, alpha=0.3)
fig.suptitle('P2 vs P1 VEM: Biofilm E(DI) Comparison', fontsize=14, fontweight='bold')
plt.tight_layout()
import os
save_dir = os.path.join(os.path.dirname(__file__), 'results')
os.makedirs(save_dir, exist_ok=True)
path = os.path.join(save_dir, 'vem_p2_comparison.png')
plt.savefig(path, dpi=150, bbox_inches='tight')
print(f" Saved: {path}")
plt.close()
def _sample_stress_along_line(vertices, elements, u, E_field, nu, y_target, y_tol,
order=1):
"""
Sample sigma_xx at element centroids near y=y_target.
Parameters
----------
order : 1 for P1 elements, 2 for P2 elements
"""
x_vals = []
sxx_vals = []
for el_id, el in enumerate(elements):
el_int = el.astype(int)
if order == 2:
n_v = len(el_int) // 2
vert_ids = el_int[:n_v]
else:
n_v = len(el_int)
vert_ids = el_int
verts = vertices[vert_ids]
centroid = verts.mean(axis=0)
if abs(centroid[1] - y_target) > y_tol:
continue
E_el = E_field[el_id] if hasattr(E_field, '__len__') else E_field
C = (E_el / (1.0 - nu**2)) * np.array([
[1.0, nu, 0.0], [nu, 1.0, 0.0], [0.0, 0.0, (1.0 - nu) / 2.0]])
# Compute average strain from boundary integral (vertex normals)
area_comp = (verts[:, 0] * np.roll(verts[:, 1], -1)
- np.roll(verts[:, 0], -1) * verts[:, 1])
area = 0.5 * abs(np.sum(area_comp))
if area < 1e-30:
continue
# Average strain via boundary integral
eps_avg = np.zeros(3) # [eps_xx, eps_yy, 2*eps_xy]
if order == 1:
# P1: trapezoidal
for i in range(n_v):
j = (i + 1) % n_v
vi_id = vert_ids[i]
vj_id = vert_ids[j]
nx = verts[j, 1] - verts[i, 1]
ny = verts[i, 0] - verts[j, 0]
ux_avg = 0.5 * (u[2*vi_id] + u[2*vj_id])
uy_avg = 0.5 * (u[2*vi_id+1] + u[2*vj_id+1])
eps_avg[0] += ux_avg * nx
eps_avg[1] += uy_avg * ny
eps_avg[2] += ux_avg * ny + uy_avg * nx
else:
# P2: Simpson with midpoints
mid_ids = el_int[n_v:]
for i in range(n_v):
j = (i + 1) % n_v
vi_id = vert_ids[i]
vj_id = vert_ids[j]
mi_id = mid_ids[i]
nx = verts[j, 1] - verts[i, 1]
ny = verts[i, 0] - verts[j, 0]
edge_len = np.sqrt(nx**2 + ny**2)
if edge_len < 1e-30:
continue
nx_u = nx / edge_len
ny_u = ny / edge_len
# Simpson weights
w1 = edge_len / 6.0
w4 = 4.0 * edge_len / 6.0
w2 = edge_len / 6.0
# eps_xx contribution: u_x * n_x
eps_avg[0] += nx_u * (w1 * u[2*vi_id] + w4 * u[2*mi_id] + w2 * u[2*vj_id])
# eps_yy contribution: u_y * n_y
eps_avg[1] += ny_u * (w1 * u[2*vi_id+1] + w4 * u[2*mi_id+1] + w2 * u[2*vj_id+1])
# 2*eps_xy: u_x*n_y + u_y*n_x
eps_avg[2] += ny_u * (w1 * u[2*vi_id] + w4 * u[2*mi_id] + w2 * u[2*vj_id])
eps_avg[2] += nx_u * (w1 * u[2*vi_id+1] + w4 * u[2*mi_id+1] + w2 * u[2*vj_id+1])
eps_avg /= area
sigma = C @ eps_avg
x_vals.append(centroid[0])