biharmonic_energy_intrinsic
biharmonic_energy_intrinsic(l_sq, F, n=None, bc='mixedfem_zero_neumann')
Constructs the biharmonic energy matrix Q such that for a per-vertex function u, the discrete biharmonic energy is u'Qu, using only intrinsic information from the mesh.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
l_sq |
(m,3) numpy array
|
squared halfedge lengths as computed by halfedge_lengths_squared |
required |
F |
(m,3) numpy int array
|
face index list of a triangle mesh |
required |
n |
int, optional (default None)
|
number of vertices in the mesh. If absent, will try to infer from F. |
None
|
bc |
string, optional (default 'mixedfem_zero_neumann')
|
Which type of discretization and boundary condition to use. Options are: {'mixedfem_zero_neumann', 'hessian', 'curved_hessian'} - 'mixedfem_zero_neumann' implements the mixed finite element discretization from Jacobson et al. 2010. "Mixed Finite Elements for Variational Surface Modeling". - 'curved_hessian' implements the Hessian energy from Stein et al. 2020 "A Smoothness Energy without Boundary Distortion for Curved Surfaces" via libigl C++ binding. |
'mixedfem_zero_neumann'
|
Returns:
Name | Type | Description |
---|---|---|
Q |
(n,n) scipy csr_matrix
|
biharmonic energy matrix |
Examples:
# Mesh in V,F
import gpytoolbox as gpy
l_sq = gpy.halfedge_lengths_squared(V,F)
Q = biharmonic_energy_intrinsic(l_sq,F)
Source code in src/gpytoolbox/biharmonic_energy_intrinsic.py
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 |
|