xref: /libCEED/julia/LibCEED.jl/test/rundevtests.jl (revision ba7b2a00d73849c86e4657fa7cfd7a31da27058f)
1using Test, LibCEED, LinearAlgebra, StaticArrays
2
3include("buildmats.jl")
4
5function checkoutput(str, fname)
6    if str != getoutput(fname)
7        write(fname, str)
8        return false
9    end
10    return true
11end
12
13@testset "LibCEED Development Tests" begin
14    @testset "Basis" begin
15        c = Ceed()
16        dim = 3
17        ncomp = 1
18        p = 4
19        q = 6
20        b1 = create_tensor_h1_lagrange_basis(c, dim, ncomp, p, q, GAUSS_LOBATTO)
21
22        @test checkoutput(showstr(b1), "b1.out")
23
24        b1d = CeedScalar[1.0 0.0; 0.5 0.5; 0.0 1.0]
25        d1d = CeedScalar[-0.5 0.5; -0.5 0.5; -0.5 0.5]
26        q1d = CeedScalar[-1.0, 0.0, 1.0]
27        w1d = CeedScalar[1/3, 4/3, 1/3]
28        q, p = size(b1d)
29        d2d = zeros(CeedScalar, 2, q*q, p*p)
30        d2d[1, :, :] = kron(b1d, d1d)
31        d2d[2, :, :] = kron(d1d, b1d)
32
33        dim2 = 2
34        b2 = create_tensor_h1_basis(c, dim2, 1, p, q, b1d, d1d, q1d, w1d)
35        @test checkoutput(showstr(b2), "b2.out")
36
37        b3 = create_h1_basis(
38            c,
39            LINE,
40            1,
41            p,
42            q,
43            b1d,
44            reshape(d1d, 1, q, p),
45            reshape(q1d, 1, q),
46            w1d,
47        )
48        @test checkoutput(showstr(b3), "b3.out")
49
50        dim = 2
51        ncomp = 1
52        p1 = 4
53        q1 = 4
54        qref1 = Array{Float64}(undef, dim, q1)
55        qweight1 = Array{Float64}(undef, q1)
56        interp1, div1 = build_mats_hdiv(qref1, qweight1)
57        b1 = create_hdiv_basis(c, QUAD, ncomp, p1, q1, interp1, div1, qref1, qweight1)
58
59        u1 = ones(Float64, p1)
60        v1 = apply(b1, u1)
61
62        for i = 1:q1
63            @test v1[i] ≈ -1.0
64            @test v1[q1+i] ≈ 1.0
65        end
66
67        p2 = 3
68        q2 = 4
69        qref2 = Array{Float64}(undef, dim, q2)
70        qweight2 = Array{Float64}(undef, q2)
71        interp2, curl2 = build_mats_hcurl(qref2, qweight2)
72        b2 = create_hcurl_basis(c, TRIANGLE, ncomp, p2, q2, interp2, curl2, qref2, qweight2)
73
74        u2 = [1.0, 2.0, 1.0]
75        v2 = apply(b2, u2)
76
77        for i = 1:q2
78            @test v2[i] ≈ 1.0
79        end
80
81        u2[1] = -1.0
82        u2[2] = 1.0
83        u2[3] = 2.0
84        v2 = apply(b2, u2)
85
86        for i = 1:q2
87            @test v2[q2+i] ≈ 1.0
88        end
89    end
90end
91