{ "cells": [ { "cell_type": "markdown", "metadata": {}, "source": [ "# Standalone libCEED examples\n", "\n", "This is a tutorial for [libCEED](https://github.com/CEED/libCEED/), the low-level API library for efficient high-order discretization methods developed by the co-design [Center for Efficient Exascale Discretizations](https://ceed.exascaleproject.org/) (CEED) of the [Exascale Computing Project](https://www.exascaleproject.org/) (ECP).\n", "\n", "While libCEED's focus is on high-order finite/spectral element method implementations, the approach is mostly algebraic and thus applicable to other discretizations in factored form, as explained in the [user manual](https://libceed.readthedocs.io/)." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Common notation\n", "\n", "For most of our examples, the spatial discretization\n", "uses high-order finite elements/spectral elements, namely, the high-order Lagrange\n", "polynomials defined over $P$ non-uniformly spaced nodes, the\n", "Gauss-Legendre-Lobatto (GLL) points, and quadrature points $\\{q_i\\}_{i=1}^Q$, with\n", "corresponding weights $\\{w_i\\}_{i=1}^Q$ (typically the ones given by Gauss\n", "or Gauss-Lobatto quadratures, that are built in the library).\n", "\n", "We discretize the domain, $\\Omega \\subset \\mathbb{R}^d$ (with $d=1,2,3$,\n", "typically) by letting $\\Omega = \\bigcup_{e=1}^{N_e}\\Omega_e$, with $N_e$\n", "disjoint elements. For most examples we use unstructured meshes for which the elements\n", "are hexahedra (although this is not a requirement in libCEED).\n", "\n", "The physical coordinates are denoted by $\\mathbf{x}=(x,y,z)\\in\\Omega_e$,\n", "while the reference coordinates are represented as\n", "$\\boldsymbol{X}=(X,Y,Z) \\equiv (X_1,X_2,X_3) \\in I=[-1,1]^3$\n", "(for $d=3$).\n", "\n", "### Ex1-Volume\n", "\n", "#### Mathematical formulation\n", "\n", "This example is located in the subdirectory `examples/ceed`. It illustrates a\n", "simple usage of libCEED to compute the volume of a given body using a matrix-free\n", "application of the mass operator. Arbitrary mesh and solution orders in 1D, 2D and 3D\n", "are supported from the same code.\n", "\n", "This example shows how to compute line/surface/volume integrals of a 1D, 2D, or 3D\n", "domain $\\Omega$ respectively, by applying the mass operator to a vector of\n", "$1$s. It computes:\n", "\n", "$$\n", " I = \\int_{\\Omega} 1 \\, dV .\n", "$$\n", "\n", "We write here the vector $u(x)\\equiv 1$ in the Galerkin approximation,\n", "and find the volume of $\\Omega$ as\n", "\n", "$$\n", " \\sum_e \\int_{\\Omega_e} v(x) 1 \\, dV\n", "$$\n", "\n", "with $v(x) \\in \\mathcal{V}_p = \\{ v \\in H^{1}(\\Omega_e) \\,|\\, v \\in P_p(\\boldsymbol{I}), e=1,\\ldots,N_e \\}$,\n", "the test functions.\n" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "#### Executing the example\n", "\n", "Clone or download libCEED by running" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "! git clone https://github.com/CEED/libCEED.git" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Then move to the libCEED folder" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "cd libCEED" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "And compile the library by running" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "! make" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Move to the examples folder " ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "cd examples/" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Then move to the standalone libCEED's examples folder" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "cd ceed/" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "And compile the examples by running" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "! make" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Now run `ex1-volume` by running" ] }, { "cell_type": "code", "execution_count": 1, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Selected options: [command line option] : \n", " Ceed specification [-c] : /cpu/self\n", " Mesh dimension [-d] : 3\n", " Mesh order [-m] : 4\n", " Solution order [-o] : 4\n", " Num. 1D quadr. pts [-q] : 6\n", " Approx. # unknowns [-s] : 262144\n", " QFunction source [-g] : gallery\n", "\n", "Mesh size: nx = 16, ny = 16, nz = 16\n", "Number of mesh nodes : 274625\n", "Number of solution nodes : 274625\n", "Computing the quadrature data for the mass operator ... done.\n", "Computing the mesh volume using the formula: vol = 1^T.M.1 ... done.\n", "Exact mesh volume : 2.3561944901923\n", "Computed mesh volume : 2.3561944901921\n", "Volume error : -2.7444713168734e-13\n" ] } ], "source": [ "! ./ex1-volume -d 3 -g" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "This example shows how to compute line/surface/volume integrals of a 1D, 2D, or 3D domain Ω respectively, by applying the mass operator to a vector of 1s. The command line option `-d` specifies the dimensionality of the domain Ω. The option `-g` specifies that the mass operator is, in this case, selected from a gallery of available built-in operators in the library." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Ex2-Surface\n", "\n", "#### Mathematical formulation\n", "\n", "This example is located in the subdirectory `examples/ceed`. It computes the\n", "surface area of a given body using matrix-free application of a diffusion operator.\n", "Arbitrary mesh and solution orders in 1D, 2D and 3D are supported from the same code.\n", "\n", "Similarly to the Ex1-Volume example, it computes:\n", "\n", "\\begin{equation}\\label{eq-ex2-surface}\\tag{eq. 1}\n", " I = \\int_{\\partial \\Omega} 1 \\, dS ,\n", "\\end{equation}\n", "\n", "by applying the divergence theorem.\n", "In particular, we select $u(\\mathbf x) = x_0 + x_1 + x_2$, for which $\\nabla u = [1, 1, 1]^T$, and thus $\\nabla u \\cdot \\hat{\\mathbf n} = 1$.\n", "\n", "Given Laplace's equation,\n", "\n", "$$\n", " \\nabla \\cdot \\nabla u = 0, \\textrm{ for } \\mathbf{x} \\in \\Omega ,\n", "$$\n", "\n", "let us multiply by a test function $v$ and integrate by parts to obtain\n", "\n", "$$\n", " \\int_\\Omega \\nabla v \\cdot \\nabla u \\, dV - \\int_{\\partial \\Omega} v \\nabla u \\cdot \\hat{\\mathbf n}\\, dS = 0 .\n", "$$\n", "\n", "Since we have chosen $u$ such that $\\nabla u \\cdot \\hat{\\mathbf n} = 1$, the boundary integrand is $v 1 \\equiv v$. Hence, similar to the previous example, we can evaluate the surface integral by applying the volumetric Laplacian as follows\n", "\n", "$$\n", " \\int_\\Omega \\nabla v \\cdot \\nabla u \\, dV \\approx \\sum_e \\int_{\\partial \\Omega_e} v(x) 1 \\, dS .\n", "$$" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Executing the example\n", "\n", "Assuming the steps above, you should be in the `examples/ceed/` subdirectory and have already compiled the example.\n", "\n", "Now run `ex2-surface` by running " ] }, { "cell_type": "code", "execution_count": 2, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Selected options: [command line option] : \n", " Ceed specification [-c] : /cpu/self\n", " Mesh dimension [-d] : 3\n", " Mesh order [-m] : 4\n", " Solution order [-o] : 4\n", " Num. 1D quadr. pts [-q] : 6\n", " Approx. # unknowns [-s] : 262144\n", " QFunction source [-g] : gallery\n", "\n", "Mesh size: nx = 16, ny = 16, nz = 16\n", "Number of mesh nodes : 274625\n", "Number of solution nodes : 274625\n", "Computing the quadrature data for the diffusion operator ... done.\n", "Computing the mesh surface area using the formula: sa = 1^T.|K.x| ... done.\n", "Exact mesh surface area : 6\n", "Computed mesh surface area : 5.9773703490853\n", "Surface area error : -0.022629650914673\n" ] } ], "source": [ "! ./ex2-surface -d 3 -g" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "This example computes the surface area of a given body using matrix-free application of a Laplace's (diffusion) operator. The command line option `-d` specifies the dimensionality of the domain Ω. The option `-g` specifies that the Laplace's operator is, in this case, selected from a gallery of available built-in operators in the library." ] } ], "metadata": { "kernelspec": { "display_name": "Python 3", "language": "python", "name": "python3" }, "language_info": { "codemirror_mode": { "name": "ipython", "version": 3 }, "file_extension": ".py", "mimetype": "text/x-python", "name": "python", "nbconvert_exporter": "python", "pygments_lexer": "ipython3", "version": "3.8.5" } }, "nbformat": 4, "nbformat_minor": 4 }