xref: /libCEED/python/tests/test-3-basis.py (revision d4cc18453651bd0f94c1a2e078b2646a92dafdcc)
1*9ba83ac0SJeremy L Thompson# Copyright (c) 2017-2026, Lawrence Livermore National Security, LLC and other CEED contributors
23d8e8822SJeremy L Thompson# All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
30ef72598Sjeremylt#
43d8e8822SJeremy L Thompson# SPDX-License-Identifier: BSD-2-Clause
50ef72598Sjeremylt#
63d8e8822SJeremy L Thompson# This file is part of CEED:  http://github.com/ceed
70ef72598Sjeremylt
80ef72598Sjeremylt# @file
90ef72598Sjeremylt# Test Ceed Basis functionality
100ef72598Sjeremylt
110ef72598Sjeremyltimport os
120ef72598Sjeremyltimport math
130ef72598Sjeremyltimport libceed
140ef72598Sjeremyltimport numpy as np
150ef72598Sjeremyltimport buildmats as bm
160ef72598Sjeremyltimport check
170ef72598Sjeremylt
1880a9ef05SNatalie BeamsTOL = libceed.EPSILON * 256
190ef72598Sjeremylt
200ef72598Sjeremylt# -------------------------------------------------------------------------------
210ef72598Sjeremylt# Utilities
220ef72598Sjeremylt# -------------------------------------------------------------------------------
230ef72598Sjeremylt
240ef72598Sjeremylt
250ef72598Sjeremyltdef eval(dim, x):
260ef72598Sjeremylt    result, center = 1, 0.1
270ef72598Sjeremylt    for d in range(dim):
280ef72598Sjeremylt        result *= math.tanh(x[d] - center)
290ef72598Sjeremylt        center += 0.1
300ef72598Sjeremylt    return result
310ef72598Sjeremylt
320ef72598Sjeremylt
330ef72598Sjeremyltdef feval(x1, x2):
340ef72598Sjeremylt    return x1 * x1 + x2 * x2 + x1 * x2 + 1
350ef72598Sjeremylt
360ef72598Sjeremylt
370ef72598Sjeremyltdef dfeval(x1, x2):
380ef72598Sjeremylt    return 2 * x1 + x2
390ef72598Sjeremylt
400ef72598Sjeremylt# -------------------------------------------------------------------------------
410ef72598Sjeremylt# Test creation and distruction of a H1Lagrange basis
420ef72598Sjeremylt# -------------------------------------------------------------------------------
430ef72598Sjeremylt
440ef72598Sjeremylt
450ef72598Sjeremyltdef test_300(ceed_resource, capsys):
460ef72598Sjeremylt    ceed = libceed.Ceed(ceed_resource)
470ef72598Sjeremylt
480ef72598Sjeremylt    b = ceed.BasisTensorH1Lagrange(1, 1, 4, 4, libceed.GAUSS_LOBATTO)
490ef72598Sjeremylt    print(b)
500ef72598Sjeremylt    del b
510ef72598Sjeremylt
520ef72598Sjeremylt    b = ceed.BasisTensorH1Lagrange(1, 1, 4, 4, libceed.GAUSS)
530ef72598Sjeremylt    print(b)
540ef72598Sjeremylt    del b
550ef72598Sjeremylt
5680a9ef05SNatalie Beams    # Only run this test in double precision
5780a9ef05SNatalie Beams    if libceed.lib.CEED_SCALAR_TYPE == libceed.SCALAR_FP64:
580ef72598Sjeremylt        stdout, stderr, ref_stdout = check.output(capsys)
590ef72598Sjeremylt        assert not stderr
600ef72598Sjeremylt        assert stdout == ref_stdout
610ef72598Sjeremylt
620ef72598Sjeremylt# -------------------------------------------------------------------------------
630ef72598Sjeremylt# Test GetNumNodes and GetNumQuadraturePoints for basis
640ef72598Sjeremylt# -------------------------------------------------------------------------------
650ef72598Sjeremylt
660ef72598Sjeremylt
670ef72598Sjeremyltdef test_306(ceed_resource):
680ef72598Sjeremylt    ceed = libceed.Ceed(ceed_resource)
690ef72598Sjeremylt
700ef72598Sjeremylt    b = ceed.BasisTensorH1Lagrange(3, 1, 4, 5, libceed.GAUSS_LOBATTO)
710ef72598Sjeremylt
720ef72598Sjeremylt    p = b.get_num_nodes()
730ef72598Sjeremylt    q = b.get_num_quadrature_points()
740ef72598Sjeremylt
750ef72598Sjeremylt    assert p == 64
760ef72598Sjeremylt    assert q == 125
770ef72598Sjeremylt
780ef72598Sjeremylt# -------------------------------------------------------------------------------
790ef72598Sjeremylt# Test interpolation in multiple dimensions
800ef72598Sjeremylt# -------------------------------------------------------------------------------
810ef72598Sjeremylt
820ef72598Sjeremylt
830ef72598Sjeremyltdef test_313(ceed_resource):
840ef72598Sjeremylt    ceed = libceed.Ceed(ceed_resource)
850ef72598Sjeremylt
860ef72598Sjeremylt    for dim in range(1, 4):
870ef72598Sjeremylt        Q = 10
880ef72598Sjeremylt        Qdim = Q**dim
890ef72598Sjeremylt        Xdim = 2**dim
9080a9ef05SNatalie Beams        x = np.empty(Xdim * dim, dtype=ceed.scalar_type())
9180a9ef05SNatalie Beams        uq = np.empty(Qdim, dtype=ceed.scalar_type())
920ef72598Sjeremylt
930ef72598Sjeremylt        for d in range(dim):
940ef72598Sjeremylt            for i in range(Xdim):
950ef72598Sjeremylt                x[d * Xdim + i] = 1 if (i %
960ef72598Sjeremylt                                        (2**(dim - d))) // (2**(dim - d - 1)) else -1
970ef72598Sjeremylt
980ef72598Sjeremylt        X = ceed.Vector(Xdim * dim)
990ef72598Sjeremylt        X.set_array(x, cmode=libceed.USE_POINTER)
1000ef72598Sjeremylt        Xq = ceed.Vector(Qdim * dim)
1010ef72598Sjeremylt        Xq.set_value(0)
1020ef72598Sjeremylt        U = ceed.Vector(Qdim)
1030ef72598Sjeremylt        U.set_value(0)
1040ef72598Sjeremylt        Uq = ceed.Vector(Qdim)
1050ef72598Sjeremylt
1060ef72598Sjeremylt        bxl = ceed.BasisTensorH1Lagrange(dim, dim, 2, Q, libceed.GAUSS_LOBATTO)
1070ef72598Sjeremylt        bul = ceed.BasisTensorH1Lagrange(dim, 1, Q, Q, libceed.GAUSS_LOBATTO)
1080ef72598Sjeremylt
1090ef72598Sjeremylt        bxl.apply(1, libceed.EVAL_INTERP, X, Xq)
1100ef72598Sjeremylt
1110ef72598Sjeremylt        with Xq.array_read() as xq:
1120ef72598Sjeremylt            for i in range(Qdim):
11380a9ef05SNatalie Beams                xx = np.empty(dim, dtype=ceed.scalar_type())
1140ef72598Sjeremylt                for d in range(dim):
1150ef72598Sjeremylt                    xx[d] = xq[d * Qdim + i]
1160ef72598Sjeremylt                uq[i] = eval(dim, xx)
1170ef72598Sjeremylt
1180ef72598Sjeremylt        Uq.set_array(uq, cmode=libceed.USE_POINTER)
1190ef72598Sjeremylt
1200ef72598Sjeremylt        # This operation is the identity because the quadrature is collocated
1210ef72598Sjeremylt        bul.T.apply(1, libceed.EVAL_INTERP, Uq, U)
1220ef72598Sjeremylt
1230ef72598Sjeremylt        bxg = ceed.BasisTensorH1Lagrange(dim, dim, 2, Q, libceed.GAUSS)
1240ef72598Sjeremylt        bug = ceed.BasisTensorH1Lagrange(dim, 1, Q, Q, libceed.GAUSS)
1250ef72598Sjeremylt
1260ef72598Sjeremylt        bxg.apply(1, libceed.EVAL_INTERP, X, Xq)
1270ef72598Sjeremylt        bug.apply(1, libceed.EVAL_INTERP, U, Uq)
1280ef72598Sjeremylt
1290ef72598Sjeremylt        with Xq.array_read() as xq, Uq.array_read() as u:
1300ef72598Sjeremylt            for i in range(Qdim):
13180a9ef05SNatalie Beams                xx = np.empty(dim, dtype=ceed.scalar_type())
1320ef72598Sjeremylt                for d in range(dim):
1330ef72598Sjeremylt                    xx[d] = xq[d * Qdim + i]
1340ef72598Sjeremylt                fx = eval(dim, xx)
1350ef72598Sjeremylt                assert math.fabs(u[i] - fx) < 1E-4
1360ef72598Sjeremylt
1370ef72598Sjeremylt# -------------------------------------------------------------------------------
1380ef72598Sjeremylt# Test grad in multiple dimensions
1390ef72598Sjeremylt# -------------------------------------------------------------------------------
1400ef72598Sjeremylt
1410ef72598Sjeremylt
1420ef72598Sjeremyltdef test_314(ceed_resource):
1430ef72598Sjeremylt    ceed = libceed.Ceed(ceed_resource)
1440ef72598Sjeremylt
1450ef72598Sjeremylt    for dim in range(1, 4):
1460ef72598Sjeremylt        P, Q = 8, 10
1470ef72598Sjeremylt        Pdim = P**dim
1480ef72598Sjeremylt        Qdim = Q**dim
1490ef72598Sjeremylt        Xdim = 2**dim
1500ef72598Sjeremylt        sum1 = sum2 = 0
15180a9ef05SNatalie Beams        x = np.empty(Xdim * dim, dtype=ceed.scalar_type())
15280a9ef05SNatalie Beams        u = np.empty(Pdim, dtype=ceed.scalar_type())
1530ef72598Sjeremylt
1540ef72598Sjeremylt        for d in range(dim):
1550ef72598Sjeremylt            for i in range(Xdim):
1560ef72598Sjeremylt                x[d * Xdim + i] = 1 if (i %
1570ef72598Sjeremylt                                        (2**(dim - d))) // (2**(dim - d - 1)) else -1
1580ef72598Sjeremylt
1590ef72598Sjeremylt        X = ceed.Vector(Xdim * dim)
1600ef72598Sjeremylt        X.set_array(x, cmode=libceed.USE_POINTER)
1610ef72598Sjeremylt        Xq = ceed.Vector(Pdim * dim)
1620ef72598Sjeremylt        Xq.set_value(0)
1630ef72598Sjeremylt        U = ceed.Vector(Pdim)
1640ef72598Sjeremylt        Uq = ceed.Vector(Qdim * dim)
1650ef72598Sjeremylt        Uq.set_value(0)
1660ef72598Sjeremylt        Ones = ceed.Vector(Qdim * dim)
1670ef72598Sjeremylt        Ones.set_value(1)
1680ef72598Sjeremylt        Gtposeones = ceed.Vector(Pdim)
1690ef72598Sjeremylt        Gtposeones.set_value(0)
1700ef72598Sjeremylt
1710ef72598Sjeremylt        # Get function values at quadrature points
1720ef72598Sjeremylt        bxl = ceed.BasisTensorH1Lagrange(dim, dim, 2, P, libceed.GAUSS_LOBATTO)
1730ef72598Sjeremylt        bxl.apply(1, libceed.EVAL_INTERP, X, Xq)
1740ef72598Sjeremylt
1750ef72598Sjeremylt        with Xq.array_read() as xq:
1760ef72598Sjeremylt            for i in range(Pdim):
17780a9ef05SNatalie Beams                xx = np.empty(dim, dtype=ceed.scalar_type())
1780ef72598Sjeremylt                for d in range(dim):
1790ef72598Sjeremylt                    xx[d] = xq[d * Pdim + i]
1800ef72598Sjeremylt                u[i] = eval(dim, xx)
1810ef72598Sjeremylt
1820ef72598Sjeremylt        U.set_array(u, cmode=libceed.USE_POINTER)
1830ef72598Sjeremylt
1840ef72598Sjeremylt        # Calculate G u at quadrature points, G' * 1 at dofs
1850ef72598Sjeremylt        bug = ceed.BasisTensorH1Lagrange(dim, 1, P, Q, libceed.GAUSS)
1860ef72598Sjeremylt        bug.apply(1, libceed.EVAL_GRAD, U, Uq)
1870ef72598Sjeremylt        bug.T.apply(1, libceed.EVAL_GRAD, Ones, Gtposeones)
1880ef72598Sjeremylt
1890ef72598Sjeremylt        # Check if 1' * G * u = u' * (G' * 1)
1900ef72598Sjeremylt        with Gtposeones.array_read() as gtposeones, Uq.array_read() as uq:
1910ef72598Sjeremylt            for i in range(Pdim):
1920ef72598Sjeremylt                sum1 += gtposeones[i] * u[i]
1930ef72598Sjeremylt            for i in range(dim * Qdim):
1940ef72598Sjeremylt                sum2 += uq[i]
1950ef72598Sjeremylt
19680a9ef05SNatalie Beams        assert math.fabs(sum1 - sum2) < 10. * TOL
1970ef72598Sjeremylt
1980ef72598Sjeremylt# -------------------------------------------------------------------------------
1990ef72598Sjeremylt# Test creation and destruction of a 2D Simplex non-tensor H1 basis
2000ef72598Sjeremylt# -------------------------------------------------------------------------------
2010ef72598Sjeremylt
2020ef72598Sjeremylt
2030ef72598Sjeremyltdef test_320(ceed_resource):
2040ef72598Sjeremylt    ceed = libceed.Ceed(ceed_resource)
2050ef72598Sjeremylt
2060ef72598Sjeremylt    P, Q, dim = 6, 4, 2
2070ef72598Sjeremylt
20880a9ef05SNatalie Beams    in_array = np.empty(P, dtype=ceed.scalar_type())
20980a9ef05SNatalie Beams    qref = np.empty(dim * Q, dtype=ceed.scalar_type())
21080a9ef05SNatalie Beams    qweight = np.empty(Q, dtype=ceed.scalar_type())
2110ef72598Sjeremylt
21280a9ef05SNatalie Beams    interp, grad = bm.buildmats(qref, qweight, libceed.scalar_types[
21380a9ef05SNatalie Beams        libceed.lib.CEED_SCALAR_TYPE])
2140ef72598Sjeremylt
2150ef72598Sjeremylt    b = ceed.BasisH1(libceed.TRIANGLE, 1, P, Q, interp, grad, qref, qweight)
2160ef72598Sjeremylt
2170ef72598Sjeremylt    print(b)
2180ef72598Sjeremylt    del b
2190ef72598Sjeremylt
2200ef72598Sjeremylt# -------------------------------------------------------------------------------
2210ef72598Sjeremylt# Test integration with a 2D Simplex non-tensor H1 basis
2220ef72598Sjeremylt# -------------------------------------------------------------------------------
2230ef72598Sjeremylt
2240ef72598Sjeremylt
2250ef72598Sjeremyltdef test_322(ceed_resource):
2260ef72598Sjeremylt    ceed = libceed.Ceed(ceed_resource)
2270ef72598Sjeremylt
2280ef72598Sjeremylt    P, Q, dim = 6, 4, 2
2290ef72598Sjeremylt
2300ef72598Sjeremylt    xr = np.array([0., 0.5, 1., 0., 0.5, 0., 0., 0.,
23180a9ef05SNatalie Beams                   0., 0.5, 0.5, 1.], dtype=ceed.scalar_type())
23280a9ef05SNatalie Beams    in_array = np.empty(P, dtype=ceed.scalar_type())
23380a9ef05SNatalie Beams    qref = np.empty(dim * Q, dtype=ceed.scalar_type())
23480a9ef05SNatalie Beams    qweight = np.empty(Q, dtype=ceed.scalar_type())
2350ef72598Sjeremylt
23680a9ef05SNatalie Beams    interp, grad = bm.buildmats(qref, qweight, libceed.scalar_types[
23780a9ef05SNatalie Beams        libceed.lib.CEED_SCALAR_TYPE])
2380ef72598Sjeremylt
2390ef72598Sjeremylt    b = ceed.BasisH1(libceed.TRIANGLE, 1, P, Q, interp, grad, qref, qweight)
2400ef72598Sjeremylt
2410ef72598Sjeremylt    # Interpolate function to quadrature points
2420ef72598Sjeremylt    for i in range(P):
2430ef72598Sjeremylt        in_array[i] = feval(xr[0 * P + i], xr[1 * P + i])
2440ef72598Sjeremylt
2450ef72598Sjeremylt    in_vec = ceed.Vector(P)
2460ef72598Sjeremylt    in_vec.set_array(in_array, cmode=libceed.USE_POINTER)
2470ef72598Sjeremylt    out_vec = ceed.Vector(Q)
2480ef72598Sjeremylt    out_vec.set_value(0)
2490ef72598Sjeremylt    weights_vec = ceed.Vector(Q)
2500ef72598Sjeremylt    weights_vec.set_value(0)
2510ef72598Sjeremylt
2520ef72598Sjeremylt    b.apply(1, libceed.EVAL_INTERP, in_vec, out_vec)
2530ef72598Sjeremylt    b.apply(1, libceed.EVAL_WEIGHT, libceed.VECTOR_NONE, weights_vec)
2540ef72598Sjeremylt
2550ef72598Sjeremylt    # Check values at quadrature points
2560ef72598Sjeremylt    with out_vec.array_read() as out_array, weights_vec.array_read() as weights_array:
2570ef72598Sjeremylt        sum = 0
2580ef72598Sjeremylt        for i in range(Q):
2590ef72598Sjeremylt            sum += out_array[i] * weights_array[i]
26080a9ef05SNatalie Beams        assert math.fabs(sum - 17. / 24.) < 10. * TOL
2610ef72598Sjeremylt
2620ef72598Sjeremylt# -------------------------------------------------------------------------------
2630ef72598Sjeremylt# Test grad with a 2D Simplex non-tensor H1 basis
2640ef72598Sjeremylt# -------------------------------------------------------------------------------
2650ef72598Sjeremylt
2660ef72598Sjeremylt
2670ef72598Sjeremyltdef test_323(ceed_resource):
2680ef72598Sjeremylt    ceed = libceed.Ceed(ceed_resource)
2690ef72598Sjeremylt
2700ef72598Sjeremylt    P, Q, dim = 6, 4, 2
2710ef72598Sjeremylt
2720ef72598Sjeremylt    xq = np.array([0.2, 0.6, 1. / 3., 0.2, 0.2, 0.2,
27380a9ef05SNatalie Beams                   1. / 3., 0.6], dtype=ceed.scalar_type())
2740ef72598Sjeremylt    xr = np.array([0., 0.5, 1., 0., 0.5, 0., 0., 0.,
27580a9ef05SNatalie Beams                   0., 0.5, 0.5, 1.], dtype=ceed.scalar_type())
27680a9ef05SNatalie Beams    in_array = np.empty(P, dtype=ceed.scalar_type())
27780a9ef05SNatalie Beams    qref = np.empty(dim * Q, dtype=ceed.scalar_type())
27880a9ef05SNatalie Beams    qweight = np.empty(Q, dtype=ceed.scalar_type())
2790ef72598Sjeremylt
28080a9ef05SNatalie Beams    interp, grad = bm.buildmats(qref, qweight, libceed.scalar_types[
28180a9ef05SNatalie Beams        libceed.lib.CEED_SCALAR_TYPE])
2820ef72598Sjeremylt
2830ef72598Sjeremylt    b = ceed.BasisH1(libceed.TRIANGLE, 1, P, Q, interp, grad, qref, qweight)
2840ef72598Sjeremylt
2850ef72598Sjeremylt    # Interpolate function to quadrature points
2860ef72598Sjeremylt    for i in range(P):
2870ef72598Sjeremylt        in_array[i] = feval(xr[0 * P + i], xr[1 * P + i])
2880ef72598Sjeremylt
2890ef72598Sjeremylt    in_vec = ceed.Vector(P)
2900ef72598Sjeremylt    in_vec.set_array(in_array, cmode=libceed.USE_POINTER)
2910ef72598Sjeremylt    out_vec = ceed.Vector(Q * dim)
2920ef72598Sjeremylt    out_vec.set_value(0)
2930ef72598Sjeremylt
2940ef72598Sjeremylt    b.apply(1, libceed.EVAL_GRAD, in_vec, out_vec)
2950ef72598Sjeremylt
2960ef72598Sjeremylt    # Check values at quadrature points
2970ef72598Sjeremylt    with out_vec.array_read() as out_array:
2980ef72598Sjeremylt        for i in range(Q):
2990ef72598Sjeremylt            value = dfeval(xq[0 * Q + i], xq[1 * Q + i])
30080a9ef05SNatalie Beams            assert math.fabs(out_array[0 * Q + i] - value) < 10. * TOL
3010ef72598Sjeremylt
3020ef72598Sjeremylt            value = dfeval(xq[1 * Q + i], xq[0 * Q + i])
30380a9ef05SNatalie Beams            assert math.fabs(out_array[1 * Q + i] - value) < 10. * TOL
3040ef72598Sjeremylt
3050ef72598Sjeremylt# -------------------------------------------------------------------------------
30697c1c57aSSebastian Grimberg# Test interpolation with a 2D Quad non-tensor H(div) basis
30797c1c57aSSebastian Grimberg# -------------------------------------------------------------------------------
30897c1c57aSSebastian Grimberg
30997c1c57aSSebastian Grimberg
31097c1c57aSSebastian Grimbergdef test_330(ceed_resource):
31197c1c57aSSebastian Grimberg    ceed = libceed.Ceed(ceed_resource)
31297c1c57aSSebastian Grimberg
31397c1c57aSSebastian Grimberg    P, Q, dim = 4, 4, 2
31497c1c57aSSebastian Grimberg
31597c1c57aSSebastian Grimberg    in_array = np.ones(P, dtype=ceed.scalar_type())
31697c1c57aSSebastian Grimberg    qref = np.empty(dim * Q, dtype=ceed.scalar_type())
31797c1c57aSSebastian Grimberg    qweight = np.empty(Q, dtype=ceed.scalar_type())
31897c1c57aSSebastian Grimberg
31997c1c57aSSebastian Grimberg    interp, div = bm.buildmatshdiv(qref, qweight, libceed.scalar_types[
32097c1c57aSSebastian Grimberg        libceed.lib.CEED_SCALAR_TYPE])
32197c1c57aSSebastian Grimberg
32297c1c57aSSebastian Grimberg    b = ceed.BasisHdiv(libceed.QUAD, 1, P, Q, interp, div, qref, qweight)
32397c1c57aSSebastian Grimberg
32497c1c57aSSebastian Grimberg    # Interpolate function to quadrature points
32597c1c57aSSebastian Grimberg    in_vec = ceed.Vector(P)
32697c1c57aSSebastian Grimberg    in_vec.set_array(in_array, cmode=libceed.USE_POINTER)
32797c1c57aSSebastian Grimberg    out_vec = ceed.Vector(Q * dim)
32897c1c57aSSebastian Grimberg    out_vec.set_value(0)
32997c1c57aSSebastian Grimberg
33097c1c57aSSebastian Grimberg    b.apply(1, libceed.EVAL_INTERP, in_vec, out_vec)
33197c1c57aSSebastian Grimberg
33297c1c57aSSebastian Grimberg    # Check values at quadrature points
33397c1c57aSSebastian Grimberg    with out_vec.array_read() as out_array:
33497c1c57aSSebastian Grimberg        for i in range(Q):
33597c1c57aSSebastian Grimberg            assert math.fabs(out_array[0 * Q + i] + 1.) < 10. * TOL
33697c1c57aSSebastian Grimberg            assert math.fabs(out_array[1 * Q + i] - 1.) < 10. * TOL
33797c1c57aSSebastian Grimberg
33897c1c57aSSebastian Grimberg# -------------------------------------------------------------------------------
33997c1c57aSSebastian Grimberg# Test interpolation with a 2D Simplex non-tensor H(curl) basis
34097c1c57aSSebastian Grimberg# -------------------------------------------------------------------------------
34197c1c57aSSebastian Grimberg
34297c1c57aSSebastian Grimberg
34397c1c57aSSebastian Grimbergdef test_340(ceed_resource):
34497c1c57aSSebastian Grimberg    ceed = libceed.Ceed(ceed_resource)
34597c1c57aSSebastian Grimberg
34697c1c57aSSebastian Grimberg    P, Q, dim = 3, 4, 2
34797c1c57aSSebastian Grimberg
34897c1c57aSSebastian Grimberg    in_array = np.empty(P, dtype=ceed.scalar_type())
34997c1c57aSSebastian Grimberg    qref = np.empty(dim * Q, dtype=ceed.scalar_type())
35097c1c57aSSebastian Grimberg    qweight = np.empty(Q, dtype=ceed.scalar_type())
35197c1c57aSSebastian Grimberg
35297c1c57aSSebastian Grimberg    interp, curl = bm.buildmatshcurl(qref, qweight, libceed.scalar_types[
35397c1c57aSSebastian Grimberg        libceed.lib.CEED_SCALAR_TYPE])
35497c1c57aSSebastian Grimberg
35597c1c57aSSebastian Grimberg    b = ceed.BasisHcurl(libceed.TRIANGLE, 1, P, Q, interp, curl, qref, qweight)
35697c1c57aSSebastian Grimberg
35797c1c57aSSebastian Grimberg    # Interpolate function to quadrature points
35897c1c57aSSebastian Grimberg    in_array[0] = 1.
35997c1c57aSSebastian Grimberg    in_array[1] = 2.
36097c1c57aSSebastian Grimberg    in_array[2] = 1.
36197c1c57aSSebastian Grimberg
36297c1c57aSSebastian Grimberg    in_vec = ceed.Vector(P)
36397c1c57aSSebastian Grimberg    in_vec.set_array(in_array, cmode=libceed.USE_POINTER)
36497c1c57aSSebastian Grimberg    out_vec = ceed.Vector(Q * dim)
36597c1c57aSSebastian Grimberg    out_vec.set_value(0)
36697c1c57aSSebastian Grimberg
36797c1c57aSSebastian Grimberg    b.apply(1, libceed.EVAL_INTERP, in_vec, out_vec)
36897c1c57aSSebastian Grimberg
36997c1c57aSSebastian Grimberg    # Check values at quadrature points
37097c1c57aSSebastian Grimberg    with out_vec.array_read() as out_array:
37197c1c57aSSebastian Grimberg        for i in range(Q):
37297c1c57aSSebastian Grimberg            assert math.fabs(out_array[0 * Q + i] - 1.) < 10. * TOL
37397c1c57aSSebastian Grimberg
37497c1c57aSSebastian Grimberg    # Interpolate function to quadrature points
37597c1c57aSSebastian Grimberg    in_array[0] = -1.
37697c1c57aSSebastian Grimberg    in_array[1] = 1.
37797c1c57aSSebastian Grimberg    in_array[2] = 2.
37897c1c57aSSebastian Grimberg
37997c1c57aSSebastian Grimberg    out_vec.set_value(0)
38097c1c57aSSebastian Grimberg
38197c1c57aSSebastian Grimberg    b.apply(1, libceed.EVAL_INTERP, in_vec, out_vec)
38297c1c57aSSebastian Grimberg
38397c1c57aSSebastian Grimberg    # Check values at quadrature points
38497c1c57aSSebastian Grimberg    with out_vec.array_read() as out_array:
38597c1c57aSSebastian Grimberg        for i in range(Q):
38697c1c57aSSebastian Grimberg            assert math.fabs(out_array[1 * Q + i] - 1.) < 10. * TOL
38797c1c57aSSebastian Grimberg
38897c1c57aSSebastian Grimberg# -------------------------------------------------------------------------------
389