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) 2011b88ddaSSebastian Grimberg- [`create_hdiv_basis`](@ref) 2111b88ddaSSebastian 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 11711b88ddaSSebastian 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. 129*b2e3f8ecSSebastian Grimberg- `qref`: Matrix of size `(dim, 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) 148*b2e3f8ecSSebastian Grimberg @assert size(qref) == (dim, 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]) 154*b2e3f8ecSSebastian Grimberg qref_rowmajor = collect(qref') 15544554ea0SWill Pazner 15644554ea0SWill Pazner ref = Ref{C.CeedBasis}() 15744554ea0SWill Pazner C.CeedBasisCreateH1( 15844554ea0SWill Pazner c[], 15944554ea0SWill Pazner topo, 16044554ea0SWill Pazner ncomp, 16144554ea0SWill Pazner nnodes, 16244554ea0SWill Pazner nqpts, 16344554ea0SWill Pazner interp_rowmajor, 16444554ea0SWill Pazner grad_rowmajor, 165*b2e3f8ecSSebastian Grimberg qref_rowmajor, 16644554ea0SWill Pazner qweight, 16744554ea0SWill Pazner ref, 16844554ea0SWill Pazner ) 16944554ea0SWill Pazner Basis(ref) 17044554ea0SWill Paznerend 17144554ea0SWill Pazner 17211b88ddaSSebastian Grimberg@doc raw""" 17311b88ddaSSebastian Grimberg create_hdiv_basis(c::Ceed, topo::Topology, ncomp, nnodes, nqpts, interp, div, qref, qweight) 17411b88ddaSSebastian Grimberg 17511b88ddaSSebastian GrimbergCreate a non tensor-product basis for H(div) discretizations 17611b88ddaSSebastian Grimberg 17711b88ddaSSebastian Grimberg# Arguments: 17811b88ddaSSebastian Grimberg- `ceed`: A [`Ceed`](@ref) object where the [`Basis`](@ref) will be created. 17911b88ddaSSebastian Grimberg- `topo`: [`Topology`](@ref) of element, e.g. hypercube, simplex, etc. 18011b88ddaSSebastian Grimberg- `ncomp`: Number of field components (1 for scalar fields). 18111b88ddaSSebastian Grimberg- `nnodes`: Total number of nodes. 18211b88ddaSSebastian Grimberg- `nqpts`: Total number of quadrature points. 18311b88ddaSSebastian Grimberg- `interp`: Matrix of size `(dim, nqpts, nnodes)` expressing the values of basis functions 18411b88ddaSSebastian Grimberg at quadrature points. 18511b88ddaSSebastian Grimberg- `div`: Array of size `(nqpts, nnodes)` expressing divergence of basis functions at 18611b88ddaSSebastian Grimberg quadrature points. 187*b2e3f8ecSSebastian Grimberg- `qref`: Matrix of size `(dim, nqpts)` holding the locations of quadrature points on the 18811b88ddaSSebastian Grimberg reference element $[-1, 1]$. 18911b88ddaSSebastian Grimberg- `qweight`: Array of length `nqpts` holding the quadrature weights on the reference 19011b88ddaSSebastian Grimberg element. 19111b88ddaSSebastian Grimberg""" 19211b88ddaSSebastian Grimbergfunction create_hdiv_basis( 19311b88ddaSSebastian Grimberg c::Ceed, 19411b88ddaSSebastian Grimberg topo::Topology, 19511b88ddaSSebastian Grimberg ncomp, 19611b88ddaSSebastian Grimberg nnodes, 19711b88ddaSSebastian Grimberg nqpts, 19811b88ddaSSebastian Grimberg interp::AbstractArray{CeedScalar}, 19911b88ddaSSebastian Grimberg div::AbstractArray{CeedScalar}, 20011b88ddaSSebastian Grimberg qref::AbstractArray{CeedScalar}, 20111b88ddaSSebastian Grimberg qweight::AbstractArray{CeedScalar}, 20211b88ddaSSebastian Grimberg) 20311b88ddaSSebastian Grimberg dim = getdimension(topo) 20411b88ddaSSebastian Grimberg @assert size(interp) == (dim, nqpts, nnodes) 20511b88ddaSSebastian Grimberg @assert size(div) == (nqpts, nnodes) 206*b2e3f8ecSSebastian Grimberg @assert size(qref) == (dim, nqpts) 20711b88ddaSSebastian Grimberg @assert length(qweight) == nqpts 20811b88ddaSSebastian Grimberg 20911b88ddaSSebastian Grimberg # Convert from Julia matrices and tensors (column-major) to row-major format 21011b88ddaSSebastian Grimberg interp_rowmajor = permutedims(interp, [3, 2, 1]) 21111b88ddaSSebastian Grimberg div_rowmajor = collect(div') 212*b2e3f8ecSSebastian Grimberg qref_rowmajor = collect(qref') 21311b88ddaSSebastian Grimberg 21411b88ddaSSebastian Grimberg ref = Ref{C.CeedBasis}() 21511b88ddaSSebastian Grimberg C.CeedBasisCreateHdiv( 21611b88ddaSSebastian Grimberg c[], 21711b88ddaSSebastian Grimberg topo, 21811b88ddaSSebastian Grimberg ncomp, 21911b88ddaSSebastian Grimberg nnodes, 22011b88ddaSSebastian Grimberg nqpts, 22111b88ddaSSebastian Grimberg interp_rowmajor, 22211b88ddaSSebastian Grimberg div_rowmajor, 223*b2e3f8ecSSebastian Grimberg qref_rowmajor, 22411b88ddaSSebastian Grimberg qweight, 22511b88ddaSSebastian Grimberg ref, 22611b88ddaSSebastian Grimberg ) 22711b88ddaSSebastian Grimberg Basis(ref) 22811b88ddaSSebastian Grimbergend 22911b88ddaSSebastian Grimberg 23011b88ddaSSebastian Grimberg@doc raw""" 23111b88ddaSSebastian Grimberg create_hcurl_basis(c::Ceed, topo::Topology, ncomp, nnodes, nqpts, interp, curl, qref, qweight) 23211b88ddaSSebastian Grimberg 23311b88ddaSSebastian GrimbergCreate a non tensor-product basis for H(curl) discretizations 23411b88ddaSSebastian Grimberg 23511b88ddaSSebastian Grimberg# Arguments: 23611b88ddaSSebastian Grimberg- `ceed`: A [`Ceed`](@ref) object where the [`Basis`](@ref) will be created. 23711b88ddaSSebastian Grimberg- `topo`: [`Topology`](@ref) of element, e.g. hypercube, simplex, etc. 23811b88ddaSSebastian Grimberg- `ncomp`: Number of field components (1 for scalar fields). 23911b88ddaSSebastian Grimberg- `nnodes`: Total number of nodes. 24011b88ddaSSebastian Grimberg- `nqpts`: Total number of quadrature points. 24111b88ddaSSebastian Grimberg- `interp`: Matrix of size `(dim, nqpts, nnodes)` expressing the values of basis functions 24211b88ddaSSebastian Grimberg at quadrature points. 24311b88ddaSSebastian Grimberg- `curl`: Matrix of size `(curlcomp, nqpts, nnodes)`, `curlcomp = 1 if dim < 3 else dim`) 24411b88ddaSSebastian Grimberg matrix expressing curl of basis functions at quadrature points. 245*b2e3f8ecSSebastian Grimberg- `qref`: Matrix of size `(dim, nqpts)` holding the locations of quadrature points on the 24611b88ddaSSebastian Grimberg reference element $[-1, 1]$. 24711b88ddaSSebastian Grimberg- `qweight`: Array of length `nqpts` holding the quadrature weights on the reference 24811b88ddaSSebastian Grimberg element. 24911b88ddaSSebastian Grimberg""" 25011b88ddaSSebastian Grimbergfunction create_hcurl_basis( 25111b88ddaSSebastian Grimberg c::Ceed, 25211b88ddaSSebastian Grimberg topo::Topology, 25311b88ddaSSebastian Grimberg ncomp, 25411b88ddaSSebastian Grimberg nnodes, 25511b88ddaSSebastian Grimberg nqpts, 25611b88ddaSSebastian Grimberg interp::AbstractArray{CeedScalar}, 25711b88ddaSSebastian Grimberg curl::AbstractArray{CeedScalar}, 25811b88ddaSSebastian Grimberg qref::AbstractArray{CeedScalar}, 25911b88ddaSSebastian Grimberg qweight::AbstractArray{CeedScalar}, 26011b88ddaSSebastian Grimberg) 26111b88ddaSSebastian Grimberg dim = getdimension(topo) 26211b88ddaSSebastian Grimberg curlcomp = dim < 3 ? 1 : dim 26311b88ddaSSebastian Grimberg @assert size(interp) == (dim, nqpts, nnodes) 26411b88ddaSSebastian Grimberg @assert size(curl) == (curlcomp, nqpts, nnodes) 265*b2e3f8ecSSebastian Grimberg @assert size(qref) == (dim, nqpts) 26611b88ddaSSebastian Grimberg @assert length(qweight) == nqpts 26711b88ddaSSebastian Grimberg 26811b88ddaSSebastian Grimberg # Convert from Julia matrices and tensors (column-major) to row-major format 26911b88ddaSSebastian Grimberg interp_rowmajor = permutedims(interp, [3, 2, 1]) 27011b88ddaSSebastian Grimberg curl_rowmajor = permutedims(curl, [3, 2, 1]) 271*b2e3f8ecSSebastian Grimberg qref_rowmajor = collect(qref') 27211b88ddaSSebastian Grimberg 27311b88ddaSSebastian Grimberg ref = Ref{C.CeedBasis}() 27411b88ddaSSebastian Grimberg C.CeedBasisCreateHcurl( 27511b88ddaSSebastian Grimberg c[], 27611b88ddaSSebastian Grimberg topo, 27711b88ddaSSebastian Grimberg ncomp, 27811b88ddaSSebastian Grimberg nnodes, 27911b88ddaSSebastian Grimberg nqpts, 28011b88ddaSSebastian Grimberg interp_rowmajor, 28111b88ddaSSebastian Grimberg curl_rowmajor, 282*b2e3f8ecSSebastian Grimberg qref_rowmajor, 28311b88ddaSSebastian Grimberg qweight, 28411b88ddaSSebastian Grimberg ref, 28511b88ddaSSebastian Grimberg ) 28611b88ddaSSebastian Grimberg Basis(ref) 28711b88ddaSSebastian Grimbergend 28811b88ddaSSebastian Grimberg 28944554ea0SWill Pazner""" 29044554ea0SWill Pazner apply!(b::Basis, nelem, tmode::TransposeMode, emode::EvalMode, u::AbstractCeedVector, v::AbstractCeedVector) 29144554ea0SWill Pazner 29244554ea0SWill PaznerApply basis evaluation from nodes to quadrature points or vice versa, storing the result in 29344554ea0SWill Paznerthe [`CeedVector`](@ref) `v`. 29444554ea0SWill Pazner 29544554ea0SWill Pazner`nelem` specifies the number of elements to apply the basis evaluation to; the backend will 29644554ea0SWill Paznerspecify the ordering in CeedElemRestrictionCreateBlocked() 29744554ea0SWill Pazner 29844554ea0SWill PaznerSet `tmode` to `CEED_NOTRANSPOSE` to evaluate from nodes to quadrature or to 29944554ea0SWill Pazner`CEED_TRANSPOSE` to apply the transpose, mapping from quadrature points to nodes. 30044554ea0SWill Pazner 30144554ea0SWill PaznerSet the [`EvalMode`](@ref) `emode` to: 30244554ea0SWill Pazner- `CEED_EVAL_NONE` to use values directly, 30344554ea0SWill Pazner- `CEED_EVAL_INTERP` to use interpolated values, 30444554ea0SWill Pazner- `CEED_EVAL_GRAD` to use gradients, 30544554ea0SWill Pazner- `CEED_EVAL_WEIGHT` to use quadrature weights. 30644554ea0SWill Pazner""" 30744554ea0SWill Paznerfunction apply!( 30844554ea0SWill Pazner b::Basis, 30944554ea0SWill Pazner nelem, 31044554ea0SWill Pazner tmode::TransposeMode, 31144554ea0SWill Pazner emode::EvalMode, 31244554ea0SWill Pazner u::AbstractCeedVector, 31344554ea0SWill Pazner v::AbstractCeedVector, 31444554ea0SWill Pazner) 31544554ea0SWill Pazner C.CeedBasisApply(b[], nelem, tmode, emode, u[], v[]) 31644554ea0SWill Paznerend 31744554ea0SWill Pazner 31844554ea0SWill Pazner""" 31944554ea0SWill Pazner apply(b::Basis, u::AbstractVector; nelem=1, tmode=NOTRANSPOSE, emode=EVAL_INTERP) 32044554ea0SWill Pazner 32144554ea0SWill PaznerPerforms the same function as the above-defined [`apply!`](@ref apply!(b::Basis, nelem, 32244554ea0SWill Paznertmode::TransposeMode, emode::EvalMode, u::AbstractCeedVector, v::AbstractCeedVector)), but 32344554ea0SWill Paznerautomatically convert from Julia arrays to [`CeedVector`](@ref) for convenience. 32444554ea0SWill Pazner 32544554ea0SWill PaznerThe result will be returned in a newly allocated array of the correct size. 32644554ea0SWill Pazner""" 32744554ea0SWill Paznerfunction apply(b::Basis, u::AbstractVector; nelem=1, tmode=NOTRANSPOSE, emode=EVAL_INTERP) 32844554ea0SWill Pazner ceed_ref = Ref{C.Ceed}() 32944554ea0SWill Pazner ccall((:CeedBasisGetCeed, C.libceed), Cint, (C.CeedBasis, Ptr{C.Ceed}), b[], ceed_ref) 33044554ea0SWill Pazner c = Ceed(ceed_ref) 33144554ea0SWill Pazner 33244554ea0SWill Pazner u_vec = CeedVector(c, u) 33344554ea0SWill Pazner 334*b2e3f8ecSSebastian Grimberg qcomp = Ref{CeedInt}() 335*b2e3f8ecSSebastian Grimberg C.CeedBasisGetNumQuadratureComponents(b[], emode, qcomp) 336*b2e3f8ecSSebastian Grimberg len_v = (tmode == TRANSPOSE) ? getnumnodes(b) : qcomp[]*getnumqpts(b) 33744554ea0SWill Pazner 33844554ea0SWill Pazner v_vec = CeedVector(c, len_v) 33944554ea0SWill Pazner 34044554ea0SWill Pazner apply!(b, nelem, tmode, emode, u_vec, v_vec) 34144554ea0SWill Pazner Vector(v_vec) 34244554ea0SWill Paznerend 34344554ea0SWill Pazner 34444554ea0SWill Pazner""" 34544554ea0SWill Pazner getdimension(b::Basis) 34644554ea0SWill Pazner 34744554ea0SWill PaznerReturn the spatial dimension of the given [`Basis`](@ref). 34844554ea0SWill Pazner""" 34944554ea0SWill Paznerfunction getdimension(b::Basis) 35044554ea0SWill Pazner dim = Ref{CeedInt}() 35144554ea0SWill Pazner C.CeedBasisGetDimension(b[], dim) 35244554ea0SWill Pazner dim[] 35344554ea0SWill Paznerend 35444554ea0SWill Pazner 35544554ea0SWill Pazner""" 35644554ea0SWill Pazner getdimension(t::Topology) 35744554ea0SWill Pazner 35844554ea0SWill PaznerReturn the spatial dimension of the given [`Topology`](@ref). 35944554ea0SWill Pazner""" 36044554ea0SWill Paznerfunction getdimension(t::Topology) 36144554ea0SWill Pazner return Int(t) >> 16 36244554ea0SWill Paznerend 36344554ea0SWill Pazner 36444554ea0SWill Pazner""" 36544554ea0SWill Pazner gettopology(b::Basis) 36644554ea0SWill Pazner 36744554ea0SWill PaznerReturn the [`Topology`](@ref) of the given [`Basis`](@ref). 36844554ea0SWill Pazner""" 36944554ea0SWill Paznerfunction gettopology(b::Basis) 37044554ea0SWill Pazner topo = Ref{Topology}() 37144554ea0SWill Pazner C.CeedBasisGetTopology(b[], topo) 37244554ea0SWill Pazner topo[] 37344554ea0SWill Paznerend 37444554ea0SWill Pazner 37544554ea0SWill Pazner""" 37644554ea0SWill Pazner getnumcomponents(b::Basis) 37744554ea0SWill Pazner 37844554ea0SWill PaznerReturn the number of components of the given [`Basis`](@ref). 37944554ea0SWill Pazner""" 38044554ea0SWill Paznerfunction getnumcomponents(b::Basis) 38144554ea0SWill Pazner ncomp = Ref{CeedInt}() 38244554ea0SWill Pazner C.CeedBasisGetNumComponents(b[], ncomp) 38344554ea0SWill Pazner ncomp[] 38444554ea0SWill Paznerend 38544554ea0SWill Pazner 38644554ea0SWill Pazner""" 38744554ea0SWill Pazner getnumnodes(b::Basis) 38844554ea0SWill Pazner 38944554ea0SWill PaznerReturn the number of nodes of the given [`Basis`](@ref). 39044554ea0SWill Pazner""" 39144554ea0SWill Paznerfunction getnumnodes(b::Basis) 39244554ea0SWill Pazner nnodes = Ref{CeedInt}() 39344554ea0SWill Pazner C.CeedBasisGetNumNodes(b[], nnodes) 39444554ea0SWill Pazner nnodes[] 39544554ea0SWill Paznerend 39644554ea0SWill Pazner 39744554ea0SWill Pazner""" 39844554ea0SWill Pazner getnumnodes1d(b::Basis) 39944554ea0SWill Pazner 40044554ea0SWill Pazner Return the number of 1D nodes of the given (tensor-product) [`Basis`](@ref). 40144554ea0SWill Pazner""" 40244554ea0SWill Paznerfunction getnumnodes1d(b::Basis) 40344554ea0SWill Pazner nnodes1d = Ref{CeedInt}() 40444554ea0SWill Pazner C.CeedBasisGetNumNodes1D(b[], nnodes1d) 40544554ea0SWill Pazner nnodes1d[] 40644554ea0SWill Paznerend 40744554ea0SWill Pazner 40844554ea0SWill Pazner""" 40944554ea0SWill Pazner getnumqpts(b::Basis) 41044554ea0SWill Pazner 41144554ea0SWill PaznerReturn the number of quadrature points of the given [`Basis`](@ref). 41244554ea0SWill Pazner""" 41344554ea0SWill Paznerfunction getnumqpts(b::Basis) 41444554ea0SWill Pazner nqpts = Ref{CeedInt}() 41544554ea0SWill Pazner C.CeedBasisGetNumQuadraturePoints(b[], nqpts) 41644554ea0SWill Pazner nqpts[] 41744554ea0SWill Paznerend 41844554ea0SWill Pazner 41944554ea0SWill Pazner""" 42044554ea0SWill Pazner getnumqpts1d(b::Basis) 42144554ea0SWill Pazner 42244554ea0SWill PaznerReturn the number of 1D quadrature points of the given (tensor-product) [`Basis`](@ref). 42344554ea0SWill Pazner""" 42444554ea0SWill Paznerfunction getnumqpts1d(b::Basis) 42544554ea0SWill Pazner nqpts1d = Ref{CeedInt}() 42644554ea0SWill Pazner C.CeedBasisGetNumQuadraturePoints1D(b[], nqpts1d) 42744554ea0SWill Pazner nqpts1d[] 42844554ea0SWill Paznerend 42944554ea0SWill Pazner 43044554ea0SWill Pazner""" 43144554ea0SWill Pazner getqref(b::Basis) 43244554ea0SWill Pazner 43344554ea0SWill PaznerGet the reference coordinates of quadrature points (in `dim` dimensions) of the given 43444554ea0SWill Pazner[`Basis`](@ref). 43544554ea0SWill Pazner""" 43644554ea0SWill Paznerfunction getqref(b::Basis) 437a9e65696SJeremy L Thompson istensor = Ref{Bool}() 438a9e65696SJeremy L Thompson C.CeedBasisIsTensor(b[], istensor) 43944554ea0SWill Pazner ref = Ref{Ptr{CeedScalar}}() 44044554ea0SWill Pazner C.CeedBasisGetQRef(b[], ref) 441a9e65696SJeremy L Thompson copy( 442*b2e3f8ecSSebastian Grimberg istensor[] ? unsafe_wrap(Array, ref[], getnumqpts1d(b)) : 443*b2e3f8ecSSebastian Grimberg unsafe_wrap(Array, ref[], (getnumqpts(b), getdimension(b)))', 444a9e65696SJeremy L Thompson ) 44544554ea0SWill Paznerend 44644554ea0SWill Pazner 44744554ea0SWill Pazner""" 44844554ea0SWill Pazner getqref(b::Basis) 44944554ea0SWill Pazner 45044554ea0SWill PaznerGet the quadrature weights of quadrature points (in `dim` dimensions) of the given 45144554ea0SWill Pazner[`Basis`](@ref). 45244554ea0SWill Pazner""" 45344554ea0SWill Paznerfunction getqweights(b::Basis) 454a9e65696SJeremy L Thompson istensor = Ref{Bool}() 455a9e65696SJeremy L Thompson C.CeedBasisIsTensor(b[], istensor) 45644554ea0SWill Pazner ref = Ref{Ptr{CeedScalar}}() 45744554ea0SWill Pazner C.CeedBasisGetQWeights(b[], ref) 458a9e65696SJeremy L Thompson copy(unsafe_wrap(Array, ref[], istensor[] ? getnumqpts1d(b) : getnumqpts(b))) 45944554ea0SWill Paznerend 46044554ea0SWill Pazner 46111b88ddaSSebastian Grimberg@doc raw""" 46244554ea0SWill Pazner getinterp(b::Basis) 46344554ea0SWill Pazner 46444554ea0SWill PaznerGet the interpolation matrix of the given [`Basis`](@ref). Returns a matrix of size 46511b88ddaSSebastian Grimberg`(getnumqpts(b), getnumnodes(b))` for a given $H^1$ basis or `(getdimension(b), 46611b88ddaSSebastian Grimberggetnumqpts(b), getnumnodes(b))` for a given vector $H(div)$ or $H(curl)$ basis. 46744554ea0SWill Pazner""" 46844554ea0SWill Paznerfunction getinterp(b::Basis) 46944554ea0SWill Pazner ref = Ref{Ptr{CeedScalar}}() 47044554ea0SWill Pazner C.CeedBasisGetInterp(b[], ref) 47144554ea0SWill Pazner q = getnumqpts(b) 47244554ea0SWill Pazner p = getnumnodes(b) 47311b88ddaSSebastian Grimberg qcomp = Ref{CeedInt}() 47411b88ddaSSebastian Grimberg C.CeedBasisGetNumQuadratureComponents(b[], C.CEED_EVAL_INTERP, qcomp) 47511b88ddaSSebastian Grimberg if qcomp[] == 1 47644554ea0SWill Pazner collect(unsafe_wrap(Array, ref[], (p, q))') 47711b88ddaSSebastian Grimberg else 47811b88ddaSSebastian Grimberg permutedims(unsafe_wrap(Array, ref[], (p, q, qcomp[])), [3, 2, 1]) 47911b88ddaSSebastian Grimberg end 48044554ea0SWill Paznerend 48144554ea0SWill Pazner 48244554ea0SWill Pazner""" 48344554ea0SWill Pazner getinterp1d(b::Basis) 48444554ea0SWill Pazner 48544554ea0SWill PaznerGet the 1D interpolation matrix of the given [`Basis`](@ref). `b` must be a tensor-product 48644554ea0SWill Paznerbasis, otherwise this function will fail. Returns a matrix of size `(getnumqpts1d(b), 48744554ea0SWill Paznergetnumnodes1d(b))`. 48844554ea0SWill Pazner""" 48944554ea0SWill Paznerfunction getinterp1d(b::Basis) 49044554ea0SWill Pazner ref = Ref{Ptr{CeedScalar}}() 49144554ea0SWill Pazner C.CeedBasisGetInterp1D(b[], ref) 49244554ea0SWill Pazner q = getnumqpts1d(b) 49344554ea0SWill Pazner p = getnumnodes1d(b) 49444554ea0SWill Pazner collect(unsafe_wrap(Array, ref[], (p, q))') 49544554ea0SWill Paznerend 49644554ea0SWill Pazner 49744554ea0SWill Pazner""" 4987574393bSJeremy L Thompson getgrad(b::Basis) 49944554ea0SWill Pazner 50044554ea0SWill PaznerGet the gradient matrix of the given [`Basis`](@ref). Returns a tensor of size 50144554ea0SWill Pazner`(getdimension(b), getnumqpts(b), getnumnodes(b))`. 50244554ea0SWill Pazner""" 50344554ea0SWill Paznerfunction getgrad(b::Basis) 50444554ea0SWill Pazner ref = Ref{Ptr{CeedScalar}}() 50544554ea0SWill Pazner C.CeedBasisGetGrad(b[], ref) 50644554ea0SWill Pazner dim = getdimension(b) 50744554ea0SWill Pazner q = getnumqpts(b) 50844554ea0SWill Pazner p = getnumnodes(b) 50944554ea0SWill Pazner permutedims(unsafe_wrap(Array, ref[], (p, q, dim)), [3, 2, 1]) 51044554ea0SWill Paznerend 51144554ea0SWill Pazner 51244554ea0SWill Pazner""" 51344554ea0SWill Pazner getgrad1d(b::Basis) 51444554ea0SWill Pazner 51544554ea0SWill PaznerGet the 1D derivative matrix of the given [`Basis`](@ref). Returns a matrix of size 51644554ea0SWill Pazner`(getnumqpts(b), getnumnodes(b))`. 51744554ea0SWill Pazner""" 51844554ea0SWill Paznerfunction getgrad1d(b::Basis) 51944554ea0SWill Pazner ref = Ref{Ptr{CeedScalar}}() 52044554ea0SWill Pazner C.CeedBasisGetGrad1D(b[], ref) 52144554ea0SWill Pazner q = getnumqpts1d(b) 52244554ea0SWill Pazner p = getnumnodes1d(b) 52344554ea0SWill Pazner collect(unsafe_wrap(Array, ref[], (p, q))') 52444554ea0SWill Paznerend 52511b88ddaSSebastian Grimberg 52611b88ddaSSebastian Grimberg""" 52711b88ddaSSebastian Grimberg getdiv(b::Basis) 52811b88ddaSSebastian Grimberg 52911b88ddaSSebastian GrimbergGet the divergence matrix of the given [`Basis`](@ref). Returns a tensor of size 53011b88ddaSSebastian Grimberg`(getnumqpts(b), getnumnodes(b))`. 53111b88ddaSSebastian Grimberg""" 53211b88ddaSSebastian Grimbergfunction getdiv(b::Basis) 53311b88ddaSSebastian Grimberg ref = Ref{Ptr{CeedScalar}}() 53411b88ddaSSebastian Grimberg C.CeedBasisGetDiv(b[], ref) 53511b88ddaSSebastian Grimberg q = getnumqpts(b) 53611b88ddaSSebastian Grimberg p = getnumnodes(b) 53711b88ddaSSebastian Grimberg collect(unsafe_wrap(Array, ref[], (p, q))') 53811b88ddaSSebastian Grimbergend 53911b88ddaSSebastian Grimberg 54011b88ddaSSebastian Grimberg""" 54111b88ddaSSebastian Grimberg getcurl(b::Basis) 54211b88ddaSSebastian Grimberg 54311b88ddaSSebastian GrimbergGet the curl matrix of the given [`Basis`](@ref). Returns a tensor of size 54411b88ddaSSebastian Grimberg`(curlcomp, getnumqpts(b), getnumnodes(b))`, `curlcomp = 1 if getdimension(b) < 3 else 54511b88ddaSSebastian Grimberggetdimension(b)`. 54611b88ddaSSebastian Grimberg""" 54711b88ddaSSebastian Grimbergfunction getcurl(b::Basis) 54811b88ddaSSebastian Grimberg ref = Ref{Ptr{CeedScalar}}() 54911b88ddaSSebastian Grimberg C.CeedBasisGetCurl(b[], ref) 55011b88ddaSSebastian Grimberg q = getnumqpts(b) 55111b88ddaSSebastian Grimberg p = getnumnodes(b) 55211b88ddaSSebastian Grimberg qcomp = Ref{CeedInt}() 55311b88ddaSSebastian Grimberg C.CeedBasisGetNumQuadratureComponents(b[], C.CEED_EVAL_CURL, qcomp) 55411b88ddaSSebastian Grimberg permutedims(unsafe_wrap(Array, ref[], (p, q, qcomp[])), [3, 2, 1]) 55511b88ddaSSebastian Grimbergend 556