Skip to content

grad_intrinsic

grad_intrinsic(l_sq, F, n=None)

Intrinsic finite element gradient matrix

Given a triangle mesh, computes the finite element gradient matrix assuming piecewise linear hat function basis.

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

Returns:

Name Type Description
G (2

Sparse FEM gradient matrix. The first m rows ar ethe gradient with respect to the edge (1,2). The m rows after that run in a pi/2 rotation counter-clockwise.

See Also

cotangent_laplacian.

Notes

Adapted from https://github.com/alecjacobson/gptoolbox/blob/master/mesh/grad_intrinsic.m

Examples:

import gpytoolbox as gpy
l = gpy.halfedge_lengths_squared(V,F)
G = gpy.grad_intrinsic(l_sq,F)
Source code in src/gpytoolbox/grad_intrinsic.py
 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
def grad_intrinsic(l_sq,F,
    n=None,):
    """Intrinsic finite element gradient matrix

    Given a triangle mesh, computes the finite element gradient matrix assuming
    piecewise linear hat function basis.

    Parameters
    ----------
    l_sq : (m,3) numpy array
        squared halfedge lengths as computed by halfedge_lengths_squared
    F : (m,3) numpy int array
        face index list of a triangle mesh
    n : int, optional (default None)
        number of vertices in the mesh.
        If absent, will try to infer from F.

    Returns
    -------
    G : (2*m,n) scipy sparse.csr_matrix
        Sparse FEM gradient matrix.
        The first m rows ar ethe gradient with respect to the edge (1,2).
        The m rows after that run in a pi/2 rotation counter-clockwise.

    See Also
    --------
    cotangent_laplacian.

    Notes
    -----
    Adapted from https://github.com/alecjacobson/gptoolbox/blob/master/mesh/grad_intrinsic.m

    Examples
    --------
    ```python
    import gpytoolbox as gpy
    l = gpy.halfedge_lengths_squared(V,F)
    G = gpy.grad_intrinsic(l_sq,F)
    ```
    """


    assert F.shape[1] == 3, "Only works on triangle meshes."
    assert l_sq.shape == F.shape
    assert np.all(l_sq>=0)

    if n==None:
        n = np.max(F)+1
    m = F.shape[0]

    l0 = np.sqrt(l_sq[:,0])[:,None]
    gx = (l_sq[:,1][:,None]-l_sq[:,0][:,None]-l_sq[:,2][:,None]) / (-2.*l0)
    gy = np.sqrt(l_sq[:,2][:,None] - gx**2)
    V2 = np.block([[gx,gy], [np.zeros((m,2))], [l0, np.zeros((m,1))]])
    F2 = np.reshape(np.arange(3*m), (m,3), order='F')
    G2 = grad(V2,F2)
    P = sp.sparse.csr_matrix((np.ones(3*m),(F2.ravel(),F.ravel())),
        shape=(3*m,n))
    G = G2*P

    return G