xref: /libCEED/julia/LibCEED.jl/src/Basis.jl (revision 11b88dda510d0aa70e79dc59ad165e2a5539c3c3)
144554ea0SWill Paznerabstract type AbstractBasis end
244554ea0SWill Pazner
344554ea0SWill Pazner"""
444554ea0SWill Pazner    BasisCollocated()
544554ea0SWill Pazner
644554ea0SWill PaznerReturns the singleton object corresponding to libCEED's `CEED_BASIS_COLLOCATED`.
744554ea0SWill Pazner"""
844554ea0SWill Paznerstruct BasisCollocated <: AbstractBasis end
944554ea0SWill PaznerBase.getindex(::BasisCollocated) = C.CEED_BASIS_COLLOCATED[]
1044554ea0SWill Pazner
1144554ea0SWill Pazner"""
1244554ea0SWill Pazner    Basis
1344554ea0SWill Pazner
1444554ea0SWill PaznerWraps a `CeedBasis` object, representing a finite element basis. A `Basis` object can be
1544554ea0SWill Paznercreated using one of:
1644554ea0SWill Pazner
1744554ea0SWill Pazner- [`create_tensor_h1_lagrange_basis`](@ref)
1844554ea0SWill Pazner- [`create_tensor_h1_basis`](@ref)
1944554ea0SWill Pazner- [`create_h1_basis`](@ref)
20*11b88ddaSSebastian Grimberg- [`create_hdiv_basis`](@ref)
21*11b88ddaSSebastian Grimberg- [`create_hcurl_basis`](@ref)
2244554ea0SWill Pazner"""
2344554ea0SWill Paznermutable struct Basis <: AbstractBasis
2444554ea0SWill Pazner    ref::RefValue{C.CeedBasis}
2544554ea0SWill Pazner    function Basis(ref)
2644554ea0SWill Pazner        obj = new(ref)
2744554ea0SWill Pazner        finalizer(obj) do x
2844554ea0SWill Pazner            # ccall(:jl_safe_printf, Cvoid, (Cstring, Cstring), "Finalizing %s.\n", repr(x))
2944554ea0SWill Pazner            destroy(x)
3044554ea0SWill Pazner        end
3144554ea0SWill Pazner        return obj
3244554ea0SWill Pazner    end
3344554ea0SWill Paznerend
3444554ea0SWill Paznerdestroy(b::Basis) = C.CeedBasisDestroy(b.ref) # COV_EXCL_LINE
3544554ea0SWill PaznerBase.getindex(b::Basis) = b.ref[]
3644554ea0SWill PaznerBase.show(io::IO, ::MIME"text/plain", b::Basis) = ceed_show(io, b, C.CeedBasisView)
3744554ea0SWill Pazner
3844554ea0SWill Pazner@doc raw"""
3944554ea0SWill Pazner    create_tensor_h1_lagrange_basis(ceed, dim, ncomp, p, q, qmode)
4044554ea0SWill Pazner
4144554ea0SWill PaznerCreate a tensor-product Lagrange basis.
4244554ea0SWill Pazner
4344554ea0SWill Pazner# Arguments:
4444554ea0SWill Pazner- `ceed`:  A [`Ceed`](@ref) object where the [`Basis`](@ref) will be created.
4544554ea0SWill Pazner- `dim`:   Topological dimension of element.
4644554ea0SWill Pazner- `ncomp`: Number of field components (1 for scalar fields).
4744554ea0SWill Pazner- `p`:     Number of Gauss-Lobatto nodes in one dimension.  The polynomial degree of the
4844554ea0SWill Pazner           resulting $Q_k$ element is $k=p-1$.
4944554ea0SWill Pazner- `q`:     Number of quadrature points in one dimension.
5044554ea0SWill Pazner- `qmode`: Distribution of the $q$ quadrature points (affects order of accuracy for the
5144554ea0SWill Pazner           quadrature).
5244554ea0SWill Pazner"""
5344554ea0SWill Paznerfunction create_tensor_h1_lagrange_basis(c::Ceed, dim, ncomp, p, q, quad_mode::QuadMode)
5444554ea0SWill Pazner    ref = Ref{C.CeedBasis}()
5544554ea0SWill Pazner    C.CeedBasisCreateTensorH1Lagrange(c[], dim, ncomp, p, q, quad_mode, ref)
5644554ea0SWill Pazner    Basis(ref)
5744554ea0SWill Paznerend
5844554ea0SWill Pazner
5944554ea0SWill Pazner@doc raw"""
6044554ea0SWill Pazner    create_tensor_h1_basis(c::Ceed, dim, ncomp, p, q, interp1d, grad1d, qref1d, qweight1d)
6144554ea0SWill Pazner
6244554ea0SWill PaznerCreate a tensor-product basis for $H^1$ discretizations.
6344554ea0SWill Pazner
6444554ea0SWill Pazner# Arguments:
6544554ea0SWill Pazner- `ceed`:      A [`Ceed`](@ref) object where the [`Basis`](@ref) will be created.
6644554ea0SWill Pazner- `dim`:       Topological dimension.
6744554ea0SWill Pazner- `ncomp`:     Number of field components (1 for scalar fields).
6844554ea0SWill Pazner- `p`:         Number of nodes in one dimension.
6944554ea0SWill Pazner- `q`:         Number of quadrature points in one dimension
7044554ea0SWill Pazner- `interp1d`:  Matrix of size `(q, p)` expressing the values of nodal basis functions at
7144554ea0SWill Pazner               quadrature points.
7244554ea0SWill Pazner- `grad1d`:    Matrix of size `(p, q)` expressing derivatives of nodal basis functions at
7344554ea0SWill Pazner               quadrature points.
7444554ea0SWill Pazner- `qref1d`:    Array of length `q` holding the locations of quadrature points on the 1D
7544554ea0SWill Pazner               reference element $[-1, 1]$.
7644554ea0SWill Pazner- `qweight1d`: Array of length `q` holding the quadrature weights on the reference element.
7744554ea0SWill Pazner"""
7844554ea0SWill Paznerfunction create_tensor_h1_basis(
7944554ea0SWill Pazner    c::Ceed,
8044554ea0SWill Pazner    dim,
8144554ea0SWill Pazner    ncomp,
8244554ea0SWill Pazner    p,
8344554ea0SWill Pazner    q,
8480a9ef05SNatalie Beams    interp1d::AbstractArray{CeedScalar},
8580a9ef05SNatalie Beams    grad1d::AbstractArray{CeedScalar},
8680a9ef05SNatalie Beams    qref1d::AbstractArray{CeedScalar},
8780a9ef05SNatalie Beams    qweight1d::AbstractArray{CeedScalar},
8844554ea0SWill Pazner)
8944554ea0SWill Pazner    @assert size(interp1d) == (q, p)
9044554ea0SWill Pazner    @assert size(grad1d) == (q, p)
9144554ea0SWill Pazner    @assert length(qref1d) == q
9244554ea0SWill Pazner    @assert length(qweight1d) == q
9344554ea0SWill Pazner
9444554ea0SWill Pazner    # Convert from Julia matrices (column-major) to row-major format
9544554ea0SWill Pazner    interp1d_rowmajor = collect(interp1d')
9644554ea0SWill Pazner    grad1d_rowmajor = collect(grad1d')
9744554ea0SWill Pazner
9844554ea0SWill Pazner    ref = Ref{C.CeedBasis}()
9944554ea0SWill Pazner    C.CeedBasisCreateTensorH1(
10044554ea0SWill Pazner        c[],
10144554ea0SWill Pazner        dim,
10244554ea0SWill Pazner        ncomp,
10344554ea0SWill Pazner        p,
10444554ea0SWill Pazner        q,
10544554ea0SWill Pazner        interp1d_rowmajor,
10644554ea0SWill Pazner        grad1d_rowmajor,
10744554ea0SWill Pazner        qref1d,
10844554ea0SWill Pazner        qweight1d,
10944554ea0SWill Pazner        ref,
11044554ea0SWill Pazner    )
11144554ea0SWill Pazner    Basis(ref)
11244554ea0SWill Paznerend
11344554ea0SWill Pazner
11444554ea0SWill Pazner@doc raw"""
11544554ea0SWill Pazner    create_h1_basis(c::Ceed, topo::Topology, ncomp, nnodes, nqpts, interp, grad, qref, qweight)
11644554ea0SWill Pazner
117*11b88ddaSSebastian GrimbergCreate a non tensor-product basis for $H^1$ discretizations
11844554ea0SWill Pazner
11944554ea0SWill Pazner# Arguments:
12044554ea0SWill Pazner- `ceed`:    A [`Ceed`](@ref) object where the [`Basis`](@ref) will be created.
12144554ea0SWill Pazner- `topo`:    [`Topology`](@ref) of element, e.g. hypercube, simplex, etc.
12244554ea0SWill Pazner- `ncomp`:   Number of field components (1 for scalar fields).
12344554ea0SWill Pazner- `nnodes`:  Total number of nodes.
12444554ea0SWill Pazner- `nqpts`:   Total number of quadrature points.
12544554ea0SWill Pazner- `interp`:  Matrix of size `(nqpts, nnodes)` expressing the values of nodal basis functions
12644554ea0SWill Pazner             at quadrature points.
12744554ea0SWill Pazner- `grad`:    Array of size `(dim, nqpts, nnodes)` expressing derivatives of nodal basis
12844554ea0SWill Pazner             functions at quadrature points.
12944554ea0SWill Pazner- `qref`:    Array of length `nqpts` holding the locations of quadrature points on the
13044554ea0SWill Pazner             reference element $[-1, 1]$.
13144554ea0SWill Pazner- `qweight`: Array of length `nqpts` holding the quadrature weights on the reference
13244554ea0SWill Pazner             element.
13344554ea0SWill Pazner"""
13444554ea0SWill Paznerfunction create_h1_basis(
13544554ea0SWill Pazner    c::Ceed,
13644554ea0SWill Pazner    topo::Topology,
13744554ea0SWill Pazner    ncomp,
13844554ea0SWill Pazner    nnodes,
13944554ea0SWill Pazner    nqpts,
14080a9ef05SNatalie Beams    interp::AbstractArray{CeedScalar},
14180a9ef05SNatalie Beams    grad::AbstractArray{CeedScalar},
14280a9ef05SNatalie Beams    qref::AbstractArray{CeedScalar},
14380a9ef05SNatalie Beams    qweight::AbstractArray{CeedScalar},
14444554ea0SWill Pazner)
14544554ea0SWill Pazner    dim = getdimension(topo)
14644554ea0SWill Pazner    @assert size(interp) == (nqpts, nnodes)
14744554ea0SWill Pazner    @assert size(grad) == (dim, nqpts, nnodes)
14844554ea0SWill Pazner    @assert length(qref) == nqpts
14944554ea0SWill Pazner    @assert length(qweight) == nqpts
15044554ea0SWill Pazner
15144554ea0SWill Pazner    # Convert from Julia matrices and tensors (column-major) to row-major format
15244554ea0SWill Pazner    interp_rowmajor = collect(interp')
15344554ea0SWill Pazner    grad_rowmajor = permutedims(grad, [3, 2, 1])
15444554ea0SWill Pazner
15544554ea0SWill Pazner    ref = Ref{C.CeedBasis}()
15644554ea0SWill Pazner    C.CeedBasisCreateH1(
15744554ea0SWill Pazner        c[],
15844554ea0SWill Pazner        topo,
15944554ea0SWill Pazner        ncomp,
16044554ea0SWill Pazner        nnodes,
16144554ea0SWill Pazner        nqpts,
16244554ea0SWill Pazner        interp_rowmajor,
16344554ea0SWill Pazner        grad_rowmajor,
16444554ea0SWill Pazner        qref,
16544554ea0SWill Pazner        qweight,
16644554ea0SWill Pazner        ref,
16744554ea0SWill Pazner    )
16844554ea0SWill Pazner    Basis(ref)
16944554ea0SWill Paznerend
17044554ea0SWill Pazner
171*11b88ddaSSebastian Grimberg@doc raw"""
172*11b88ddaSSebastian Grimberg    create_hdiv_basis(c::Ceed, topo::Topology, ncomp, nnodes, nqpts, interp, div, qref, qweight)
173*11b88ddaSSebastian Grimberg
174*11b88ddaSSebastian GrimbergCreate a non tensor-product basis for H(div) discretizations
175*11b88ddaSSebastian Grimberg
176*11b88ddaSSebastian Grimberg# Arguments:
177*11b88ddaSSebastian Grimberg- `ceed`:    A [`Ceed`](@ref) object where the [`Basis`](@ref) will be created.
178*11b88ddaSSebastian Grimberg- `topo`:    [`Topology`](@ref) of element, e.g. hypercube, simplex, etc.
179*11b88ddaSSebastian Grimberg- `ncomp`:   Number of field components (1 for scalar fields).
180*11b88ddaSSebastian Grimberg- `nnodes`:  Total number of nodes.
181*11b88ddaSSebastian Grimberg- `nqpts`:   Total number of quadrature points.
182*11b88ddaSSebastian Grimberg- `interp`:  Matrix of size `(dim, nqpts, nnodes)` expressing the values of basis functions
183*11b88ddaSSebastian Grimberg             at quadrature points.
184*11b88ddaSSebastian Grimberg- `div`:     Array of size `(nqpts, nnodes)` expressing divergence of basis functions at
185*11b88ddaSSebastian Grimberg             quadrature points.
186*11b88ddaSSebastian Grimberg- `qref`:    Array of length `nqpts` holding the locations of quadrature points on the
187*11b88ddaSSebastian Grimberg             reference element $[-1, 1]$.
188*11b88ddaSSebastian Grimberg- `qweight`: Array of length `nqpts` holding the quadrature weights on the reference
189*11b88ddaSSebastian Grimberg             element.
190*11b88ddaSSebastian Grimberg"""
191*11b88ddaSSebastian Grimbergfunction create_hdiv_basis(
192*11b88ddaSSebastian Grimberg    c::Ceed,
193*11b88ddaSSebastian Grimberg    topo::Topology,
194*11b88ddaSSebastian Grimberg    ncomp,
195*11b88ddaSSebastian Grimberg    nnodes,
196*11b88ddaSSebastian Grimberg    nqpts,
197*11b88ddaSSebastian Grimberg    interp::AbstractArray{CeedScalar},
198*11b88ddaSSebastian Grimberg    div::AbstractArray{CeedScalar},
199*11b88ddaSSebastian Grimberg    qref::AbstractArray{CeedScalar},
200*11b88ddaSSebastian Grimberg    qweight::AbstractArray{CeedScalar},
201*11b88ddaSSebastian Grimberg)
202*11b88ddaSSebastian Grimberg    dim = getdimension(topo)
203*11b88ddaSSebastian Grimberg    @assert size(interp) == (dim, nqpts, nnodes)
204*11b88ddaSSebastian Grimberg    @assert size(div) == (nqpts, nnodes)
205*11b88ddaSSebastian Grimberg    @assert length(qref) == nqpts
206*11b88ddaSSebastian Grimberg    @assert length(qweight) == nqpts
207*11b88ddaSSebastian Grimberg
208*11b88ddaSSebastian Grimberg    # Convert from Julia matrices and tensors (column-major) to row-major format
209*11b88ddaSSebastian Grimberg    interp_rowmajor = permutedims(interp, [3, 2, 1])
210*11b88ddaSSebastian Grimberg    div_rowmajor = collect(div')
211*11b88ddaSSebastian Grimberg
212*11b88ddaSSebastian Grimberg    ref = Ref{C.CeedBasis}()
213*11b88ddaSSebastian Grimberg    C.CeedBasisCreateHdiv(
214*11b88ddaSSebastian Grimberg        c[],
215*11b88ddaSSebastian Grimberg        topo,
216*11b88ddaSSebastian Grimberg        ncomp,
217*11b88ddaSSebastian Grimberg        nnodes,
218*11b88ddaSSebastian Grimberg        nqpts,
219*11b88ddaSSebastian Grimberg        interp_rowmajor,
220*11b88ddaSSebastian Grimberg        div_rowmajor,
221*11b88ddaSSebastian Grimberg        qref,
222*11b88ddaSSebastian Grimberg        qweight,
223*11b88ddaSSebastian Grimberg        ref,
224*11b88ddaSSebastian Grimberg    )
225*11b88ddaSSebastian Grimberg    Basis(ref)
226*11b88ddaSSebastian Grimbergend
227*11b88ddaSSebastian Grimberg
228*11b88ddaSSebastian Grimberg@doc raw"""
229*11b88ddaSSebastian Grimberg    create_hcurl_basis(c::Ceed, topo::Topology, ncomp, nnodes, nqpts, interp, curl, qref, qweight)
230*11b88ddaSSebastian Grimberg
231*11b88ddaSSebastian GrimbergCreate a non tensor-product basis for H(curl) discretizations
232*11b88ddaSSebastian Grimberg
233*11b88ddaSSebastian Grimberg# Arguments:
234*11b88ddaSSebastian Grimberg- `ceed`:    A [`Ceed`](@ref) object where the [`Basis`](@ref) will be created.
235*11b88ddaSSebastian Grimberg- `topo`:    [`Topology`](@ref) of element, e.g. hypercube, simplex, etc.
236*11b88ddaSSebastian Grimberg- `ncomp`:   Number of field components (1 for scalar fields).
237*11b88ddaSSebastian Grimberg- `nnodes`:  Total number of nodes.
238*11b88ddaSSebastian Grimberg- `nqpts`:   Total number of quadrature points.
239*11b88ddaSSebastian Grimberg- `interp`:  Matrix of size `(dim, nqpts, nnodes)` expressing the values of basis functions
240*11b88ddaSSebastian Grimberg             at quadrature points.
241*11b88ddaSSebastian Grimberg- `curl`:    Matrix of size `(curlcomp, nqpts, nnodes)`, `curlcomp = 1 if dim < 3 else dim`)
242*11b88ddaSSebastian Grimberg             matrix expressing curl of basis functions at quadrature points.
243*11b88ddaSSebastian Grimberg- `qref`:    Array of length `nqpts` holding the locations of quadrature points on the
244*11b88ddaSSebastian Grimberg             reference element $[-1, 1]$.
245*11b88ddaSSebastian Grimberg- `qweight`: Array of length `nqpts` holding the quadrature weights on the reference
246*11b88ddaSSebastian Grimberg             element.
247*11b88ddaSSebastian Grimberg"""
248*11b88ddaSSebastian Grimbergfunction create_hcurl_basis(
249*11b88ddaSSebastian Grimberg    c::Ceed,
250*11b88ddaSSebastian Grimberg    topo::Topology,
251*11b88ddaSSebastian Grimberg    ncomp,
252*11b88ddaSSebastian Grimberg    nnodes,
253*11b88ddaSSebastian Grimberg    nqpts,
254*11b88ddaSSebastian Grimberg    interp::AbstractArray{CeedScalar},
255*11b88ddaSSebastian Grimberg    curl::AbstractArray{CeedScalar},
256*11b88ddaSSebastian Grimberg    qref::AbstractArray{CeedScalar},
257*11b88ddaSSebastian Grimberg    qweight::AbstractArray{CeedScalar},
258*11b88ddaSSebastian Grimberg)
259*11b88ddaSSebastian Grimberg    dim = getdimension(topo)
260*11b88ddaSSebastian Grimberg    curlcomp = dim < 3 ? 1 : dim
261*11b88ddaSSebastian Grimberg    @assert size(interp) == (dim, nqpts, nnodes)
262*11b88ddaSSebastian Grimberg    @assert size(curl) == (curlcomp, nqpts, nnodes)
263*11b88ddaSSebastian Grimberg    @assert length(qref) == nqpts
264*11b88ddaSSebastian Grimberg    @assert length(qweight) == nqpts
265*11b88ddaSSebastian Grimberg
266*11b88ddaSSebastian Grimberg    # Convert from Julia matrices and tensors (column-major) to row-major format
267*11b88ddaSSebastian Grimberg    interp_rowmajor = permutedims(interp, [3, 2, 1])
268*11b88ddaSSebastian Grimberg    curl_rowmajor = permutedims(curl, [3, 2, 1])
269*11b88ddaSSebastian Grimberg
270*11b88ddaSSebastian Grimberg    ref = Ref{C.CeedBasis}()
271*11b88ddaSSebastian Grimberg    C.CeedBasisCreateHcurl(
272*11b88ddaSSebastian Grimberg        c[],
273*11b88ddaSSebastian Grimberg        topo,
274*11b88ddaSSebastian Grimberg        ncomp,
275*11b88ddaSSebastian Grimberg        nnodes,
276*11b88ddaSSebastian Grimberg        nqpts,
277*11b88ddaSSebastian Grimberg        interp_rowmajor,
278*11b88ddaSSebastian Grimberg        curl_rowmajor,
279*11b88ddaSSebastian Grimberg        qref,
280*11b88ddaSSebastian Grimberg        qweight,
281*11b88ddaSSebastian Grimberg        ref,
282*11b88ddaSSebastian Grimberg    )
283*11b88ddaSSebastian Grimberg    Basis(ref)
284*11b88ddaSSebastian Grimbergend
285*11b88ddaSSebastian Grimberg
28644554ea0SWill Pazner"""
28744554ea0SWill Pazner    apply!(b::Basis, nelem, tmode::TransposeMode, emode::EvalMode, u::AbstractCeedVector, v::AbstractCeedVector)
28844554ea0SWill Pazner
28944554ea0SWill PaznerApply basis evaluation from nodes to quadrature points or vice versa, storing the result in
29044554ea0SWill Paznerthe [`CeedVector`](@ref) `v`.
29144554ea0SWill Pazner
29244554ea0SWill Pazner`nelem` specifies the number of elements to apply the basis evaluation to; the backend will
29344554ea0SWill Paznerspecify the ordering in CeedElemRestrictionCreateBlocked()
29444554ea0SWill Pazner
29544554ea0SWill PaznerSet `tmode` to `CEED_NOTRANSPOSE` to evaluate from nodes to quadrature or to
29644554ea0SWill Pazner`CEED_TRANSPOSE` to apply the transpose, mapping from quadrature points to nodes.
29744554ea0SWill Pazner
29844554ea0SWill PaznerSet the [`EvalMode`](@ref) `emode` to:
29944554ea0SWill Pazner- `CEED_EVAL_NONE` to use values directly,
30044554ea0SWill Pazner- `CEED_EVAL_INTERP` to use interpolated values,
30144554ea0SWill Pazner- `CEED_EVAL_GRAD` to use gradients,
30244554ea0SWill Pazner- `CEED_EVAL_WEIGHT` to use quadrature weights.
30344554ea0SWill Pazner"""
30444554ea0SWill Paznerfunction apply!(
30544554ea0SWill Pazner    b::Basis,
30644554ea0SWill Pazner    nelem,
30744554ea0SWill Pazner    tmode::TransposeMode,
30844554ea0SWill Pazner    emode::EvalMode,
30944554ea0SWill Pazner    u::AbstractCeedVector,
31044554ea0SWill Pazner    v::AbstractCeedVector,
31144554ea0SWill Pazner)
31244554ea0SWill Pazner    C.CeedBasisApply(b[], nelem, tmode, emode, u[], v[])
31344554ea0SWill Paznerend
31444554ea0SWill Pazner
31544554ea0SWill Pazner"""
31644554ea0SWill Pazner    apply(b::Basis, u::AbstractVector; nelem=1, tmode=NOTRANSPOSE, emode=EVAL_INTERP)
31744554ea0SWill Pazner
31844554ea0SWill PaznerPerforms the same function as the above-defined [`apply!`](@ref apply!(b::Basis, nelem,
31944554ea0SWill Paznertmode::TransposeMode, emode::EvalMode, u::AbstractCeedVector, v::AbstractCeedVector)), but
32044554ea0SWill Paznerautomatically convert from Julia arrays to [`CeedVector`](@ref) for convenience.
32144554ea0SWill Pazner
32244554ea0SWill PaznerThe result will be returned in a newly allocated array of the correct size.
32344554ea0SWill Pazner"""
32444554ea0SWill Paznerfunction apply(b::Basis, u::AbstractVector; nelem=1, tmode=NOTRANSPOSE, emode=EVAL_INTERP)
32544554ea0SWill Pazner    ceed_ref = Ref{C.Ceed}()
32644554ea0SWill Pazner    ccall((:CeedBasisGetCeed, C.libceed), Cint, (C.CeedBasis, Ptr{C.Ceed}), b[], ceed_ref)
32744554ea0SWill Pazner    c = Ceed(ceed_ref)
32844554ea0SWill Pazner
32944554ea0SWill Pazner    u_vec = CeedVector(c, u)
33044554ea0SWill Pazner
33144554ea0SWill Pazner    len_v = (tmode == TRANSPOSE) ? getnumnodes(b) : getnumqpts(b)
33244554ea0SWill Pazner    if emode == EVAL_GRAD
33344554ea0SWill Pazner        len_v *= getdimension(b)
33444554ea0SWill Pazner    end
33544554ea0SWill Pazner
33644554ea0SWill Pazner    v_vec = CeedVector(c, len_v)
33744554ea0SWill Pazner
33844554ea0SWill Pazner    apply!(b, nelem, tmode, emode, u_vec, v_vec)
33944554ea0SWill Pazner    Vector(v_vec)
34044554ea0SWill Paznerend
34144554ea0SWill Pazner
34244554ea0SWill Pazner"""
34344554ea0SWill Pazner    getdimension(b::Basis)
34444554ea0SWill Pazner
34544554ea0SWill PaznerReturn the spatial dimension of the given [`Basis`](@ref).
34644554ea0SWill Pazner"""
34744554ea0SWill Paznerfunction getdimension(b::Basis)
34844554ea0SWill Pazner    dim = Ref{CeedInt}()
34944554ea0SWill Pazner    C.CeedBasisGetDimension(b[], dim)
35044554ea0SWill Pazner    dim[]
35144554ea0SWill Paznerend
35244554ea0SWill Pazner
35344554ea0SWill Pazner"""
35444554ea0SWill Pazner    getdimension(t::Topology)
35544554ea0SWill Pazner
35644554ea0SWill PaznerReturn the spatial dimension of the given [`Topology`](@ref).
35744554ea0SWill Pazner"""
35844554ea0SWill Paznerfunction getdimension(t::Topology)
35944554ea0SWill Pazner    return Int(t) >> 16
36044554ea0SWill Paznerend
36144554ea0SWill Pazner
36244554ea0SWill Pazner"""
36344554ea0SWill Pazner    gettopology(b::Basis)
36444554ea0SWill Pazner
36544554ea0SWill PaznerReturn the [`Topology`](@ref) of the given [`Basis`](@ref).
36644554ea0SWill Pazner"""
36744554ea0SWill Paznerfunction gettopology(b::Basis)
36844554ea0SWill Pazner    topo = Ref{Topology}()
36944554ea0SWill Pazner    C.CeedBasisGetTopology(b[], topo)
37044554ea0SWill Pazner    topo[]
37144554ea0SWill Paznerend
37244554ea0SWill Pazner
37344554ea0SWill Pazner"""
37444554ea0SWill Pazner    getnumcomponents(b::Basis)
37544554ea0SWill Pazner
37644554ea0SWill PaznerReturn the number of components of the given [`Basis`](@ref).
37744554ea0SWill Pazner"""
37844554ea0SWill Paznerfunction getnumcomponents(b::Basis)
37944554ea0SWill Pazner    ncomp = Ref{CeedInt}()
38044554ea0SWill Pazner    C.CeedBasisGetNumComponents(b[], ncomp)
38144554ea0SWill Pazner    ncomp[]
38244554ea0SWill Paznerend
38344554ea0SWill Pazner
38444554ea0SWill Pazner"""
38544554ea0SWill Pazner    getnumnodes(b::Basis)
38644554ea0SWill Pazner
38744554ea0SWill PaznerReturn the number of nodes of the given [`Basis`](@ref).
38844554ea0SWill Pazner"""
38944554ea0SWill Paznerfunction getnumnodes(b::Basis)
39044554ea0SWill Pazner    nnodes = Ref{CeedInt}()
39144554ea0SWill Pazner    C.CeedBasisGetNumNodes(b[], nnodes)
39244554ea0SWill Pazner    nnodes[]
39344554ea0SWill Paznerend
39444554ea0SWill Pazner
39544554ea0SWill Pazner"""
39644554ea0SWill Pazner    getnumnodes1d(b::Basis)
39744554ea0SWill Pazner
39844554ea0SWill Pazner    Return the number of 1D nodes of the given (tensor-product) [`Basis`](@ref).
39944554ea0SWill Pazner"""
40044554ea0SWill Paznerfunction getnumnodes1d(b::Basis)
40144554ea0SWill Pazner    nnodes1d = Ref{CeedInt}()
40244554ea0SWill Pazner    C.CeedBasisGetNumNodes1D(b[], nnodes1d)
40344554ea0SWill Pazner    nnodes1d[]
40444554ea0SWill Paznerend
40544554ea0SWill Pazner
40644554ea0SWill Pazner"""
40744554ea0SWill Pazner    getnumqpts(b::Basis)
40844554ea0SWill Pazner
40944554ea0SWill PaznerReturn the number of quadrature points of the given [`Basis`](@ref).
41044554ea0SWill Pazner"""
41144554ea0SWill Paznerfunction getnumqpts(b::Basis)
41244554ea0SWill Pazner    nqpts = Ref{CeedInt}()
41344554ea0SWill Pazner    C.CeedBasisGetNumQuadraturePoints(b[], nqpts)
41444554ea0SWill Pazner    nqpts[]
41544554ea0SWill Paznerend
41644554ea0SWill Pazner
41744554ea0SWill Pazner"""
41844554ea0SWill Pazner    getnumqpts1d(b::Basis)
41944554ea0SWill Pazner
42044554ea0SWill PaznerReturn the number of 1D quadrature points of the given (tensor-product) [`Basis`](@ref).
42144554ea0SWill Pazner"""
42244554ea0SWill Paznerfunction getnumqpts1d(b::Basis)
42344554ea0SWill Pazner    nqpts1d = Ref{CeedInt}()
42444554ea0SWill Pazner    C.CeedBasisGetNumQuadraturePoints1D(b[], nqpts1d)
42544554ea0SWill Pazner    nqpts1d[]
42644554ea0SWill Paznerend
42744554ea0SWill Pazner
42844554ea0SWill Pazner"""
42944554ea0SWill Pazner    getqref(b::Basis)
43044554ea0SWill Pazner
43144554ea0SWill PaznerGet the reference coordinates of quadrature points (in `dim` dimensions) of the given
43244554ea0SWill Pazner[`Basis`](@ref).
43344554ea0SWill Pazner"""
43444554ea0SWill Paznerfunction getqref(b::Basis)
435a9e65696SJeremy L Thompson    istensor = Ref{Bool}()
436a9e65696SJeremy L Thompson    C.CeedBasisIsTensor(b[], istensor)
43744554ea0SWill Pazner    ref = Ref{Ptr{CeedScalar}}()
43844554ea0SWill Pazner    C.CeedBasisGetQRef(b[], ref)
439a9e65696SJeremy L Thompson    copy(
440a9e65696SJeremy L Thompson        unsafe_wrap(
441a9e65696SJeremy L Thompson            Array,
442a9e65696SJeremy L Thompson            ref[],
443a9e65696SJeremy L Thompson            istensor[] ? getnumqpts1d(b) : (getnumqpts(b)*getdimension(b)),
444a9e65696SJeremy L Thompson        ),
445a9e65696SJeremy L Thompson    )
44644554ea0SWill Paznerend
44744554ea0SWill Pazner
44844554ea0SWill Pazner"""
44944554ea0SWill Pazner    getqref(b::Basis)
45044554ea0SWill Pazner
45144554ea0SWill PaznerGet the quadrature weights of quadrature points (in `dim` dimensions) of the given
45244554ea0SWill Pazner[`Basis`](@ref).
45344554ea0SWill Pazner"""
45444554ea0SWill Paznerfunction getqweights(b::Basis)
455a9e65696SJeremy L Thompson    istensor = Ref{Bool}()
456a9e65696SJeremy L Thompson    C.CeedBasisIsTensor(b[], istensor)
45744554ea0SWill Pazner    ref = Ref{Ptr{CeedScalar}}()
45844554ea0SWill Pazner    C.CeedBasisGetQWeights(b[], ref)
459a9e65696SJeremy L Thompson    copy(unsafe_wrap(Array, ref[], istensor[] ? getnumqpts1d(b) : getnumqpts(b)))
46044554ea0SWill Paznerend
46144554ea0SWill Pazner
462*11b88ddaSSebastian Grimberg@doc raw"""
46344554ea0SWill Pazner    getinterp(b::Basis)
46444554ea0SWill Pazner
46544554ea0SWill PaznerGet the interpolation matrix of the given [`Basis`](@ref). Returns a matrix of size
466*11b88ddaSSebastian Grimberg`(getnumqpts(b), getnumnodes(b))` for a given $H^1$ basis or `(getdimension(b),
467*11b88ddaSSebastian Grimberggetnumqpts(b), getnumnodes(b))` for a given vector $H(div)$ or $H(curl)$ basis.
46844554ea0SWill Pazner"""
46944554ea0SWill Paznerfunction getinterp(b::Basis)
47044554ea0SWill Pazner    ref = Ref{Ptr{CeedScalar}}()
47144554ea0SWill Pazner    C.CeedBasisGetInterp(b[], ref)
47244554ea0SWill Pazner    q = getnumqpts(b)
47344554ea0SWill Pazner    p = getnumnodes(b)
474*11b88ddaSSebastian Grimberg    qcomp = Ref{CeedInt}()
475*11b88ddaSSebastian Grimberg    C.CeedBasisGetNumQuadratureComponents(b[], C.CEED_EVAL_INTERP, qcomp)
476*11b88ddaSSebastian Grimberg    if qcomp[] == 1
47744554ea0SWill Pazner        collect(unsafe_wrap(Array, ref[], (p, q))')
478*11b88ddaSSebastian Grimberg    else
479*11b88ddaSSebastian Grimberg        permutedims(unsafe_wrap(Array, ref[], (p, q, qcomp[])), [3, 2, 1])
480*11b88ddaSSebastian Grimberg    end
48144554ea0SWill Paznerend
48244554ea0SWill Pazner
48344554ea0SWill Pazner"""
48444554ea0SWill Pazner    getinterp1d(b::Basis)
48544554ea0SWill Pazner
48644554ea0SWill PaznerGet the 1D interpolation matrix of the given [`Basis`](@ref). `b` must be a tensor-product
48744554ea0SWill Paznerbasis, otherwise this function will fail. Returns a matrix of size `(getnumqpts1d(b),
48844554ea0SWill Paznergetnumnodes1d(b))`.
48944554ea0SWill Pazner"""
49044554ea0SWill Paznerfunction getinterp1d(b::Basis)
49144554ea0SWill Pazner    ref = Ref{Ptr{CeedScalar}}()
49244554ea0SWill Pazner    C.CeedBasisGetInterp1D(b[], ref)
49344554ea0SWill Pazner    q = getnumqpts1d(b)
49444554ea0SWill Pazner    p = getnumnodes1d(b)
49544554ea0SWill Pazner    collect(unsafe_wrap(Array, ref[], (p, q))')
49644554ea0SWill Paznerend
49744554ea0SWill Pazner
49844554ea0SWill Pazner"""
4997574393bSJeremy L Thompson    getgrad(b::Basis)
50044554ea0SWill Pazner
50144554ea0SWill PaznerGet the gradient matrix of the given [`Basis`](@ref). Returns a tensor of size
50244554ea0SWill Pazner`(getdimension(b), getnumqpts(b), getnumnodes(b))`.
50344554ea0SWill Pazner"""
50444554ea0SWill Paznerfunction getgrad(b::Basis)
50544554ea0SWill Pazner    ref = Ref{Ptr{CeedScalar}}()
50644554ea0SWill Pazner    C.CeedBasisGetGrad(b[], ref)
50744554ea0SWill Pazner    dim = getdimension(b)
50844554ea0SWill Pazner    q = getnumqpts(b)
50944554ea0SWill Pazner    p = getnumnodes(b)
51044554ea0SWill Pazner    permutedims(unsafe_wrap(Array, ref[], (p, q, dim)), [3, 2, 1])
51144554ea0SWill Paznerend
51244554ea0SWill Pazner
51344554ea0SWill Pazner"""
51444554ea0SWill Pazner    getgrad1d(b::Basis)
51544554ea0SWill Pazner
51644554ea0SWill PaznerGet the 1D derivative matrix of the given [`Basis`](@ref). Returns a matrix of size
51744554ea0SWill Pazner`(getnumqpts(b), getnumnodes(b))`.
51844554ea0SWill Pazner"""
51944554ea0SWill Paznerfunction getgrad1d(b::Basis)
52044554ea0SWill Pazner    ref = Ref{Ptr{CeedScalar}}()
52144554ea0SWill Pazner    C.CeedBasisGetGrad1D(b[], ref)
52244554ea0SWill Pazner    q = getnumqpts1d(b)
52344554ea0SWill Pazner    p = getnumnodes1d(b)
52444554ea0SWill Pazner    collect(unsafe_wrap(Array, ref[], (p, q))')
52544554ea0SWill Paznerend
526*11b88ddaSSebastian Grimberg
527*11b88ddaSSebastian Grimberg"""
528*11b88ddaSSebastian Grimberg    getdiv(b::Basis)
529*11b88ddaSSebastian Grimberg
530*11b88ddaSSebastian GrimbergGet the divergence matrix of the given [`Basis`](@ref). Returns a tensor of size
531*11b88ddaSSebastian Grimberg`(getnumqpts(b), getnumnodes(b))`.
532*11b88ddaSSebastian Grimberg"""
533*11b88ddaSSebastian Grimbergfunction getdiv(b::Basis)
534*11b88ddaSSebastian Grimberg    ref = Ref{Ptr{CeedScalar}}()
535*11b88ddaSSebastian Grimberg    C.CeedBasisGetDiv(b[], ref)
536*11b88ddaSSebastian Grimberg    q = getnumqpts(b)
537*11b88ddaSSebastian Grimberg    p = getnumnodes(b)
538*11b88ddaSSebastian Grimberg    collect(unsafe_wrap(Array, ref[], (p, q))')
539*11b88ddaSSebastian Grimbergend
540*11b88ddaSSebastian Grimberg
541*11b88ddaSSebastian Grimberg"""
542*11b88ddaSSebastian Grimberg    getcurl(b::Basis)
543*11b88ddaSSebastian Grimberg
544*11b88ddaSSebastian GrimbergGet the curl matrix of the given [`Basis`](@ref). Returns a tensor of size
545*11b88ddaSSebastian Grimberg`(curlcomp, getnumqpts(b), getnumnodes(b))`, `curlcomp = 1 if getdimension(b) < 3 else
546*11b88ddaSSebastian Grimberggetdimension(b)`.
547*11b88ddaSSebastian Grimberg"""
548*11b88ddaSSebastian Grimbergfunction getcurl(b::Basis)
549*11b88ddaSSebastian Grimberg    ref = Ref{Ptr{CeedScalar}}()
550*11b88ddaSSebastian Grimberg    C.CeedBasisGetCurl(b[], ref)
551*11b88ddaSSebastian Grimberg    q = getnumqpts(b)
552*11b88ddaSSebastian Grimberg    p = getnumnodes(b)
553*11b88ddaSSebastian Grimberg    qcomp = Ref{CeedInt}()
554*11b88ddaSSebastian Grimberg    C.CeedBasisGetNumQuadratureComponents(b[], C.CEED_EVAL_CURL, qcomp)
555*11b88ddaSSebastian Grimberg    permutedims(unsafe_wrap(Array, ref[], (p, q, qcomp[])), [3, 2, 1])
556*11b88ddaSSebastian Grimbergend
557