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) 2044554ea0SWill Pazner""" 2144554ea0SWill Paznermutable struct Basis <: AbstractBasis 2244554ea0SWill Pazner ref::RefValue{C.CeedBasis} 2344554ea0SWill Pazner function Basis(ref) 2444554ea0SWill Pazner obj = new(ref) 2544554ea0SWill Pazner finalizer(obj) do x 2644554ea0SWill Pazner # ccall(:jl_safe_printf, Cvoid, (Cstring, Cstring), "Finalizing %s.\n", repr(x)) 2744554ea0SWill Pazner destroy(x) 2844554ea0SWill Pazner end 2944554ea0SWill Pazner return obj 3044554ea0SWill Pazner end 3144554ea0SWill Paznerend 3244554ea0SWill Paznerdestroy(b::Basis) = C.CeedBasisDestroy(b.ref) # COV_EXCL_LINE 3344554ea0SWill PaznerBase.getindex(b::Basis) = b.ref[] 3444554ea0SWill PaznerBase.show(io::IO, ::MIME"text/plain", b::Basis) = ceed_show(io, b, C.CeedBasisView) 3544554ea0SWill Pazner 3644554ea0SWill Pazner@doc raw""" 3744554ea0SWill Pazner create_tensor_h1_lagrange_basis(ceed, dim, ncomp, p, q, qmode) 3844554ea0SWill Pazner 3944554ea0SWill PaznerCreate a tensor-product Lagrange basis. 4044554ea0SWill Pazner 4144554ea0SWill Pazner# Arguments: 4244554ea0SWill Pazner- `ceed`: A [`Ceed`](@ref) object where the [`Basis`](@ref) will be created. 4344554ea0SWill Pazner- `dim`: Topological dimension of element. 4444554ea0SWill Pazner- `ncomp`: Number of field components (1 for scalar fields). 4544554ea0SWill Pazner- `p`: Number of Gauss-Lobatto nodes in one dimension. The polynomial degree of the 4644554ea0SWill Pazner resulting $Q_k$ element is $k=p-1$. 4744554ea0SWill Pazner- `q`: Number of quadrature points in one dimension. 4844554ea0SWill Pazner- `qmode`: Distribution of the $q$ quadrature points (affects order of accuracy for the 4944554ea0SWill Pazner quadrature). 5044554ea0SWill Pazner""" 5144554ea0SWill Paznerfunction create_tensor_h1_lagrange_basis(c::Ceed, dim, ncomp, p, q, quad_mode::QuadMode) 5244554ea0SWill Pazner ref = Ref{C.CeedBasis}() 5344554ea0SWill Pazner C.CeedBasisCreateTensorH1Lagrange(c[], dim, ncomp, p, q, quad_mode, ref) 5444554ea0SWill Pazner Basis(ref) 5544554ea0SWill Paznerend 5644554ea0SWill Pazner 5744554ea0SWill Pazner@doc raw""" 5844554ea0SWill Pazner create_tensor_h1_basis(c::Ceed, dim, ncomp, p, q, interp1d, grad1d, qref1d, qweight1d) 5944554ea0SWill Pazner 6044554ea0SWill PaznerCreate a tensor-product basis for $H^1$ discretizations. 6144554ea0SWill Pazner 6244554ea0SWill Pazner# Arguments: 6344554ea0SWill Pazner- `ceed`: A [`Ceed`](@ref) object where the [`Basis`](@ref) will be created. 6444554ea0SWill Pazner- `dim`: Topological dimension. 6544554ea0SWill Pazner- `ncomp`: Number of field components (1 for scalar fields). 6644554ea0SWill Pazner- `p`: Number of nodes in one dimension. 6744554ea0SWill Pazner- `q`: Number of quadrature points in one dimension 6844554ea0SWill Pazner- `interp1d`: Matrix of size `(q, p)` expressing the values of nodal basis functions at 6944554ea0SWill Pazner quadrature points. 7044554ea0SWill Pazner- `grad1d`: Matrix of size `(p, q)` expressing derivatives of nodal basis functions at 7144554ea0SWill Pazner quadrature points. 7244554ea0SWill Pazner- `qref1d`: Array of length `q` holding the locations of quadrature points on the 1D 7344554ea0SWill Pazner reference element $[-1, 1]$. 7444554ea0SWill Pazner- `qweight1d`: Array of length `q` holding the quadrature weights on the reference element. 7544554ea0SWill Pazner""" 7644554ea0SWill Paznerfunction create_tensor_h1_basis( 7744554ea0SWill Pazner c::Ceed, 7844554ea0SWill Pazner dim, 7944554ea0SWill Pazner ncomp, 8044554ea0SWill Pazner p, 8144554ea0SWill Pazner q, 8280a9ef05SNatalie Beams interp1d::AbstractArray{CeedScalar}, 8380a9ef05SNatalie Beams grad1d::AbstractArray{CeedScalar}, 8480a9ef05SNatalie Beams qref1d::AbstractArray{CeedScalar}, 8580a9ef05SNatalie Beams qweight1d::AbstractArray{CeedScalar}, 8644554ea0SWill Pazner) 8744554ea0SWill Pazner @assert size(interp1d) == (q, p) 8844554ea0SWill Pazner @assert size(grad1d) == (q, p) 8944554ea0SWill Pazner @assert length(qref1d) == q 9044554ea0SWill Pazner @assert length(qweight1d) == q 9144554ea0SWill Pazner 9244554ea0SWill Pazner # Convert from Julia matrices (column-major) to row-major format 9344554ea0SWill Pazner interp1d_rowmajor = collect(interp1d') 9444554ea0SWill Pazner grad1d_rowmajor = collect(grad1d') 9544554ea0SWill Pazner 9644554ea0SWill Pazner ref = Ref{C.CeedBasis}() 9744554ea0SWill Pazner C.CeedBasisCreateTensorH1( 9844554ea0SWill Pazner c[], 9944554ea0SWill Pazner dim, 10044554ea0SWill Pazner ncomp, 10144554ea0SWill Pazner p, 10244554ea0SWill Pazner q, 10344554ea0SWill Pazner interp1d_rowmajor, 10444554ea0SWill Pazner grad1d_rowmajor, 10544554ea0SWill Pazner qref1d, 10644554ea0SWill Pazner qweight1d, 10744554ea0SWill Pazner ref, 10844554ea0SWill Pazner ) 10944554ea0SWill Pazner Basis(ref) 11044554ea0SWill Paznerend 11144554ea0SWill Pazner 11244554ea0SWill Pazner@doc raw""" 11344554ea0SWill Pazner create_h1_basis(c::Ceed, topo::Topology, ncomp, nnodes, nqpts, interp, grad, qref, qweight) 11444554ea0SWill Pazner 11544554ea0SWill PaznerCreate a non tensor-product basis for H^1 discretizations 11644554ea0SWill Pazner 11744554ea0SWill Pazner# Arguments: 11844554ea0SWill Pazner- `ceed`: A [`Ceed`](@ref) object where the [`Basis`](@ref) will be created. 11944554ea0SWill Pazner- `topo`: [`Topology`](@ref) of element, e.g. hypercube, simplex, etc. 12044554ea0SWill Pazner- `ncomp`: Number of field components (1 for scalar fields). 12144554ea0SWill Pazner- `nnodes`: Total number of nodes. 12244554ea0SWill Pazner- `nqpts`: Total number of quadrature points. 12344554ea0SWill Pazner- `interp`: Matrix of size `(nqpts, nnodes)` expressing the values of nodal basis functions 12444554ea0SWill Pazner at quadrature points. 12544554ea0SWill Pazner- `grad`: Array of size `(dim, nqpts, nnodes)` expressing derivatives of nodal basis 12644554ea0SWill Pazner functions at quadrature points. 12744554ea0SWill Pazner- `qref`: Array of length `nqpts` holding the locations of quadrature points on the 12844554ea0SWill Pazner reference element $[-1, 1]$. 12944554ea0SWill Pazner- `qweight`: Array of length `nqpts` holding the quadrature weights on the reference 13044554ea0SWill Pazner element. 13144554ea0SWill Pazner""" 13244554ea0SWill Paznerfunction create_h1_basis( 13344554ea0SWill Pazner c::Ceed, 13444554ea0SWill Pazner topo::Topology, 13544554ea0SWill Pazner ncomp, 13644554ea0SWill Pazner nnodes, 13744554ea0SWill Pazner nqpts, 13880a9ef05SNatalie Beams interp::AbstractArray{CeedScalar}, 13980a9ef05SNatalie Beams grad::AbstractArray{CeedScalar}, 14080a9ef05SNatalie Beams qref::AbstractArray{CeedScalar}, 14180a9ef05SNatalie Beams qweight::AbstractArray{CeedScalar}, 14244554ea0SWill Pazner) 14344554ea0SWill Pazner dim = getdimension(topo) 14444554ea0SWill Pazner @assert size(interp) == (nqpts, nnodes) 14544554ea0SWill Pazner @assert size(grad) == (dim, nqpts, nnodes) 14644554ea0SWill Pazner @assert length(qref) == nqpts 14744554ea0SWill Pazner @assert length(qweight) == nqpts 14844554ea0SWill Pazner 14944554ea0SWill Pazner # Convert from Julia matrices and tensors (column-major) to row-major format 15044554ea0SWill Pazner interp_rowmajor = collect(interp') 15144554ea0SWill Pazner grad_rowmajor = permutedims(grad, [3, 2, 1]) 15244554ea0SWill Pazner 15344554ea0SWill Pazner ref = Ref{C.CeedBasis}() 15444554ea0SWill Pazner C.CeedBasisCreateH1( 15544554ea0SWill Pazner c[], 15644554ea0SWill Pazner topo, 15744554ea0SWill Pazner ncomp, 15844554ea0SWill Pazner nnodes, 15944554ea0SWill Pazner nqpts, 16044554ea0SWill Pazner interp_rowmajor, 16144554ea0SWill Pazner grad_rowmajor, 16244554ea0SWill Pazner qref, 16344554ea0SWill Pazner qweight, 16444554ea0SWill Pazner ref, 16544554ea0SWill Pazner ) 16644554ea0SWill Pazner Basis(ref) 16744554ea0SWill Paznerend 16844554ea0SWill Pazner 16944554ea0SWill Pazner""" 17044554ea0SWill Pazner apply!(b::Basis, nelem, tmode::TransposeMode, emode::EvalMode, u::AbstractCeedVector, v::AbstractCeedVector) 17144554ea0SWill Pazner 17244554ea0SWill PaznerApply basis evaluation from nodes to quadrature points or vice versa, storing the result in 17344554ea0SWill Paznerthe [`CeedVector`](@ref) `v`. 17444554ea0SWill Pazner 17544554ea0SWill Pazner`nelem` specifies the number of elements to apply the basis evaluation to; the backend will 17644554ea0SWill Paznerspecify the ordering in CeedElemRestrictionCreateBlocked() 17744554ea0SWill Pazner 17844554ea0SWill PaznerSet `tmode` to `CEED_NOTRANSPOSE` to evaluate from nodes to quadrature or to 17944554ea0SWill Pazner`CEED_TRANSPOSE` to apply the transpose, mapping from quadrature points to nodes. 18044554ea0SWill Pazner 18144554ea0SWill PaznerSet the [`EvalMode`](@ref) `emode` to: 18244554ea0SWill Pazner- `CEED_EVAL_NONE` to use values directly, 18344554ea0SWill Pazner- `CEED_EVAL_INTERP` to use interpolated values, 18444554ea0SWill Pazner- `CEED_EVAL_GRAD` to use gradients, 18544554ea0SWill Pazner- `CEED_EVAL_WEIGHT` to use quadrature weights. 18644554ea0SWill Pazner""" 18744554ea0SWill Paznerfunction apply!( 18844554ea0SWill Pazner b::Basis, 18944554ea0SWill Pazner nelem, 19044554ea0SWill Pazner tmode::TransposeMode, 19144554ea0SWill Pazner emode::EvalMode, 19244554ea0SWill Pazner u::AbstractCeedVector, 19344554ea0SWill Pazner v::AbstractCeedVector, 19444554ea0SWill Pazner) 19544554ea0SWill Pazner C.CeedBasisApply(b[], nelem, tmode, emode, u[], v[]) 19644554ea0SWill Paznerend 19744554ea0SWill Pazner 19844554ea0SWill Pazner""" 19944554ea0SWill Pazner apply(b::Basis, u::AbstractVector; nelem=1, tmode=NOTRANSPOSE, emode=EVAL_INTERP) 20044554ea0SWill Pazner 20144554ea0SWill PaznerPerforms the same function as the above-defined [`apply!`](@ref apply!(b::Basis, nelem, 20244554ea0SWill Paznertmode::TransposeMode, emode::EvalMode, u::AbstractCeedVector, v::AbstractCeedVector)), but 20344554ea0SWill Paznerautomatically convert from Julia arrays to [`CeedVector`](@ref) for convenience. 20444554ea0SWill Pazner 20544554ea0SWill PaznerThe result will be returned in a newly allocated array of the correct size. 20644554ea0SWill Pazner""" 20744554ea0SWill Paznerfunction apply(b::Basis, u::AbstractVector; nelem=1, tmode=NOTRANSPOSE, emode=EVAL_INTERP) 20844554ea0SWill Pazner ceed_ref = Ref{C.Ceed}() 20944554ea0SWill Pazner ccall((:CeedBasisGetCeed, C.libceed), Cint, (C.CeedBasis, Ptr{C.Ceed}), b[], ceed_ref) 21044554ea0SWill Pazner c = Ceed(ceed_ref) 21144554ea0SWill Pazner 21244554ea0SWill Pazner u_vec = CeedVector(c, u) 21344554ea0SWill Pazner 21444554ea0SWill Pazner len_v = (tmode == TRANSPOSE) ? getnumnodes(b) : getnumqpts(b) 21544554ea0SWill Pazner if emode == EVAL_GRAD 21644554ea0SWill Pazner len_v *= getdimension(b) 21744554ea0SWill Pazner end 21844554ea0SWill Pazner 21944554ea0SWill Pazner v_vec = CeedVector(c, len_v) 22044554ea0SWill Pazner 22144554ea0SWill Pazner apply!(b, nelem, tmode, emode, u_vec, v_vec) 22244554ea0SWill Pazner Vector(v_vec) 22344554ea0SWill Paznerend 22444554ea0SWill Pazner 22544554ea0SWill Pazner""" 22644554ea0SWill Pazner getdimension(b::Basis) 22744554ea0SWill Pazner 22844554ea0SWill PaznerReturn the spatial dimension of the given [`Basis`](@ref). 22944554ea0SWill Pazner""" 23044554ea0SWill Paznerfunction getdimension(b::Basis) 23144554ea0SWill Pazner dim = Ref{CeedInt}() 23244554ea0SWill Pazner C.CeedBasisGetDimension(b[], dim) 23344554ea0SWill Pazner dim[] 23444554ea0SWill Paznerend 23544554ea0SWill Pazner 23644554ea0SWill Pazner""" 23744554ea0SWill Pazner getdimension(t::Topology) 23844554ea0SWill Pazner 23944554ea0SWill PaznerReturn the spatial dimension of the given [`Topology`](@ref). 24044554ea0SWill Pazner""" 24144554ea0SWill Paznerfunction getdimension(t::Topology) 24244554ea0SWill Pazner return Int(t) >> 16 24344554ea0SWill Paznerend 24444554ea0SWill Pazner 24544554ea0SWill Pazner""" 24644554ea0SWill Pazner gettopology(b::Basis) 24744554ea0SWill Pazner 24844554ea0SWill PaznerReturn the [`Topology`](@ref) of the given [`Basis`](@ref). 24944554ea0SWill Pazner""" 25044554ea0SWill Paznerfunction gettopology(b::Basis) 25144554ea0SWill Pazner topo = Ref{Topology}() 25244554ea0SWill Pazner C.CeedBasisGetTopology(b[], topo) 25344554ea0SWill Pazner topo[] 25444554ea0SWill Paznerend 25544554ea0SWill Pazner 25644554ea0SWill Pazner""" 25744554ea0SWill Pazner getnumcomponents(b::Basis) 25844554ea0SWill Pazner 25944554ea0SWill PaznerReturn the number of components of the given [`Basis`](@ref). 26044554ea0SWill Pazner""" 26144554ea0SWill Paznerfunction getnumcomponents(b::Basis) 26244554ea0SWill Pazner ncomp = Ref{CeedInt}() 26344554ea0SWill Pazner C.CeedBasisGetNumComponents(b[], ncomp) 26444554ea0SWill Pazner ncomp[] 26544554ea0SWill Paznerend 26644554ea0SWill Pazner 26744554ea0SWill Pazner""" 26844554ea0SWill Pazner getnumnodes(b::Basis) 26944554ea0SWill Pazner 27044554ea0SWill PaznerReturn the number of nodes of the given [`Basis`](@ref). 27144554ea0SWill Pazner""" 27244554ea0SWill Paznerfunction getnumnodes(b::Basis) 27344554ea0SWill Pazner nnodes = Ref{CeedInt}() 27444554ea0SWill Pazner C.CeedBasisGetNumNodes(b[], nnodes) 27544554ea0SWill Pazner nnodes[] 27644554ea0SWill Paznerend 27744554ea0SWill Pazner 27844554ea0SWill Pazner""" 27944554ea0SWill Pazner getnumnodes1d(b::Basis) 28044554ea0SWill Pazner 28144554ea0SWill Pazner Return the number of 1D nodes of the given (tensor-product) [`Basis`](@ref). 28244554ea0SWill Pazner""" 28344554ea0SWill Paznerfunction getnumnodes1d(b::Basis) 28444554ea0SWill Pazner nnodes1d = Ref{CeedInt}() 28544554ea0SWill Pazner C.CeedBasisGetNumNodes1D(b[], nnodes1d) 28644554ea0SWill Pazner nnodes1d[] 28744554ea0SWill Paznerend 28844554ea0SWill Pazner 28944554ea0SWill Pazner""" 29044554ea0SWill Pazner getnumqpts(b::Basis) 29144554ea0SWill Pazner 29244554ea0SWill PaznerReturn the number of quadrature points of the given [`Basis`](@ref). 29344554ea0SWill Pazner""" 29444554ea0SWill Paznerfunction getnumqpts(b::Basis) 29544554ea0SWill Pazner nqpts = Ref{CeedInt}() 29644554ea0SWill Pazner C.CeedBasisGetNumQuadraturePoints(b[], nqpts) 29744554ea0SWill Pazner nqpts[] 29844554ea0SWill Paznerend 29944554ea0SWill Pazner 30044554ea0SWill Pazner""" 30144554ea0SWill Pazner getnumqpts1d(b::Basis) 30244554ea0SWill Pazner 30344554ea0SWill PaznerReturn the number of 1D quadrature points of the given (tensor-product) [`Basis`](@ref). 30444554ea0SWill Pazner""" 30544554ea0SWill Paznerfunction getnumqpts1d(b::Basis) 30644554ea0SWill Pazner nqpts1d = Ref{CeedInt}() 30744554ea0SWill Pazner C.CeedBasisGetNumQuadraturePoints1D(b[], nqpts1d) 30844554ea0SWill Pazner nqpts1d[] 30944554ea0SWill Paznerend 31044554ea0SWill Pazner 31144554ea0SWill Pazner""" 31244554ea0SWill Pazner getqref(b::Basis) 31344554ea0SWill Pazner 31444554ea0SWill PaznerGet the reference coordinates of quadrature points (in `dim` dimensions) of the given 31544554ea0SWill Pazner[`Basis`](@ref). 31644554ea0SWill Pazner""" 31744554ea0SWill Paznerfunction getqref(b::Basis) 318*a9e65696SJeremy L Thompson istensor = Ref{Bool}() 319*a9e65696SJeremy L Thompson C.CeedBasisIsTensor(b[], istensor) 32044554ea0SWill Pazner ref = Ref{Ptr{CeedScalar}}() 32144554ea0SWill Pazner C.CeedBasisGetQRef(b[], ref) 322*a9e65696SJeremy L Thompson copy( 323*a9e65696SJeremy L Thompson unsafe_wrap( 324*a9e65696SJeremy L Thompson Array, 325*a9e65696SJeremy L Thompson ref[], 326*a9e65696SJeremy L Thompson istensor[] ? getnumqpts1d(b) : (getnumqpts(b)*getdimension(b)), 327*a9e65696SJeremy L Thompson ), 328*a9e65696SJeremy L Thompson ) 32944554ea0SWill Paznerend 33044554ea0SWill Pazner 33144554ea0SWill Pazner""" 33244554ea0SWill Pazner getqref(b::Basis) 33344554ea0SWill Pazner 33444554ea0SWill PaznerGet the quadrature weights of quadrature points (in `dim` dimensions) of the given 33544554ea0SWill Pazner[`Basis`](@ref). 33644554ea0SWill Pazner""" 33744554ea0SWill Paznerfunction getqweights(b::Basis) 338*a9e65696SJeremy L Thompson istensor = Ref{Bool}() 339*a9e65696SJeremy L Thompson C.CeedBasisIsTensor(b[], istensor) 34044554ea0SWill Pazner ref = Ref{Ptr{CeedScalar}}() 34144554ea0SWill Pazner C.CeedBasisGetQWeights(b[], ref) 342*a9e65696SJeremy L Thompson copy(unsafe_wrap(Array, ref[], istensor[] ? getnumqpts1d(b) : getnumqpts(b))) 34344554ea0SWill Paznerend 34444554ea0SWill Pazner 34544554ea0SWill Pazner""" 34644554ea0SWill Pazner getinterp(b::Basis) 34744554ea0SWill Pazner 34844554ea0SWill PaznerGet the interpolation matrix of the given [`Basis`](@ref). Returns a matrix of size 34944554ea0SWill Pazner`(getnumqpts(b), getnumnodes(b))`. 35044554ea0SWill Pazner""" 35144554ea0SWill Paznerfunction getinterp(b::Basis) 35244554ea0SWill Pazner ref = Ref{Ptr{CeedScalar}}() 35344554ea0SWill Pazner C.CeedBasisGetInterp(b[], ref) 35444554ea0SWill Pazner q = getnumqpts(b) 35544554ea0SWill Pazner p = getnumnodes(b) 35644554ea0SWill Pazner collect(unsafe_wrap(Array, ref[], (p, q))') 35744554ea0SWill Paznerend 35844554ea0SWill Pazner 35944554ea0SWill Pazner""" 36044554ea0SWill Pazner getinterp1d(b::Basis) 36144554ea0SWill Pazner 36244554ea0SWill PaznerGet the 1D interpolation matrix of the given [`Basis`](@ref). `b` must be a tensor-product 36344554ea0SWill Paznerbasis, otherwise this function will fail. Returns a matrix of size `(getnumqpts1d(b), 36444554ea0SWill Paznergetnumnodes1d(b))`. 36544554ea0SWill Pazner""" 36644554ea0SWill Paznerfunction getinterp1d(b::Basis) 36744554ea0SWill Pazner ref = Ref{Ptr{CeedScalar}}() 36844554ea0SWill Pazner C.CeedBasisGetInterp1D(b[], ref) 36944554ea0SWill Pazner q = getnumqpts1d(b) 37044554ea0SWill Pazner p = getnumnodes1d(b) 37144554ea0SWill Pazner collect(unsafe_wrap(Array, ref[], (p, q))') 37244554ea0SWill Paznerend 37344554ea0SWill Pazner 37444554ea0SWill Pazner""" 37544554ea0SWill Pazner getgad(b::Basis) 37644554ea0SWill Pazner 37744554ea0SWill PaznerGet the gradient matrix of the given [`Basis`](@ref). Returns a tensor of size 37844554ea0SWill Pazner`(getdimension(b), getnumqpts(b), getnumnodes(b))`. 37944554ea0SWill Pazner""" 38044554ea0SWill Paznerfunction getgrad(b::Basis) 38144554ea0SWill Pazner ref = Ref{Ptr{CeedScalar}}() 38244554ea0SWill Pazner C.CeedBasisGetGrad(b[], ref) 38344554ea0SWill Pazner dim = getdimension(b) 38444554ea0SWill Pazner q = getnumqpts(b) 38544554ea0SWill Pazner p = getnumnodes(b) 38644554ea0SWill Pazner permutedims(unsafe_wrap(Array, ref[], (p, q, dim)), [3, 2, 1]) 38744554ea0SWill Paznerend 38844554ea0SWill Pazner 38944554ea0SWill Pazner""" 39044554ea0SWill Pazner getgrad1d(b::Basis) 39144554ea0SWill Pazner 39244554ea0SWill PaznerGet the 1D derivative matrix of the given [`Basis`](@ref). Returns a matrix of size 39344554ea0SWill Pazner`(getnumqpts(b), getnumnodes(b))`. 39444554ea0SWill Pazner""" 39544554ea0SWill Paznerfunction getgrad1d(b::Basis) 39644554ea0SWill Pazner ref = Ref{Ptr{CeedScalar}}() 39744554ea0SWill Pazner C.CeedBasisGetGrad1D(b[], ref) 39844554ea0SWill Pazner q = getnumqpts1d(b) 39944554ea0SWill Pazner p = getnumnodes1d(b) 40044554ea0SWill Pazner collect(unsafe_wrap(Array, ref[], (p, q))') 40144554ea0SWill Paznerend 402