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