1*44554ea0SWill Paznerabstract type AbstractBasis end 2*44554ea0SWill Pazner 3*44554ea0SWill Pazner""" 4*44554ea0SWill Pazner BasisCollocated() 5*44554ea0SWill Pazner 6*44554ea0SWill PaznerReturns the singleton object corresponding to libCEED's `CEED_BASIS_COLLOCATED`. 7*44554ea0SWill Pazner""" 8*44554ea0SWill Paznerstruct BasisCollocated <: AbstractBasis end 9*44554ea0SWill PaznerBase.getindex(::BasisCollocated) = C.CEED_BASIS_COLLOCATED[] 10*44554ea0SWill Pazner 11*44554ea0SWill Pazner""" 12*44554ea0SWill Pazner Basis 13*44554ea0SWill Pazner 14*44554ea0SWill PaznerWraps a `CeedBasis` object, representing a finite element basis. A `Basis` object can be 15*44554ea0SWill Paznercreated using one of: 16*44554ea0SWill Pazner 17*44554ea0SWill Pazner- [`create_tensor_h1_lagrange_basis`](@ref) 18*44554ea0SWill Pazner- [`create_tensor_h1_basis`](@ref) 19*44554ea0SWill Pazner- [`create_h1_basis`](@ref) 20*44554ea0SWill Pazner""" 21*44554ea0SWill Paznermutable struct Basis <: AbstractBasis 22*44554ea0SWill Pazner ref::RefValue{C.CeedBasis} 23*44554ea0SWill Pazner function Basis(ref) 24*44554ea0SWill Pazner obj = new(ref) 25*44554ea0SWill Pazner finalizer(obj) do x 26*44554ea0SWill Pazner # ccall(:jl_safe_printf, Cvoid, (Cstring, Cstring), "Finalizing %s.\n", repr(x)) 27*44554ea0SWill Pazner destroy(x) 28*44554ea0SWill Pazner end 29*44554ea0SWill Pazner return obj 30*44554ea0SWill Pazner end 31*44554ea0SWill Paznerend 32*44554ea0SWill Paznerdestroy(b::Basis) = C.CeedBasisDestroy(b.ref) # COV_EXCL_LINE 33*44554ea0SWill PaznerBase.getindex(b::Basis) = b.ref[] 34*44554ea0SWill PaznerBase.show(io::IO, ::MIME"text/plain", b::Basis) = ceed_show(io, b, C.CeedBasisView) 35*44554ea0SWill Pazner 36*44554ea0SWill Pazner@doc raw""" 37*44554ea0SWill Pazner create_tensor_h1_lagrange_basis(ceed, dim, ncomp, p, q, qmode) 38*44554ea0SWill Pazner 39*44554ea0SWill PaznerCreate a tensor-product Lagrange basis. 40*44554ea0SWill Pazner 41*44554ea0SWill Pazner# Arguments: 42*44554ea0SWill Pazner- `ceed`: A [`Ceed`](@ref) object where the [`Basis`](@ref) will be created. 43*44554ea0SWill Pazner- `dim`: Topological dimension of element. 44*44554ea0SWill Pazner- `ncomp`: Number of field components (1 for scalar fields). 45*44554ea0SWill Pazner- `p`: Number of Gauss-Lobatto nodes in one dimension. The polynomial degree of the 46*44554ea0SWill Pazner resulting $Q_k$ element is $k=p-1$. 47*44554ea0SWill Pazner- `q`: Number of quadrature points in one dimension. 48*44554ea0SWill Pazner- `qmode`: Distribution of the $q$ quadrature points (affects order of accuracy for the 49*44554ea0SWill Pazner quadrature). 50*44554ea0SWill Pazner""" 51*44554ea0SWill Paznerfunction create_tensor_h1_lagrange_basis(c::Ceed, dim, ncomp, p, q, quad_mode::QuadMode) 52*44554ea0SWill Pazner ref = Ref{C.CeedBasis}() 53*44554ea0SWill Pazner C.CeedBasisCreateTensorH1Lagrange(c[], dim, ncomp, p, q, quad_mode, ref) 54*44554ea0SWill Pazner Basis(ref) 55*44554ea0SWill Paznerend 56*44554ea0SWill Pazner 57*44554ea0SWill Pazner@doc raw""" 58*44554ea0SWill Pazner create_tensor_h1_basis(c::Ceed, dim, ncomp, p, q, interp1d, grad1d, qref1d, qweight1d) 59*44554ea0SWill Pazner 60*44554ea0SWill PaznerCreate a tensor-product basis for $H^1$ discretizations. 61*44554ea0SWill Pazner 62*44554ea0SWill Pazner# Arguments: 63*44554ea0SWill Pazner- `ceed`: A [`Ceed`](@ref) object where the [`Basis`](@ref) will be created. 64*44554ea0SWill Pazner- `dim`: Topological dimension. 65*44554ea0SWill Pazner- `ncomp`: Number of field components (1 for scalar fields). 66*44554ea0SWill Pazner- `p`: Number of nodes in one dimension. 67*44554ea0SWill Pazner- `q`: Number of quadrature points in one dimension 68*44554ea0SWill Pazner- `interp1d`: Matrix of size `(q, p)` expressing the values of nodal basis functions at 69*44554ea0SWill Pazner quadrature points. 70*44554ea0SWill Pazner- `grad1d`: Matrix of size `(p, q)` expressing derivatives of nodal basis functions at 71*44554ea0SWill Pazner quadrature points. 72*44554ea0SWill Pazner- `qref1d`: Array of length `q` holding the locations of quadrature points on the 1D 73*44554ea0SWill Pazner reference element $[-1, 1]$. 74*44554ea0SWill Pazner- `qweight1d`: Array of length `q` holding the quadrature weights on the reference element. 75*44554ea0SWill Pazner""" 76*44554ea0SWill Paznerfunction create_tensor_h1_basis( 77*44554ea0SWill Pazner c::Ceed, 78*44554ea0SWill Pazner dim, 79*44554ea0SWill Pazner ncomp, 80*44554ea0SWill Pazner p, 81*44554ea0SWill Pazner q, 82*44554ea0SWill Pazner interp1d, 83*44554ea0SWill Pazner grad1d, 84*44554ea0SWill Pazner qref1d, 85*44554ea0SWill Pazner qweight1d, 86*44554ea0SWill Pazner) 87*44554ea0SWill Pazner @assert size(interp1d) == (q, p) 88*44554ea0SWill Pazner @assert size(grad1d) == (q, p) 89*44554ea0SWill Pazner @assert length(qref1d) == q 90*44554ea0SWill Pazner @assert length(qweight1d) == q 91*44554ea0SWill Pazner 92*44554ea0SWill Pazner # Convert from Julia matrices (column-major) to row-major format 93*44554ea0SWill Pazner interp1d_rowmajor = collect(interp1d') 94*44554ea0SWill Pazner grad1d_rowmajor = collect(grad1d') 95*44554ea0SWill Pazner 96*44554ea0SWill Pazner ref = Ref{C.CeedBasis}() 97*44554ea0SWill Pazner C.CeedBasisCreateTensorH1( 98*44554ea0SWill Pazner c[], 99*44554ea0SWill Pazner dim, 100*44554ea0SWill Pazner ncomp, 101*44554ea0SWill Pazner p, 102*44554ea0SWill Pazner q, 103*44554ea0SWill Pazner interp1d_rowmajor, 104*44554ea0SWill Pazner grad1d_rowmajor, 105*44554ea0SWill Pazner qref1d, 106*44554ea0SWill Pazner qweight1d, 107*44554ea0SWill Pazner ref, 108*44554ea0SWill Pazner ) 109*44554ea0SWill Pazner Basis(ref) 110*44554ea0SWill Paznerend 111*44554ea0SWill Pazner 112*44554ea0SWill Pazner@doc raw""" 113*44554ea0SWill Pazner create_h1_basis(c::Ceed, topo::Topology, ncomp, nnodes, nqpts, interp, grad, qref, qweight) 114*44554ea0SWill Pazner 115*44554ea0SWill PaznerCreate a non tensor-product basis for H^1 discretizations 116*44554ea0SWill Pazner 117*44554ea0SWill Pazner# Arguments: 118*44554ea0SWill Pazner- `ceed`: A [`Ceed`](@ref) object where the [`Basis`](@ref) will be created. 119*44554ea0SWill Pazner- `topo`: [`Topology`](@ref) of element, e.g. hypercube, simplex, etc. 120*44554ea0SWill Pazner- `ncomp`: Number of field components (1 for scalar fields). 121*44554ea0SWill Pazner- `nnodes`: Total number of nodes. 122*44554ea0SWill Pazner- `nqpts`: Total number of quadrature points. 123*44554ea0SWill Pazner- `interp`: Matrix of size `(nqpts, nnodes)` expressing the values of nodal basis functions 124*44554ea0SWill Pazner at quadrature points. 125*44554ea0SWill Pazner- `grad`: Array of size `(dim, nqpts, nnodes)` expressing derivatives of nodal basis 126*44554ea0SWill Pazner functions at quadrature points. 127*44554ea0SWill Pazner- `qref`: Array of length `nqpts` holding the locations of quadrature points on the 128*44554ea0SWill Pazner reference element $[-1, 1]$. 129*44554ea0SWill Pazner- `qweight`: Array of length `nqpts` holding the quadrature weights on the reference 130*44554ea0SWill Pazner element. 131*44554ea0SWill Pazner""" 132*44554ea0SWill Paznerfunction create_h1_basis( 133*44554ea0SWill Pazner c::Ceed, 134*44554ea0SWill Pazner topo::Topology, 135*44554ea0SWill Pazner ncomp, 136*44554ea0SWill Pazner nnodes, 137*44554ea0SWill Pazner nqpts, 138*44554ea0SWill Pazner interp, 139*44554ea0SWill Pazner grad, 140*44554ea0SWill Pazner qref, 141*44554ea0SWill Pazner qweight, 142*44554ea0SWill Pazner) 143*44554ea0SWill Pazner dim = getdimension(topo) 144*44554ea0SWill Pazner @assert size(interp) == (nqpts, nnodes) 145*44554ea0SWill Pazner @assert size(grad) == (dim, nqpts, nnodes) 146*44554ea0SWill Pazner @assert length(qref) == nqpts 147*44554ea0SWill Pazner @assert length(qweight) == nqpts 148*44554ea0SWill Pazner 149*44554ea0SWill Pazner # Convert from Julia matrices and tensors (column-major) to row-major format 150*44554ea0SWill Pazner interp_rowmajor = collect(interp') 151*44554ea0SWill Pazner grad_rowmajor = permutedims(grad, [3, 2, 1]) 152*44554ea0SWill Pazner 153*44554ea0SWill Pazner ref = Ref{C.CeedBasis}() 154*44554ea0SWill Pazner C.CeedBasisCreateH1( 155*44554ea0SWill Pazner c[], 156*44554ea0SWill Pazner topo, 157*44554ea0SWill Pazner ncomp, 158*44554ea0SWill Pazner nnodes, 159*44554ea0SWill Pazner nqpts, 160*44554ea0SWill Pazner interp_rowmajor, 161*44554ea0SWill Pazner grad_rowmajor, 162*44554ea0SWill Pazner qref, 163*44554ea0SWill Pazner qweight, 164*44554ea0SWill Pazner ref, 165*44554ea0SWill Pazner ) 166*44554ea0SWill Pazner Basis(ref) 167*44554ea0SWill Paznerend 168*44554ea0SWill Pazner 169*44554ea0SWill Pazner""" 170*44554ea0SWill Pazner apply!(b::Basis, nelem, tmode::TransposeMode, emode::EvalMode, u::AbstractCeedVector, v::AbstractCeedVector) 171*44554ea0SWill Pazner 172*44554ea0SWill PaznerApply basis evaluation from nodes to quadrature points or vice versa, storing the result in 173*44554ea0SWill Paznerthe [`CeedVector`](@ref) `v`. 174*44554ea0SWill Pazner 175*44554ea0SWill Pazner`nelem` specifies the number of elements to apply the basis evaluation to; the backend will 176*44554ea0SWill Paznerspecify the ordering in CeedElemRestrictionCreateBlocked() 177*44554ea0SWill Pazner 178*44554ea0SWill PaznerSet `tmode` to `CEED_NOTRANSPOSE` to evaluate from nodes to quadrature or to 179*44554ea0SWill Pazner`CEED_TRANSPOSE` to apply the transpose, mapping from quadrature points to nodes. 180*44554ea0SWill Pazner 181*44554ea0SWill PaznerSet the [`EvalMode`](@ref) `emode` to: 182*44554ea0SWill Pazner- `CEED_EVAL_NONE` to use values directly, 183*44554ea0SWill Pazner- `CEED_EVAL_INTERP` to use interpolated values, 184*44554ea0SWill Pazner- `CEED_EVAL_GRAD` to use gradients, 185*44554ea0SWill Pazner- `CEED_EVAL_WEIGHT` to use quadrature weights. 186*44554ea0SWill Pazner""" 187*44554ea0SWill Paznerfunction apply!( 188*44554ea0SWill Pazner b::Basis, 189*44554ea0SWill Pazner nelem, 190*44554ea0SWill Pazner tmode::TransposeMode, 191*44554ea0SWill Pazner emode::EvalMode, 192*44554ea0SWill Pazner u::AbstractCeedVector, 193*44554ea0SWill Pazner v::AbstractCeedVector, 194*44554ea0SWill Pazner) 195*44554ea0SWill Pazner C.CeedBasisApply(b[], nelem, tmode, emode, u[], v[]) 196*44554ea0SWill Paznerend 197*44554ea0SWill Pazner 198*44554ea0SWill Pazner""" 199*44554ea0SWill Pazner apply(b::Basis, u::AbstractVector; nelem=1, tmode=NOTRANSPOSE, emode=EVAL_INTERP) 200*44554ea0SWill Pazner 201*44554ea0SWill PaznerPerforms the same function as the above-defined [`apply!`](@ref apply!(b::Basis, nelem, 202*44554ea0SWill Paznertmode::TransposeMode, emode::EvalMode, u::AbstractCeedVector, v::AbstractCeedVector)), but 203*44554ea0SWill Paznerautomatically convert from Julia arrays to [`CeedVector`](@ref) for convenience. 204*44554ea0SWill Pazner 205*44554ea0SWill PaznerThe result will be returned in a newly allocated array of the correct size. 206*44554ea0SWill Pazner""" 207*44554ea0SWill Paznerfunction apply(b::Basis, u::AbstractVector; nelem=1, tmode=NOTRANSPOSE, emode=EVAL_INTERP) 208*44554ea0SWill Pazner ceed_ref = Ref{C.Ceed}() 209*44554ea0SWill Pazner ccall((:CeedBasisGetCeed, C.libceed), Cint, (C.CeedBasis, Ptr{C.Ceed}), b[], ceed_ref) 210*44554ea0SWill Pazner c = Ceed(ceed_ref) 211*44554ea0SWill Pazner 212*44554ea0SWill Pazner u_vec = CeedVector(c, u) 213*44554ea0SWill Pazner 214*44554ea0SWill Pazner len_v = (tmode == TRANSPOSE) ? getnumnodes(b) : getnumqpts(b) 215*44554ea0SWill Pazner if emode == EVAL_GRAD 216*44554ea0SWill Pazner len_v *= getdimension(b) 217*44554ea0SWill Pazner end 218*44554ea0SWill Pazner 219*44554ea0SWill Pazner v_vec = CeedVector(c, len_v) 220*44554ea0SWill Pazner 221*44554ea0SWill Pazner apply!(b, nelem, tmode, emode, u_vec, v_vec) 222*44554ea0SWill Pazner Vector(v_vec) 223*44554ea0SWill Paznerend 224*44554ea0SWill Pazner 225*44554ea0SWill Pazner""" 226*44554ea0SWill Pazner getdimension(b::Basis) 227*44554ea0SWill Pazner 228*44554ea0SWill PaznerReturn the spatial dimension of the given [`Basis`](@ref). 229*44554ea0SWill Pazner""" 230*44554ea0SWill Paznerfunction getdimension(b::Basis) 231*44554ea0SWill Pazner dim = Ref{CeedInt}() 232*44554ea0SWill Pazner C.CeedBasisGetDimension(b[], dim) 233*44554ea0SWill Pazner dim[] 234*44554ea0SWill Paznerend 235*44554ea0SWill Pazner 236*44554ea0SWill Pazner""" 237*44554ea0SWill Pazner getdimension(t::Topology) 238*44554ea0SWill Pazner 239*44554ea0SWill PaznerReturn the spatial dimension of the given [`Topology`](@ref). 240*44554ea0SWill Pazner""" 241*44554ea0SWill Paznerfunction getdimension(t::Topology) 242*44554ea0SWill Pazner return Int(t) >> 16 243*44554ea0SWill Paznerend 244*44554ea0SWill Pazner 245*44554ea0SWill Pazner""" 246*44554ea0SWill Pazner gettopology(b::Basis) 247*44554ea0SWill Pazner 248*44554ea0SWill PaznerReturn the [`Topology`](@ref) of the given [`Basis`](@ref). 249*44554ea0SWill Pazner""" 250*44554ea0SWill Paznerfunction gettopology(b::Basis) 251*44554ea0SWill Pazner topo = Ref{Topology}() 252*44554ea0SWill Pazner C.CeedBasisGetTopology(b[], topo) 253*44554ea0SWill Pazner topo[] 254*44554ea0SWill Paznerend 255*44554ea0SWill Pazner 256*44554ea0SWill Pazner""" 257*44554ea0SWill Pazner getnumcomponents(b::Basis) 258*44554ea0SWill Pazner 259*44554ea0SWill PaznerReturn the number of components of the given [`Basis`](@ref). 260*44554ea0SWill Pazner""" 261*44554ea0SWill Paznerfunction getnumcomponents(b::Basis) 262*44554ea0SWill Pazner ncomp = Ref{CeedInt}() 263*44554ea0SWill Pazner C.CeedBasisGetNumComponents(b[], ncomp) 264*44554ea0SWill Pazner ncomp[] 265*44554ea0SWill Paznerend 266*44554ea0SWill Pazner 267*44554ea0SWill Pazner""" 268*44554ea0SWill Pazner getnumnodes(b::Basis) 269*44554ea0SWill Pazner 270*44554ea0SWill PaznerReturn the number of nodes of the given [`Basis`](@ref). 271*44554ea0SWill Pazner""" 272*44554ea0SWill Paznerfunction getnumnodes(b::Basis) 273*44554ea0SWill Pazner nnodes = Ref{CeedInt}() 274*44554ea0SWill Pazner C.CeedBasisGetNumNodes(b[], nnodes) 275*44554ea0SWill Pazner nnodes[] 276*44554ea0SWill Paznerend 277*44554ea0SWill Pazner 278*44554ea0SWill Pazner""" 279*44554ea0SWill Pazner getnumnodes1d(b::Basis) 280*44554ea0SWill Pazner 281*44554ea0SWill Pazner Return the number of 1D nodes of the given (tensor-product) [`Basis`](@ref). 282*44554ea0SWill Pazner""" 283*44554ea0SWill Paznerfunction getnumnodes1d(b::Basis) 284*44554ea0SWill Pazner nnodes1d = Ref{CeedInt}() 285*44554ea0SWill Pazner C.CeedBasisGetNumNodes1D(b[], nnodes1d) 286*44554ea0SWill Pazner nnodes1d[] 287*44554ea0SWill Paznerend 288*44554ea0SWill Pazner 289*44554ea0SWill Pazner""" 290*44554ea0SWill Pazner getnumqpts(b::Basis) 291*44554ea0SWill Pazner 292*44554ea0SWill PaznerReturn the number of quadrature points of the given [`Basis`](@ref). 293*44554ea0SWill Pazner""" 294*44554ea0SWill Paznerfunction getnumqpts(b::Basis) 295*44554ea0SWill Pazner nqpts = Ref{CeedInt}() 296*44554ea0SWill Pazner C.CeedBasisGetNumQuadraturePoints(b[], nqpts) 297*44554ea0SWill Pazner nqpts[] 298*44554ea0SWill Paznerend 299*44554ea0SWill Pazner 300*44554ea0SWill Pazner""" 301*44554ea0SWill Pazner getnumqpts1d(b::Basis) 302*44554ea0SWill Pazner 303*44554ea0SWill PaznerReturn the number of 1D quadrature points of the given (tensor-product) [`Basis`](@ref). 304*44554ea0SWill Pazner""" 305*44554ea0SWill Paznerfunction getnumqpts1d(b::Basis) 306*44554ea0SWill Pazner nqpts1d = Ref{CeedInt}() 307*44554ea0SWill Pazner C.CeedBasisGetNumQuadraturePoints1D(b[], nqpts1d) 308*44554ea0SWill Pazner nqpts1d[] 309*44554ea0SWill Paznerend 310*44554ea0SWill Pazner 311*44554ea0SWill Pazner""" 312*44554ea0SWill Pazner getqref(b::Basis) 313*44554ea0SWill Pazner 314*44554ea0SWill PaznerGet the reference coordinates of quadrature points (in `dim` dimensions) of the given 315*44554ea0SWill Pazner[`Basis`](@ref). 316*44554ea0SWill Pazner""" 317*44554ea0SWill Paznerfunction getqref(b::Basis) 318*44554ea0SWill Pazner ref = Ref{Ptr{CeedScalar}}() 319*44554ea0SWill Pazner C.CeedBasisGetQRef(b[], ref) 320*44554ea0SWill Pazner copy(unsafe_wrap(Array, ref[], getnumqpts(b))) 321*44554ea0SWill Paznerend 322*44554ea0SWill Pazner 323*44554ea0SWill Pazner""" 324*44554ea0SWill Pazner getqref(b::Basis) 325*44554ea0SWill Pazner 326*44554ea0SWill PaznerGet the quadrature weights of quadrature points (in `dim` dimensions) of the given 327*44554ea0SWill Pazner[`Basis`](@ref). 328*44554ea0SWill Pazner""" 329*44554ea0SWill Paznerfunction getqweights(b::Basis) 330*44554ea0SWill Pazner ref = Ref{Ptr{CeedScalar}}() 331*44554ea0SWill Pazner C.CeedBasisGetQWeights(b[], ref) 332*44554ea0SWill Pazner copy(unsafe_wrap(Array, ref[], getnumqpts(b))) 333*44554ea0SWill Paznerend 334*44554ea0SWill Pazner 335*44554ea0SWill Pazner""" 336*44554ea0SWill Pazner getinterp(b::Basis) 337*44554ea0SWill Pazner 338*44554ea0SWill PaznerGet the interpolation matrix of the given [`Basis`](@ref). Returns a matrix of size 339*44554ea0SWill Pazner`(getnumqpts(b), getnumnodes(b))`. 340*44554ea0SWill Pazner""" 341*44554ea0SWill Paznerfunction getinterp(b::Basis) 342*44554ea0SWill Pazner ref = Ref{Ptr{CeedScalar}}() 343*44554ea0SWill Pazner C.CeedBasisGetInterp(b[], ref) 344*44554ea0SWill Pazner q = getnumqpts(b) 345*44554ea0SWill Pazner p = getnumnodes(b) 346*44554ea0SWill Pazner collect(unsafe_wrap(Array, ref[], (p, q))') 347*44554ea0SWill Paznerend 348*44554ea0SWill Pazner 349*44554ea0SWill Pazner""" 350*44554ea0SWill Pazner getinterp1d(b::Basis) 351*44554ea0SWill Pazner 352*44554ea0SWill PaznerGet the 1D interpolation matrix of the given [`Basis`](@ref). `b` must be a tensor-product 353*44554ea0SWill Paznerbasis, otherwise this function will fail. Returns a matrix of size `(getnumqpts1d(b), 354*44554ea0SWill Paznergetnumnodes1d(b))`. 355*44554ea0SWill Pazner""" 356*44554ea0SWill Paznerfunction getinterp1d(b::Basis) 357*44554ea0SWill Pazner ref = Ref{Ptr{CeedScalar}}() 358*44554ea0SWill Pazner C.CeedBasisGetInterp1D(b[], ref) 359*44554ea0SWill Pazner q = getnumqpts1d(b) 360*44554ea0SWill Pazner p = getnumnodes1d(b) 361*44554ea0SWill Pazner collect(unsafe_wrap(Array, ref[], (p, q))') 362*44554ea0SWill Paznerend 363*44554ea0SWill Pazner 364*44554ea0SWill Pazner""" 365*44554ea0SWill Pazner getgad(b::Basis) 366*44554ea0SWill Pazner 367*44554ea0SWill PaznerGet the gradient matrix of the given [`Basis`](@ref). Returns a tensor of size 368*44554ea0SWill Pazner`(getdimension(b), getnumqpts(b), getnumnodes(b))`. 369*44554ea0SWill Pazner""" 370*44554ea0SWill Paznerfunction getgrad(b::Basis) 371*44554ea0SWill Pazner ref = Ref{Ptr{CeedScalar}}() 372*44554ea0SWill Pazner C.CeedBasisGetGrad(b[], ref) 373*44554ea0SWill Pazner dim = getdimension(b) 374*44554ea0SWill Pazner q = getnumqpts(b) 375*44554ea0SWill Pazner p = getnumnodes(b) 376*44554ea0SWill Pazner permutedims(unsafe_wrap(Array, ref[], (p, q, dim)), [3, 2, 1]) 377*44554ea0SWill Paznerend 378*44554ea0SWill Pazner 379*44554ea0SWill Pazner""" 380*44554ea0SWill Pazner getgrad1d(b::Basis) 381*44554ea0SWill Pazner 382*44554ea0SWill PaznerGet the 1D derivative matrix of the given [`Basis`](@ref). Returns a matrix of size 383*44554ea0SWill Pazner`(getnumqpts(b), getnumnodes(b))`. 384*44554ea0SWill Pazner""" 385*44554ea0SWill Paznerfunction getgrad1d(b::Basis) 386*44554ea0SWill Pazner ref = Ref{Ptr{CeedScalar}}() 387*44554ea0SWill Pazner C.CeedBasisGetGrad1D(b[], ref) 388*44554ea0SWill Pazner q = getnumqpts1d(b) 389*44554ea0SWill Pazner p = getnumnodes1d(b) 390*44554ea0SWill Pazner collect(unsafe_wrap(Array, ref[], (p, q))') 391*44554ea0SWill Paznerend 392