-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathBasis.html
4 lines (3 loc) · 27.8 KB
/
Basis.html
1
2
3
4
<!DOCTYPE html>
<html lang="en"><head><meta charset="UTF-8"/><meta name="viewport" content="width=device-width, initial-scale=1.0"/><title>Basis · LibCEED.jl Docs</title><link href="https://fonts.googleapis.com/css?family=Lato|Roboto+Mono" rel="stylesheet" type="text/css"/><link href="https://cdnjs.cloudflare.com/ajax/libs/font-awesome/5.15.0/css/fontawesome.min.css" rel="stylesheet" type="text/css"/><link href="https://cdnjs.cloudflare.com/ajax/libs/font-awesome/5.15.0/css/solid.min.css" rel="stylesheet" type="text/css"/><link href="https://cdnjs.cloudflare.com/ajax/libs/font-awesome/5.15.0/css/brands.min.css" rel="stylesheet" type="text/css"/><link href="https://cdnjs.cloudflare.com/ajax/libs/KaTeX/0.11.1/katex.min.css" rel="stylesheet" type="text/css"/><script>documenterBaseURL="."</script><script src="https://cdnjs.cloudflare.com/ajax/libs/require.js/2.3.6/require.min.js" data-main="assets/documenter.js"></script><script src="siteinfo.js"></script><script src="../versions.js"></script><link class="docs-theme-link" rel="stylesheet" type="text/css" href="assets/themes/documenter-dark.css" data-theme-name="documenter-dark"/><link class="docs-theme-link" rel="stylesheet" type="text/css" href="assets/themes/documenter-light.css" data-theme-name="documenter-light" data-theme-primary/><script src="assets/themeswap.js"></script></head><body><div id="documenter"><nav class="docs-sidebar"><div class="docs-package-name"><span class="docs-autofit">LibCEED.jl Docs</span></div><form class="docs-search" action="search.html"><input class="docs-search-query" id="documenter-search-query" name="q" type="text" placeholder="Search docs"/></form><ul class="docs-menu"><li><a class="tocitem" href="index.html">Home</a></li><li><span class="tocitem">Ceed Objects</span><ul><li><a class="tocitem" href="Ceed.html">Ceed</a></li><li><a class="tocitem" href="CeedVector.html">CeedVector</a></li><li><a class="tocitem" href="ElemRestriction.html">ElemRestriction</a></li><li class="is-active"><a class="tocitem" href="Basis.html">Basis</a></li><li><a class="tocitem" href="QFunction.html">QFunction</a></li><li><a class="tocitem" href="Operator.html">Operator</a></li></ul></li><li><span class="tocitem">Utilities</span><ul><li><a class="tocitem" href="Misc.html">Linear Algebra</a></li><li><a class="tocitem" href="Globals.html">Constants and Enumerations</a></li><li><a class="tocitem" href="Quadrature.html">Quadrature</a></li></ul></li><li><a class="tocitem" href="LibCEED.html">Library configuration</a></li><li><a class="tocitem" href="C.html">Low-level C interface</a></li><li><a class="tocitem" href="UserQFunctions.html">Defining User Q-Functions</a></li><li><a class="tocitem" href="Examples.html">Examples</a></li></ul><div class="docs-version-selector field has-addons"><div class="control"><span class="docs-label button is-static is-size-7">Version</span></div><div class="docs-selector control is-expanded"><div class="select is-fullwidth is-size-7"><select id="documenter-version-selector"></select></div></div></div></nav><div class="docs-main"><header class="docs-navbar"><nav class="breadcrumb"><ul class="is-hidden-mobile"><li><a class="is-disabled">Ceed Objects</a></li><li class="is-active"><a href="Basis.html">Basis</a></li></ul><ul class="is-hidden-tablet"><li class="is-active"><a href="Basis.html">Basis</a></li></ul></nav><div class="docs-right"><a class="docs-edit-link" href="https://github.com/CEED/libCEED/blob/master/julia/LibCEED.jl/docs/src/Basis.md" title="Edit on GitHub"><span class="docs-icon fab"></span><span class="docs-label is-hidden-touch">Edit on GitHub</span></a><a class="docs-settings-button fas fa-cog" id="documenter-settings-button" href="#" title="Settings"></a><a class="docs-sidebar-button fa fa-bars is-hidden-desktop" id="documenter-sidebar-button" href="#"></a></div></header><article class="content" id="documenter-page"><h1 id="Basis"><a class="docs-heading-anchor" href="#Basis">Basis</a><a id="Basis-1"></a><a class="docs-heading-anchor-permalink" href="#Basis" title="Permalink"></a></h1><div class="admonition is-info"><header class="admonition-header">Column-major vs. row-major storage</header><div class="admonition-body"><p>libCEED internally uses row-major (C convention) storage of matrices, while Julia uses column-major (Fortran convention) storage.</p><p>LibCEED.jl will typically handle the conversion between these formats by transposing or permuting the dimensions of the input and output matrices and tensors.</p></div></div><article class="docstring"><header><a class="docstring-binding" id="LibCEED.Basis" href="#LibCEED.Basis"><code>LibCEED.Basis</code></a> — <span class="docstring-category">Type</span></header><section><div><pre><code class="language-julia">Basis</code></pre><p>Wraps a <code>CeedBasis</code> object, representing a finite element basis. A <code>Basis</code> object can be created using one of:</p><ul><li><a href="Basis.html#LibCEED.create_tensor_h1_lagrange_basis"><code>create_tensor_h1_lagrange_basis</code></a></li><li><a href="Basis.html#LibCEED.create_tensor_h1_basis"><code>create_tensor_h1_basis</code></a></li><li><a href="Basis.html#LibCEED.create_h1_basis"><code>create_h1_basis</code></a></li><li><a href="Basis.html#LibCEED.create_hdiv_basis"><code>create_hdiv_basis</code></a></li><li><a href="Basis.html#LibCEED.create_hcurl_basis"><code>create_hcurl_basis</code></a></li></ul></div><a class="docs-sourcelink" target="_blank" href="https://github.com/CEED/libCEED/blob/4018a20a98d451fac24765d3ddb936861647ce8d/julia/LibCEED.jl/src/Basis.jl#L11-L22">source</a></section></article><div class="admonition is-warning"><header class="admonition-header">Missing docstring.</header><div class="admonition-body"><p>Missing docstring for <code>BasisCollocated</code>. Check Documenter's build log for details.</p></div></div><article class="docstring"><header><a class="docstring-binding" id="LibCEED.create_tensor_h1_lagrange_basis" href="#LibCEED.create_tensor_h1_lagrange_basis"><code>LibCEED.create_tensor_h1_lagrange_basis</code></a> — <span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">create_tensor_h1_lagrange_basis(ceed, dim, ncomp, p, q, qmode)</code></pre><p>Create a tensor-product Lagrange basis.</p><p><strong>Arguments:</strong></p><ul><li><code>ceed</code>: A <a href="Ceed.html#LibCEED.Ceed"><code>Ceed</code></a> object where the <a href="Basis.html#LibCEED.Basis"><code>Basis</code></a> will be created.</li><li><code>dim</code>: Topological dimension of element.</li><li><code>ncomp</code>: Number of field components (1 for scalar fields).</li><li><code>p</code>: Number of Gauss-Lobatto nodes in one dimension. The polynomial degree of the resulting <span>$Q_k$</span> element is <span>$k=p-1$</span>.</li><li><code>q</code>: Number of quadrature points in one dimension.</li><li><code>qmode</code>: Distribution of the <span>$q$</span> quadrature points (affects order of accuracy for the quadrature).</li></ul></div><a class="docs-sourcelink" target="_blank" href="https://github.com/CEED/libCEED/blob/4018a20a98d451fac24765d3ddb936861647ce8d/julia/LibCEED.jl/src/Basis.jl#L38-L52">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="LibCEED.create_tensor_h1_basis" href="#LibCEED.create_tensor_h1_basis"><code>LibCEED.create_tensor_h1_basis</code></a> — <span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">create_tensor_h1_basis(c::Ceed, dim, ncomp, p, q, interp1d, grad1d, qref1d, qweight1d)</code></pre><p>Create a tensor-product basis for <span>$H^1$</span> discretizations.</p><p><strong>Arguments:</strong></p><ul><li><code>ceed</code>: A <a href="Ceed.html#LibCEED.Ceed"><code>Ceed</code></a> object where the <a href="Basis.html#LibCEED.Basis"><code>Basis</code></a> will be created.</li><li><code>dim</code>: Topological dimension.</li><li><code>ncomp</code>: Number of field components (1 for scalar fields).</li><li><code>p</code>: Number of nodes in one dimension.</li><li><code>q</code>: Number of quadrature points in one dimension</li><li><code>interp1d</code>: Matrix of size <code>(q, p)</code> expressing the values of nodal basis functions at quadrature points.</li><li><code>grad1d</code>: Matrix of size <code>(p, q)</code> expressing derivatives of nodal basis functions at quadrature points.</li><li><code>qref1d</code>: Array of length <code>q</code> holding the locations of quadrature points on the 1D reference element <span>$[-1, 1]$</span>.</li><li><code>qweight1d</code>: Array of length <code>q</code> holding the quadrature weights on the reference element.</li></ul></div><a class="docs-sourcelink" target="_blank" href="https://github.com/CEED/libCEED/blob/4018a20a98d451fac24765d3ddb936861647ce8d/julia/LibCEED.jl/src/Basis.jl#L59-L77">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="LibCEED.create_h1_basis" href="#LibCEED.create_h1_basis"><code>LibCEED.create_h1_basis</code></a> — <span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">create_h1_basis(c::Ceed, topo::Topology, ncomp, nnodes, nqpts, interp, grad, qref, qweight)</code></pre><p>Create a non tensor-product basis for <span>$H^1$</span> discretizations</p><p><strong>Arguments:</strong></p><ul><li><code>ceed</code>: A <a href="Ceed.html#LibCEED.Ceed"><code>Ceed</code></a> object where the <a href="Basis.html#LibCEED.Basis"><code>Basis</code></a> will be created.</li><li><code>topo</code>: <a href="Globals.html#LibCEED.Topology"><code>Topology</code></a> of element, e.g. hypercube, simplex, etc.</li><li><code>ncomp</code>: Number of field components (1 for scalar fields).</li><li><code>nnodes</code>: Total number of nodes.</li><li><code>nqpts</code>: Total number of quadrature points.</li><li><code>interp</code>: Matrix of size <code>(nqpts, nnodes)</code> expressing the values of nodal basis functions at quadrature points.</li><li><code>grad</code>: Array of size <code>(dim, nqpts, nnodes)</code> expressing derivatives of nodal basis functions at quadrature points.</li><li><code>qref</code>: Matrix of size <code>(dim, nqpts)</code> holding the locations of quadrature points on the reference element <span>$[-1, 1]$</span>.</li><li><code>qweight</code>: Array of length <code>nqpts</code> holding the quadrature weights on the reference element.</li></ul></div><a class="docs-sourcelink" target="_blank" href="https://github.com/CEED/libCEED/blob/4018a20a98d451fac24765d3ddb936861647ce8d/julia/LibCEED.jl/src/Basis.jl#L114-L133">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="LibCEED.create_hdiv_basis" href="#LibCEED.create_hdiv_basis"><code>LibCEED.create_hdiv_basis</code></a> — <span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">create_hdiv_basis(c::Ceed, topo::Topology, ncomp, nnodes, nqpts, interp, div, qref, qweight)</code></pre><p>Create a non tensor-product basis for H(div) discretizations</p><p><strong>Arguments:</strong></p><ul><li><code>ceed</code>: A <a href="Ceed.html#LibCEED.Ceed"><code>Ceed</code></a> object where the <a href="Basis.html#LibCEED.Basis"><code>Basis</code></a> will be created.</li><li><code>topo</code>: <a href="Globals.html#LibCEED.Topology"><code>Topology</code></a> of element, e.g. hypercube, simplex, etc.</li><li><code>ncomp</code>: Number of field components (1 for scalar fields).</li><li><code>nnodes</code>: Total number of nodes.</li><li><code>nqpts</code>: Total number of quadrature points.</li><li><code>interp</code>: Matrix of size <code>(dim, nqpts, nnodes)</code> expressing the values of basis functions at quadrature points.</li><li><code>div</code>: Array of size <code>(nqpts, nnodes)</code> expressing divergence of basis functions at quadrature points.</li><li><code>qref</code>: Matrix of size <code>(dim, nqpts)</code> holding the locations of quadrature points on the reference element <span>$[-1, 1]$</span>.</li><li><code>qweight</code>: Array of length <code>nqpts</code> holding the quadrature weights on the reference element.</li></ul></div><a class="docs-sourcelink" target="_blank" href="https://github.com/CEED/libCEED/blob/4018a20a98d451fac24765d3ddb936861647ce8d/julia/LibCEED.jl/src/Basis.jl#L172-L191">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="LibCEED.create_hcurl_basis" href="#LibCEED.create_hcurl_basis"><code>LibCEED.create_hcurl_basis</code></a> — <span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">create_hcurl_basis(c::Ceed, topo::Topology, ncomp, nnodes, nqpts, interp, curl, qref, qweight)</code></pre><p>Create a non tensor-product basis for H(curl) discretizations</p><p><strong>Arguments:</strong></p><ul><li><code>ceed</code>: A <a href="Ceed.html#LibCEED.Ceed"><code>Ceed</code></a> object where the <a href="Basis.html#LibCEED.Basis"><code>Basis</code></a> will be created.</li><li><code>topo</code>: <a href="Globals.html#LibCEED.Topology"><code>Topology</code></a> of element, e.g. hypercube, simplex, etc.</li><li><code>ncomp</code>: Number of field components (1 for scalar fields).</li><li><code>nnodes</code>: Total number of nodes.</li><li><code>nqpts</code>: Total number of quadrature points.</li><li><code>interp</code>: Matrix of size <code>(dim, nqpts, nnodes)</code> expressing the values of basis functions at quadrature points.</li><li><code>curl</code>: Matrix of size <code>(curlcomp, nqpts, nnodes)</code>, <code>curlcomp = 1 if dim < 3 else dim</code>) matrix expressing curl of basis functions at quadrature points.</li><li><code>qref</code>: Matrix of size <code>(dim, nqpts)</code> holding the locations of quadrature points on the reference element <span>$[-1, 1]$</span>.</li><li><code>qweight</code>: Array of length <code>nqpts</code> holding the quadrature weights on the reference element.</li></ul></div><a class="docs-sourcelink" target="_blank" href="https://github.com/CEED/libCEED/blob/4018a20a98d451fac24765d3ddb936861647ce8d/julia/LibCEED.jl/src/Basis.jl#L230-L249">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="LibCEED.apply!-Tuple{Basis, Any, LibCEED.C.CeedTransposeMode, LibCEED.C.CeedEvalMode, LibCEED.AbstractCeedVector, LibCEED.AbstractCeedVector}" href="#LibCEED.apply!-Tuple{Basis, Any, LibCEED.C.CeedTransposeMode, LibCEED.C.CeedEvalMode, LibCEED.AbstractCeedVector, LibCEED.AbstractCeedVector}"><code>LibCEED.apply!</code></a> — <span class="docstring-category">Method</span></header><section><div><pre><code class="language-julia">apply!(b::Basis, nelem, tmode::TransposeMode, emode::EvalMode, u::AbstractCeedVector, v::AbstractCeedVector)</code></pre><p>Apply basis evaluation from nodes to quadrature points or vice versa, storing the result in the <a href="CeedVector.html#LibCEED.CeedVector"><code>CeedVector</code></a> <code>v</code>.</p><p><code>nelem</code> specifies the number of elements to apply the basis evaluation to; the backend will specify the ordering in CeedElemRestrictionCreateBlocked()</p><p>Set <code>tmode</code> to <code>CEED_NOTRANSPOSE</code> to evaluate from nodes to quadrature or to <code>CEED_TRANSPOSE</code> to apply the transpose, mapping from quadrature points to nodes.</p><p>Set the <a href="Globals.html#LibCEED.EvalMode"><code>EvalMode</code></a> <code>emode</code> to:</p><ul><li><code>CEED_EVAL_NONE</code> to use values directly,</li><li><code>CEED_EVAL_INTERP</code> to use interpolated values,</li><li><code>CEED_EVAL_GRAD</code> to use gradients,</li><li><code>CEED_EVAL_WEIGHT</code> to use quadrature weights.</li></ul></div><a class="docs-sourcelink" target="_blank" href="https://github.com/CEED/libCEED/blob/4018a20a98d451fac24765d3ddb936861647ce8d/julia/LibCEED.jl/src/Basis.jl#L289-L306">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="LibCEED.apply-Tuple{Basis, AbstractVector}" href="#LibCEED.apply-Tuple{Basis, AbstractVector}"><code>LibCEED.apply</code></a> — <span class="docstring-category">Method</span></header><section><div><pre><code class="language-julia">apply(b::Basis, u::AbstractVector; nelem=1, tmode=NOTRANSPOSE, emode=EVAL_INTERP)</code></pre><p>Performs the same function as the above-defined <a href="Basis.html#LibCEED.apply!-Tuple{Basis, Any, LibCEED.C.CeedTransposeMode, LibCEED.C.CeedEvalMode, LibCEED.AbstractCeedVector, LibCEED.AbstractCeedVector}"><code>apply!</code></a>, but automatically convert from Julia arrays to <a href="CeedVector.html#LibCEED.CeedVector"><code>CeedVector</code></a> for convenience.</p><p>The result will be returned in a newly allocated array of the correct size.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/CEED/libCEED/blob/4018a20a98d451fac24765d3ddb936861647ce8d/julia/LibCEED.jl/src/Basis.jl#L318-L326">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="LibCEED.getdimension" href="#LibCEED.getdimension"><code>LibCEED.getdimension</code></a> — <span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">getdimension(b::Basis)</code></pre><p>Return the spatial dimension of the given <a href="Basis.html#LibCEED.Basis"><code>Basis</code></a>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/CEED/libCEED/blob/4018a20a98d451fac24765d3ddb936861647ce8d/julia/LibCEED.jl/src/Basis.jl#L344-L348">source</a></section><section><div><pre><code class="language-none">getdimension(t::Topology)</code></pre><p>Return the spatial dimension of the given <a href="Globals.html#LibCEED.Topology"><code>Topology</code></a>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/CEED/libCEED/blob/4018a20a98d451fac24765d3ddb936861647ce8d/julia/LibCEED.jl/src/Basis.jl#L355-L359">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="LibCEED.gettopology" href="#LibCEED.gettopology"><code>LibCEED.gettopology</code></a> — <span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">gettopology(b::Basis)</code></pre><p>Return the <a href="Globals.html#LibCEED.Topology"><code>Topology</code></a> of the given <a href="Basis.html#LibCEED.Basis"><code>Basis</code></a>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/CEED/libCEED/blob/4018a20a98d451fac24765d3ddb936861647ce8d/julia/LibCEED.jl/src/Basis.jl#L364-L368">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="LibCEED.getnumcomponents" href="#LibCEED.getnumcomponents"><code>LibCEED.getnumcomponents</code></a> — <span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">getnumcomponents(b::Basis)</code></pre><p>Return the number of components of the given <a href="Basis.html#LibCEED.Basis"><code>Basis</code></a>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/CEED/libCEED/blob/4018a20a98d451fac24765d3ddb936861647ce8d/julia/LibCEED.jl/src/Basis.jl#L375-L379">source</a></section><section><div><pre><code class="language-julia">getnumcomponents(r::ElemRestriction)</code></pre><p>Get the number of components in the elements of an <a href="ElemRestriction.html#LibCEED.ElemRestriction"><code>ElemRestriction</code></a>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/CEED/libCEED/blob/4018a20a98d451fac24765d3ddb936861647ce8d/julia/LibCEED.jl/src/ElemRestriction.jl#L393-L397">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="LibCEED.getnumnodes" href="#LibCEED.getnumnodes"><code>LibCEED.getnumnodes</code></a> — <span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">getnumnodes(b::Basis)</code></pre><p>Return the number of nodes of the given <a href="Basis.html#LibCEED.Basis"><code>Basis</code></a>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/CEED/libCEED/blob/4018a20a98d451fac24765d3ddb936861647ce8d/julia/LibCEED.jl/src/Basis.jl#L386-L390">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="LibCEED.getnumnodes1d" href="#LibCEED.getnumnodes1d"><code>LibCEED.getnumnodes1d</code></a> — <span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">getnumnodes1d(b::Basis)
Return the number of 1D nodes of the given (tensor-product) [`Basis`](@ref).</code></pre></div><a class="docs-sourcelink" target="_blank" href="https://github.com/CEED/libCEED/blob/4018a20a98d451fac24765d3ddb936861647ce8d/julia/LibCEED.jl/src/Basis.jl#L397-L401">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="LibCEED.getnumqpts" href="#LibCEED.getnumqpts"><code>LibCEED.getnumqpts</code></a> — <span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">getnumqpts(b::Basis)</code></pre><p>Return the number of quadrature points of the given <a href="Basis.html#LibCEED.Basis"><code>Basis</code></a>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/CEED/libCEED/blob/4018a20a98d451fac24765d3ddb936861647ce8d/julia/LibCEED.jl/src/Basis.jl#L408-L412">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="LibCEED.getnumqpts1d" href="#LibCEED.getnumqpts1d"><code>LibCEED.getnumqpts1d</code></a> — <span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">getnumqpts1d(b::Basis)</code></pre><p>Return the number of 1D quadrature points of the given (tensor-product) <a href="Basis.html#LibCEED.Basis"><code>Basis</code></a>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/CEED/libCEED/blob/4018a20a98d451fac24765d3ddb936861647ce8d/julia/LibCEED.jl/src/Basis.jl#L419-L423">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="LibCEED.getqref" href="#LibCEED.getqref"><code>LibCEED.getqref</code></a> — <span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">getqref(b::Basis)</code></pre><p>Get the reference coordinates of quadrature points (in <code>dim</code> dimensions) of the given <a href="Basis.html#LibCEED.Basis"><code>Basis</code></a>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/CEED/libCEED/blob/4018a20a98d451fac24765d3ddb936861647ce8d/julia/LibCEED.jl/src/Basis.jl#L430-L435">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="LibCEED.getqweights" href="#LibCEED.getqweights"><code>LibCEED.getqweights</code></a> — <span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">getqref(b::Basis)</code></pre><p>Get the quadrature weights of quadrature points (in <code>dim</code> dimensions) of the given <a href="Basis.html#LibCEED.Basis"><code>Basis</code></a>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/CEED/libCEED/blob/4018a20a98d451fac24765d3ddb936861647ce8d/julia/LibCEED.jl/src/Basis.jl#L447-L452">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="LibCEED.getinterp" href="#LibCEED.getinterp"><code>LibCEED.getinterp</code></a> — <span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">getinterp(b::Basis)</code></pre><p>Get the interpolation matrix of the given <a href="Basis.html#LibCEED.Basis"><code>Basis</code></a>. Returns a matrix of size <code>(getnumqpts(b), getnumnodes(b))</code> for a given <span>$H^1$</span> basis or <code>(getdimension(b), getnumqpts(b), getnumnodes(b))</code> for a given vector <span>$H(div)$</span> or <span>$H(curl)$</span> basis.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/CEED/libCEED/blob/4018a20a98d451fac24765d3ddb936861647ce8d/julia/LibCEED.jl/src/Basis.jl#L461-L467">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="LibCEED.getinterp1d" href="#LibCEED.getinterp1d"><code>LibCEED.getinterp1d</code></a> — <span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">getinterp1d(b::Basis)</code></pre><p>Get the 1D interpolation matrix of the given <a href="Basis.html#LibCEED.Basis"><code>Basis</code></a>. <code>b</code> must be a tensor-product basis, otherwise this function will fail. Returns a matrix of size <code>(getnumqpts1d(b), getnumnodes1d(b))</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/CEED/libCEED/blob/4018a20a98d451fac24765d3ddb936861647ce8d/julia/LibCEED.jl/src/Basis.jl#L482-L488">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="LibCEED.getgrad" href="#LibCEED.getgrad"><code>LibCEED.getgrad</code></a> — <span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">getgrad(b::Basis)</code></pre><p>Get the gradient matrix of the given <a href="Basis.html#LibCEED.Basis"><code>Basis</code></a>. Returns a tensor of size <code>(getdimension(b), getnumqpts(b), getnumnodes(b))</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/CEED/libCEED/blob/4018a20a98d451fac24765d3ddb936861647ce8d/julia/LibCEED.jl/src/Basis.jl#L497-L502">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="LibCEED.getgrad1d" href="#LibCEED.getgrad1d"><code>LibCEED.getgrad1d</code></a> — <span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">getgrad1d(b::Basis)</code></pre><p>Get the 1D derivative matrix of the given <a href="Basis.html#LibCEED.Basis"><code>Basis</code></a>. Returns a matrix of size <code>(getnumqpts(b), getnumnodes(b))</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/CEED/libCEED/blob/4018a20a98d451fac24765d3ddb936861647ce8d/julia/LibCEED.jl/src/Basis.jl#L512-L517">source</a></section></article><div class="admonition is-warning"><header class="admonition-header">Missing docstring.</header><div class="admonition-body"><p>Missing docstring for <code>getdiv</code>. Check Documenter's build log for details.</p></div></div><div class="admonition is-warning"><header class="admonition-header">Missing docstring.</header><div class="admonition-body"><p>Missing docstring for <code>getcurl</code>. Check Documenter's build log for details.</p></div></div></article><nav class="docs-footer"><a class="docs-footer-prevpage" href="ElemRestriction.html">« ElemRestriction</a><a class="docs-footer-nextpage" href="QFunction.html">QFunction »</a><div class="flexbox-break"></div><p class="footer-message">Powered by <a href="https://github.com/JuliaDocs/Documenter.jl">Documenter.jl</a> and the <a href="https://julialang.org/">Julia Programming Language</a>.</p></nav></div><div class="modal" id="documenter-settings"><div class="modal-background"></div><div class="modal-card"><header class="modal-card-head"><p class="modal-card-title">Settings</p><button class="delete"></button></header><section class="modal-card-body"><p><label class="label">Theme</label><div class="select"><select id="documenter-themepicker"><option value="documenter-light">documenter-light</option><option value="documenter-dark">documenter-dark</option></select></div></p><hr/><p>This document was generated with <a href="https://github.com/JuliaDocs/Documenter.jl">Documenter.jl</a> on <span class="colophon-date" title="Wednesday 1 November 2023 06:18">Wednesday 1 November 2023</span>. Using Julia version 1.9.3.</p></section><footer class="modal-card-foot"></footer></div></div></div></body></html>