Lines Matching +full:format +full:- +full:c
9 Base.getindex(::BasisNone) = C.CEED_BASIS_NONE[]
17 - [`create_tensor_h1_lagrange_basis`](@ref)
18 - [`create_tensor_h1_basis`](@ref)
19 - [`create_h1_basis`](@ref)
20 - [`create_hdiv_basis`](@ref)
21 - [`create_hcurl_basis`](@ref)
24 ref::RefValue{C.CeedBasis}
34 destroy(b::Basis) = C.CeedBasisDestroy(b.ref) # COV_EXCL_LINE
36 Base.show(io::IO, ::MIME"text/plain", b::Basis) = ceed_show(io, b, C.CeedBasisView)
41 Create a tensor-product Lagrange basis.
44 - `ceed`: A [`Ceed`](@ref) object where the [`Basis`](@ref) will be created.
45 - `dim`: Topological dimension of element.
46 - `ncomp`: Number of field components (1 for scalar fields).
47 - `p`: Number of Gauss-Lobatto nodes in one dimension. The polynomial degree of the
48 resulting $Q_k$ element is $k=p-1$.
49 - `q`: Number of quadrature points in one dimension.
50 - `qmode`: Distribution of the $q$ quadrature points (affects order of accuracy for the
53 function create_tensor_h1_lagrange_basis(c::Ceed, dim, ncomp, p, q, quad_mode::QuadMode)
54 ref = Ref{C.CeedBasis}()
55 C.CeedBasisCreateTensorH1Lagrange(c[], dim, ncomp, p, q, quad_mode, ref)
60 create_tensor_h1_basis(c::Ceed, dim, ncomp, p, q, interp1d, grad1d, qref1d, qweight1d)
62 Create a tensor-product basis for $H^1$ discretizations.
65 - `ceed`: A [`Ceed`](@ref) object where the [`Basis`](@ref) will be created.
66 - `dim`: Topological dimension.
67 - `ncomp`: Number of field components (1 for scalar fields).
68 - `p`: Number of nodes in one dimension.
69 - `q`: Number of quadrature points in one dimension
70 - `interp1d`: Matrix of size `(q, p)` expressing the values of nodal basis functions at
72 - `grad1d`: Matrix of size `(p, q)` expressing derivatives of nodal basis functions at
74 - `qref1d`: Array of length `q` holding the locations of quadrature points on the 1D
75 reference element $[-1, 1]$.
76 - `qweight1d`: Array of length `q` holding the quadrature weights on the reference element.
79 c::Ceed,
94 # Convert from Julia matrices (column-major) to row-major format
98 ref = Ref{C.CeedBasis}()
99 C.CeedBasisCreateTensorH1(
100 c[],
115 create_h1_basis(c::Ceed, topo::Topology, ncomp, nnodes, nqpts, interp, grad, qref, qweight)
117 Create a non tensor-product basis for $H^1$ discretizations
120 - `ceed`: A [`Ceed`](@ref) object where the [`Basis`](@ref) will be created.
121 - `topo`: [`Topology`](@ref) of element, e.g. hypercube, simplex, etc.
122 - `ncomp`: Number of field components (1 for scalar fields).
123 - `nnodes`: Total number of nodes.
124 - `nqpts`: Total number of quadrature points.
125 - `interp`: Matrix of size `(nqpts, nnodes)` expressing the values of nodal basis functions
127 - `grad`: Array of size `(dim, nqpts, nnodes)` expressing derivatives of nodal basis
129 - `qref`: Matrix of size `(dim, nqpts)` holding the locations of quadrature points on the
130 reference element $[-1, 1]$.
131 - `qweight`: Array of length `nqpts` holding the quadrature weights on the reference
135 c::Ceed,
151 # Convert from Julia matrices and tensors (column-major) to row-major format
156 ref = Ref{C.CeedBasis}()
157 C.CeedBasisCreateH1(
158 c[],
173 create_hdiv_basis(c::Ceed, topo::Topology, ncomp, nnodes, nqpts, interp, div, qref, qweight)
175 Create a non tensor-product basis for H(div) discretizations
178 - `ceed`: A [`Ceed`](@ref) object where the [`Basis`](@ref) will be created.
179 - `topo`: [`Topology`](@ref) of element, e.g. hypercube, simplex, etc.
180 - `ncomp`: Number of field components (1 for scalar fields).
181 - `nnodes`: Total number of nodes.
182 - `nqpts`: Total number of quadrature points.
183 - `interp`: Matrix of size `(dim, nqpts, nnodes)` expressing the values of basis functions
185 - `div`: Array of size `(nqpts, nnodes)` expressing divergence of basis functions at
187 - `qref`: Matrix of size `(dim, nqpts)` holding the locations of quadrature points on the
188 reference element $[-1, 1]$.
189 - `qweight`: Array of length `nqpts` holding the quadrature weights on the reference
193 c::Ceed,
209 # Convert from Julia matrices and tensors (column-major) to row-major format
214 ref = Ref{C.CeedBasis}()
215 C.CeedBasisCreateHdiv(
216 c[],
231 create_hcurl_basis(c::Ceed, topo::Topology, ncomp, nnodes, nqpts, interp, curl, qref, qweight)
233 Create a non tensor-product basis for H(curl) discretizations
236 - `ceed`: A [`Ceed`](@ref) object where the [`Basis`](@ref) will be created.
237 - `topo`: [`Topology`](@ref) of element, e.g. hypercube, simplex, etc.
238 - `ncomp`: Number of field components (1 for scalar fields).
239 - `nnodes`: Total number of nodes.
240 - `nqpts`: Total number of quadrature points.
241 - `interp`: Matrix of size `(dim, nqpts, nnodes)` expressing the values of basis functions
243 - `curl`: Matrix of size `(curlcomp, nqpts, nnodes)`, `curlcomp = 1 if dim < 3 else dim`)
245 - `qref`: Matrix of size `(dim, nqpts)` holding the locations of quadrature points on the
246 reference element $[-1, 1]$.
247 - `qweight`: Array of length `nqpts` holding the quadrature weights on the reference
251 c::Ceed,
268 # Convert from Julia matrices and tensors (column-major) to row-major format
273 ref = Ref{C.CeedBasis}()
274 C.CeedBasisCreateHcurl(
275 c[],
302 - `CEED_EVAL_NONE` to use values directly,
303 - `CEED_EVAL_INTERP` to use interpolated values,
304 - `CEED_EVAL_GRAD` to use gradients,
305 - `CEED_EVAL_WEIGHT` to use quadrature weights.
315 C.CeedBasisApply(b[], nelem, tmode, emode, u[], v[])
321 Performs the same function as the above-defined [`apply!`](@ref apply!(b::Basis, nelem,
328 ceed_ref = Ref{C.Ceed}()
329 ccall((:CeedBasisGetCeed, C.libceed), Cint, (C.CeedBasis, Ptr{C.Ceed}), b[], ceed_ref)
330 c = Ceed(ceed_ref)
332 u_vec = CeedVector(c, u)
335 C.CeedBasisGetNumQuadratureComponents(b[], emode, qcomp)
338 v_vec = CeedVector(c, len_v)
351 C.CeedBasisGetDimension(b[], dim)
371 C.CeedBasisGetTopology(b[], topo)
382 C.CeedBasisGetNumComponents(b[], ncomp)
393 C.CeedBasisGetNumNodes(b[], nnodes)
400 Return the number of 1D nodes of the given (tensor-product) [`Basis`](@ref).
404 C.CeedBasisGetNumNodes1D(b[], nnodes1d)
415 C.CeedBasisGetNumQuadraturePoints(b[], nqpts)
422 Return the number of 1D quadrature points of the given (tensor-product) [`Basis`](@ref).
426 C.CeedBasisGetNumQuadraturePoints1D(b[], nqpts1d)
438 C.CeedBasisIsTensor(b[], istensor)
440 C.CeedBasisGetQRef(b[], ref)
455 C.CeedBasisIsTensor(b[], istensor)
457 C.CeedBasisGetQWeights(b[], ref)
470 C.CeedBasisGetInterp(b[], ref)
474 C.CeedBasisGetNumQuadratureComponents(b[], C.CEED_EVAL_INTERP, qcomp)
485 Get the 1D interpolation matrix of the given [`Basis`](@ref). `b` must be a tensor-product
491 C.CeedBasisGetInterp1D(b[], ref)
505 C.CeedBasisGetGrad(b[], ref)
520 C.CeedBasisGetGrad1D(b[], ref)
534 C.CeedBasisGetDiv(b[], ref)
549 C.CeedBasisGetCurl(b[], ref)
553 C.CeedBasisGetNumQuadratureComponents(b[], C.CEED_EVAL_CURL, qcomp)