Skip to content

biharmonic_energy

biharmonic_energy(V, F, 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.

Parameters:

Name Type Description Default
V (n,d) numpy array

vertex list of a triangle mesh

required
F (m,3) numpy int array

face index list of a triangle mesh

required
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". - 'hessian' implements the Hessian energy from Stein et al. 2018. "Natural Boundary Conditions for Smoothing in Geometry Processing". - '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
Q = gpy.biharmonic_energy(V,F)
Source code in src/gpytoolbox/biharmonic_energy.py
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
def biharmonic_energy(V,F,
    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.

    Parameters
    ----------
    V : (n,d) numpy array
        vertex list of a triangle mesh
    F : (m,3) numpy int array
        face index list of a triangle mesh
    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".
         - 'hessian' implements the Hessian energy from Stein et al. 2018.
         "Natural Boundary Conditions for Smoothing in Geometry Processing".
         - 'curved_hessian' implements the Hessian energy from Stein et al.
         2020 "A Smoothness Energy without Boundary Distortion for Curved Surfaces"
         via libigl C++ binding.

    Returns
    -------
    Q : (n,n) scipy csr_matrix
        biharmonic energy matrix

    Examples
    --------
    ```python
    # Mesh in V,F
    import gpytoolbox as gpy
    Q = gpy.biharmonic_energy(V,F)
    ```
    """

    if bc=='hessian':
        return _hessian_energy(V, F)
    else:
        l_sq = halfedge_lengths_squared(V,F)
        return biharmonic_energy_intrinsic(l_sq,F,n=V.shape[0],bc=bc)