1edab6123Sjeremylt{ 2edab6123Sjeremylt "cells": [ 3edab6123Sjeremylt { 4edab6123Sjeremylt "cell_type": "markdown", 5edab6123Sjeremylt "metadata": {}, 6edab6123Sjeremylt "source": [ 7edab6123Sjeremylt "# libCEED for Python examples\n", 8edab6123Sjeremylt "\n", 9edab6123Sjeremylt "This is a tutorial to illustrate the main feautures of the Python interface 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", 10edab6123Sjeremylt "\n", 1113964f07SJed Brown "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.org/)." 12edab6123Sjeremylt ] 13edab6123Sjeremylt }, 14edab6123Sjeremylt { 15edab6123Sjeremylt "cell_type": "markdown", 16edab6123Sjeremylt "metadata": {}, 17edab6123Sjeremylt "source": [ 18edab6123Sjeremylt "## Setting up libCEED for Python\n", 19edab6123Sjeremylt "\n", 20edab6123Sjeremylt "Install libCEED for Python by running" 21edab6123Sjeremylt ] 22edab6123Sjeremylt }, 23edab6123Sjeremylt { 24edab6123Sjeremylt "cell_type": "code", 25edab6123Sjeremylt "execution_count": null, 26edab6123Sjeremylt "metadata": {}, 27edab6123Sjeremylt "outputs": [], 28edab6123Sjeremylt "source": [ 29edab6123Sjeremylt "! python -m pip install libceed" 30edab6123Sjeremylt ] 31edab6123Sjeremylt }, 32edab6123Sjeremylt { 33edab6123Sjeremylt "cell_type": "markdown", 34edab6123Sjeremylt "metadata": {}, 35edab6123Sjeremylt "source": [ 36edab6123Sjeremylt "## CeedBasis\n", 37edab6123Sjeremylt "\n", 3813964f07SJed Brown "Here we show some basic examples to illustrate the `libceed.Basis` class. In libCEED, a `libceed.Basis` defines the finite element basis and associated quadrature rule (see [the API documentation](https://libceed.org/en/latest/libCEEDapi.html#finite-element-operator-decomposition))." 39edab6123Sjeremylt ] 40edab6123Sjeremylt }, 41edab6123Sjeremylt { 42edab6123Sjeremylt "cell_type": "markdown", 43edab6123Sjeremylt "metadata": {}, 44edab6123Sjeremylt "source": [ 45edab6123Sjeremylt "First we declare some auxiliary functions needed in the following examples" 46edab6123Sjeremylt ] 47edab6123Sjeremylt }, 48edab6123Sjeremylt { 49edab6123Sjeremylt "cell_type": "code", 50edab6123Sjeremylt "execution_count": null, 51edab6123Sjeremylt "metadata": {}, 52edab6123Sjeremylt "outputs": [], 53edab6123Sjeremylt "source": [ 54edab6123Sjeremylt "%matplotlib inline\n", 55edab6123Sjeremylt "import numpy as np\n", 56edab6123Sjeremylt "import matplotlib.pyplot as plt\n", 57edab6123Sjeremylt "plt.style.use('ggplot')\n", 58edab6123Sjeremylt "\n", 59edab6123Sjeremylt "def eval(dim, x):\n", 60edab6123Sjeremylt " result, center = 1, 0.1\n", 61edab6123Sjeremylt " for d in range(dim):\n", 62edab6123Sjeremylt " result *= np.tanh(x[d] - center)\n", 63edab6123Sjeremylt " center += 0.1\n", 64edab6123Sjeremylt " return result\n", 65edab6123Sjeremylt "\n", 66*235e4399SJeremy L Thompson "def feval(x_1, x_2):\n", 67*235e4399SJeremy L Thompson " return x_1*x_1 + x_2*x_2 + x_1*x_2 + 1\n", 68edab6123Sjeremylt "\n", 69*235e4399SJeremy L Thompson "def dfeval(x_1, x_2):\n", 70*235e4399SJeremy L Thompson " return 2*x_1 + x_2" 71edab6123Sjeremylt ] 72edab6123Sjeremylt }, 73edab6123Sjeremylt { 74edab6123Sjeremylt "cell_type": "markdown", 75edab6123Sjeremylt "metadata": {}, 76edab6123Sjeremylt "source": [ 77edab6123Sjeremylt "## $H^1$ Lagrange bases in 1D\n", 78edab6123Sjeremylt "\n", 79edab6123Sjeremylt "The Lagrange interpolation nodes are at the Gauss-Lobatto points, so interpolation to Gauss-Lobatto quadrature points is the identity." 80edab6123Sjeremylt ] 81edab6123Sjeremylt }, 82edab6123Sjeremylt { 83edab6123Sjeremylt "cell_type": "code", 84edab6123Sjeremylt "execution_count": null, 85edab6123Sjeremylt "metadata": {}, 86edab6123Sjeremylt "outputs": [], 87edab6123Sjeremylt "source": [ 88edab6123Sjeremylt "import libceed\n", 89edab6123Sjeremylt "\n", 90edab6123Sjeremylt "ceed = libceed.Ceed()\n", 91edab6123Sjeremylt "\n", 92edab6123Sjeremylt "b = ceed.BasisTensorH1Lagrange(\n", 93edab6123Sjeremylt " dim=1, # topological dimension\n", 94edab6123Sjeremylt " ncomp=1, # number of components\n", 95edab6123Sjeremylt " P=4, # number of basis functions (nodes) per dimension\n", 96edab6123Sjeremylt " Q=4, # number of quadrature points per dimension\n", 97edab6123Sjeremylt " qmode=libceed.GAUSS_LOBATTO)\n", 98edab6123Sjeremylt "print(b)" 99edab6123Sjeremylt ] 100edab6123Sjeremylt }, 101edab6123Sjeremylt { 102edab6123Sjeremylt "cell_type": "markdown", 103edab6123Sjeremylt "metadata": {}, 104edab6123Sjeremylt "source": [ 105edab6123Sjeremylt "Although a `libceed.Basis` is fully discrete, we can use the Lagrange construction to extend the basis to continuous functions by applying `EVAL_INTERP` to the identity. This is the Vandermonde matrix of the continuous basis." 106edab6123Sjeremylt ] 107edab6123Sjeremylt }, 108edab6123Sjeremylt { 109edab6123Sjeremylt "cell_type": "code", 110edab6123Sjeremylt "execution_count": null, 111edab6123Sjeremylt "metadata": {}, 112edab6123Sjeremylt "outputs": [], 113edab6123Sjeremylt "source": [ 114edab6123Sjeremylt "P = b.get_num_nodes()\n", 115*235e4399SJeremy L Thompson "Q_viz = 50\n", 116*235e4399SJeremy L Thompson "basis_viz = ceed.BasisTensorH1Lagrange(1, 1, P, Q_viz, libceed.GAUSS_LOBATTO)\n", 117edab6123Sjeremylt "\n", 118edab6123Sjeremylt "# Construct P \"elements\" with one node activated\n", 119edab6123Sjeremylt "I = ceed.Vector(P * P)\n", 120*235e4399SJeremy L Thompson "with I.array_write(P, P) as x:\n", 121edab6123Sjeremylt " x[...] = np.eye(P)\n", 122edab6123Sjeremylt "\n", 123*235e4399SJeremy L Thompson "basis_fns = ceed.Vector(P * Q_viz)\n", 124*235e4399SJeremy L Thompson "basis_viz.apply(4, libceed.EVAL_INTERP, I, basis_fns)\n", 125edab6123Sjeremylt "\n", 126*235e4399SJeremy L Thompson "qpts_viz, _ = ceed.lobatto_quadrature(Q_viz)\n", 127*235e4399SJeremy L Thompson "with basis_fns.array_read(Q_viz, P) as B_array:\n", 128*235e4399SJeremy L Thompson " plt.plot(qpts_viz, B_array)\n", 129edab6123Sjeremylt "\n", 130edab6123Sjeremylt "# Mark tho Lobatto nodes\n", 131*235e4399SJeremy L Thompson "nodes, _ = ceed.lobatto_quadrature(P)\n", 132*235e4399SJeremy L Thompson "plt.plot(nodes, 0*nodes, 'ok');" 133edab6123Sjeremylt ] 134edab6123Sjeremylt }, 135edab6123Sjeremylt { 136edab6123Sjeremylt "cell_type": "markdown", 137edab6123Sjeremylt "metadata": {}, 138edab6123Sjeremylt "source": [ 139edab6123Sjeremylt "In contrast, the Gauss quadrature points are not collocated, and thus all basis functions are generally nonzero at every quadrature point." 140edab6123Sjeremylt ] 141edab6123Sjeremylt }, 142edab6123Sjeremylt { 143edab6123Sjeremylt "cell_type": "code", 144edab6123Sjeremylt "execution_count": null, 145edab6123Sjeremylt "metadata": {}, 146edab6123Sjeremylt "outputs": [], 147edab6123Sjeremylt "source": [ 148edab6123Sjeremylt "b = ceed.BasisTensorH1Lagrange(1, 1, 4, 4, libceed.GAUSS)\n", 149edab6123Sjeremylt "print(b)\n", 150edab6123Sjeremylt "\n", 151*235e4399SJeremy L Thompson "with basis_fns.array_read(Q_viz, P) as B_array:\n", 152*235e4399SJeremy L Thompson " plt.plot(qpts_viz, B_array)\n", 153edab6123Sjeremylt "# Mark tho Gauss quadrature points\n", 154*235e4399SJeremy L Thompson "qpts, _ = ceed.gauss_quadrature(P)\n", 155*235e4399SJeremy L Thompson "plt.plot(qpts, 0*qpts, 'ok');" 156edab6123Sjeremylt ] 157edab6123Sjeremylt }, 158edab6123Sjeremylt { 159edab6123Sjeremylt "cell_type": "markdown", 160edab6123Sjeremylt "metadata": {}, 161edab6123Sjeremylt "source": [ 162edab6123Sjeremylt "Although the underlying functions are not an intrinsic property of a `libceed.Basis` in libCEED, the sizes are.\n", 163edab6123Sjeremylt "Here, we create a 3D tensor product element with more quadrature points than Lagrange interpolation nodes." 164edab6123Sjeremylt ] 165edab6123Sjeremylt }, 166edab6123Sjeremylt { 167edab6123Sjeremylt "cell_type": "code", 168edab6123Sjeremylt "execution_count": null, 169edab6123Sjeremylt "metadata": {}, 170edab6123Sjeremylt "outputs": [], 171edab6123Sjeremylt "source": [ 172edab6123Sjeremylt "b = ceed.BasisTensorH1Lagrange(3, 1, 4, 5, libceed.GAUSS_LOBATTO)\n", 173edab6123Sjeremylt "\n", 174edab6123Sjeremylt "p = b.get_num_nodes()\n", 175edab6123Sjeremylt "print('p =', p)\n", 176edab6123Sjeremylt "\n", 177edab6123Sjeremylt "q = b.get_num_quadrature_points()\n", 178edab6123Sjeremylt "print('q =', q)" 179edab6123Sjeremylt ] 180edab6123Sjeremylt }, 181edab6123Sjeremylt { 182edab6123Sjeremylt "cell_type": "markdown", 183edab6123Sjeremylt "metadata": {}, 184edab6123Sjeremylt "source": [ 185edab6123Sjeremylt "* In the following example, we demonstrate the application of an interpolatory basis in multiple dimensions" 186edab6123Sjeremylt ] 187edab6123Sjeremylt }, 188edab6123Sjeremylt { 189edab6123Sjeremylt "cell_type": "code", 190edab6123Sjeremylt "execution_count": null, 191edab6123Sjeremylt "metadata": {}, 192edab6123Sjeremylt "outputs": [], 193edab6123Sjeremylt "source": [ 194edab6123Sjeremylt "for dim in range(1, 4):\n", 195edab6123Sjeremylt " Q = 4\n", 196*235e4399SJeremy L Thompson " Q_dim = Q**dim\n", 197*235e4399SJeremy L Thompson " X_dim = 2**dim\n", 198*235e4399SJeremy L Thompson " x = np.empty(X_dim*dim, dtype=\"float64\")\n", 199*235e4399SJeremy L Thompson " u_array = np.empty(Q_dim, dtype=\"float64\")\n", 200edab6123Sjeremylt "\n", 201edab6123Sjeremylt " for d in range(dim):\n", 202*235e4399SJeremy L Thompson " for i in range(X_dim):\n", 203*235e4399SJeremy L Thompson " x[d*X_dim + i] = 1 if (i % (2**(dim-d))) // (2**(dim-d-1)) else -1\n", 204edab6123Sjeremylt "\n", 205*235e4399SJeremy L Thompson " X = ceed.Vector(X_dim*dim)\n", 206edab6123Sjeremylt " X.set_array(x, cmode=libceed.USE_POINTER)\n", 207*235e4399SJeremy L Thompson " X_q = ceed.Vector(Q_dim*dim)\n", 208*235e4399SJeremy L Thompson " X_q.set_value(0)\n", 209*235e4399SJeremy L Thompson " U = ceed.Vector(Q_dim)\n", 210edab6123Sjeremylt " U.set_value(0)\n", 211*235e4399SJeremy L Thompson " U_q = ceed.Vector(Q_dim)\n", 212edab6123Sjeremylt "\n", 213*235e4399SJeremy L Thompson " basis_x_lobatto = ceed.BasisTensorH1Lagrange(dim, dim, 2, Q, libceed.GAUSS_LOBATTO)\n", 214*235e4399SJeremy L Thompson " basis_u_lobatto = ceed.BasisTensorH1Lagrange(dim, 1, Q, Q, libceed.GAUSS_LOBATTO)\n", 215edab6123Sjeremylt "\n", 216*235e4399SJeremy L Thompson " basis_x_lobatto.apply(1, libceed.EVAL_INTERP, X, X_q)\n", 217edab6123Sjeremylt "\n", 218*235e4399SJeremy L Thompson " with X_q.array_read() as x_array:\n", 219*235e4399SJeremy L Thompson " for i in range(Q_dim):\n", 220*235e4399SJeremy L Thompson " x = np.empty(dim, dtype=\"float64\")\n", 221edab6123Sjeremylt " for d in range(dim):\n", 222*235e4399SJeremy L Thompson " x[d] = x_array[d*Q_dim + i]\n", 223*235e4399SJeremy L Thompson " u_array[i] = eval(dim, x)\n", 224edab6123Sjeremylt "\n", 225*235e4399SJeremy L Thompson " U_q.set_array(u_array, cmode=libceed.USE_POINTER)\n", 226edab6123Sjeremylt "\n", 227edab6123Sjeremylt " # This operation is the identity because the quadrature is collocated\n", 228*235e4399SJeremy L Thompson " basis_u_lobatto.T.apply(1, libceed.EVAL_INTERP, U_q, U)\n", 229edab6123Sjeremylt "\n", 230*235e4399SJeremy L Thompson " basis_x_gauss = ceed.BasisTensorH1Lagrange(dim, dim, 2, Q, libceed.GAUSS)\n", 231*235e4399SJeremy L Thompson " basis_u_gauss = ceed.BasisTensorH1Lagrange(dim, 1, Q, Q, libceed.GAUSS)\n", 232edab6123Sjeremylt "\n", 233*235e4399SJeremy L Thompson " basis_x_gauss.apply(1, libceed.EVAL_INTERP, X, X_q)\n", 234*235e4399SJeremy L Thompson " basis_u_gauss.apply(1, libceed.EVAL_INTERP, U, U_q)\n", 235edab6123Sjeremylt "\n", 236*235e4399SJeremy L Thompson " with X_q.array_read() as x_array, U_q.array_read() as u_array:\n", 237edab6123Sjeremylt " if dim == 2:\n", 238edab6123Sjeremylt " # Default ordering is contiguous in x direction, but\n", 239edab6123Sjeremylt " # pyplot expects meshgrid convention, which is transposed.\n", 240*235e4399SJeremy L Thompson " x, y = x_array.reshape(2, Q, Q).transpose(0, 2, 1)\n", 241*235e4399SJeremy L Thompson " plt.scatter(x, y, c=np.array(u_array).reshape(Q, Q))\n", 242edab6123Sjeremylt " plt.xlim(-1, 1)\n", 243edab6123Sjeremylt " plt.ylim(-1, 1)\n", 244edab6123Sjeremylt " plt.colorbar(label='u')" 245edab6123Sjeremylt ] 246edab6123Sjeremylt }, 247edab6123Sjeremylt { 248edab6123Sjeremylt "cell_type": "markdown", 249edab6123Sjeremylt "metadata": {}, 250edab6123Sjeremylt "source": [ 251edab6123Sjeremylt "* In the following example, we demonstrate the application of the gradient of the shape functions in multiple dimensions" 252edab6123Sjeremylt ] 253edab6123Sjeremylt }, 254edab6123Sjeremylt { 255edab6123Sjeremylt "cell_type": "code", 256edab6123Sjeremylt "execution_count": null, 257edab6123Sjeremylt "metadata": {}, 258edab6123Sjeremylt "outputs": [], 259edab6123Sjeremylt "source": [ 260edab6123Sjeremylt "for dim in range (1, 4):\n", 261edab6123Sjeremylt " P, Q = 8, 10\n", 262*235e4399SJeremy L Thompson " P_dim = P**dim\n", 263*235e4399SJeremy L Thompson " Q_dim = Q**dim\n", 264*235e4399SJeremy L Thompson " X_dim = 2**dim\n", 265*235e4399SJeremy L Thompson " sum_1 = sum_2 = 0\n", 266*235e4399SJeremy L Thompson " x_array = np.empty(X_dim*dim, dtype=\"float64\")\n", 267*235e4399SJeremy L Thompson " u_array = np.empty(P_dim, dtype=\"float64\")\n", 268edab6123Sjeremylt "\n", 269edab6123Sjeremylt " for d in range(dim):\n", 270*235e4399SJeremy L Thompson " for i in range(X_dim):\n", 271*235e4399SJeremy L Thompson " x_array[d*X_dim + i] = 1 if (i % (2**(dim-d))) // (2**(dim-d-1)) else -1\n", 272edab6123Sjeremylt "\n", 273*235e4399SJeremy L Thompson " X = ceed.Vector(X_dim*dim)\n", 274*235e4399SJeremy L Thompson " X.set_array(x_array, cmode=libceed.USE_POINTER)\n", 275*235e4399SJeremy L Thompson " X_q = ceed.Vector(P_dim*dim)\n", 276*235e4399SJeremy L Thompson " X_q.set_value(0)\n", 277*235e4399SJeremy L Thompson " U = ceed.Vector(P_dim)\n", 278*235e4399SJeremy L Thompson " U_q = ceed.Vector(Q_dim*dim)\n", 279*235e4399SJeremy L Thompson " U_q.set_value(0)\n", 280*235e4399SJeremy L Thompson " Ones = ceed.Vector(Q_dim*dim)\n", 281edab6123Sjeremylt " Ones.set_value(1)\n", 282*235e4399SJeremy L Thompson " G_transpose_ones = ceed.Vector(P_dim)\n", 283*235e4399SJeremy L Thompson " G_transpose_ones.set_value(0)\n", 284edab6123Sjeremylt "\n", 285edab6123Sjeremylt " # Get function values at quadrature points\n", 286*235e4399SJeremy L Thompson " basis_x_lobatto = ceed.BasisTensorH1Lagrange(dim, dim, 2, P, libceed.GAUSS_LOBATTO)\n", 287*235e4399SJeremy L Thompson " basis_x_lobatto.apply(1, libceed.EVAL_INTERP, X, X_q)\n", 288edab6123Sjeremylt "\n", 289*235e4399SJeremy L Thompson " with X_q.array_read() as x_array:\n", 290*235e4399SJeremy L Thompson " for i in range(P_dim):\n", 291*235e4399SJeremy L Thompson " x = np.empty(dim, dtype=\"float64\")\n", 292edab6123Sjeremylt " for d in range(dim):\n", 293*235e4399SJeremy L Thompson " x[d] = x_array[d*P_dim + i]\n", 294*235e4399SJeremy L Thompson " u_array[i] = eval(dim, x)\n", 295edab6123Sjeremylt "\n", 296*235e4399SJeremy L Thompson " U.set_array(u_array, cmode=libceed.USE_POINTER)\n", 297edab6123Sjeremylt "\n", 298edab6123Sjeremylt " # Calculate G u at quadrature points, G' * 1 at dofs\n", 299*235e4399SJeremy L Thompson " basis_u_gauss = ceed.BasisTensorH1Lagrange(dim, 1, P, Q, libceed.GAUSS)\n", 300*235e4399SJeremy L Thompson " basis_u_gauss.apply(1, libceed.EVAL_GRAD, U, U_q)\n", 301*235e4399SJeremy L Thompson " basis_u_gauss.T.apply(1, libceed.EVAL_GRAD, Ones, G_transpose_ones)\n", 302edab6123Sjeremylt "\n", 303edab6123Sjeremylt " # Check if 1' * G * u = u' * (G' * 1)\n", 304*235e4399SJeremy L Thompson " with G_transpose_ones.array_read() as g_array, U_q.array_read() as uq_array:\n", 305*235e4399SJeremy L Thompson " for i in range(P_dim):\n", 306*235e4399SJeremy L Thompson " sum_1 += g_array[i]*u_array[i]\n", 307*235e4399SJeremy L Thompson " for i in range(dim*Q_dim):\n", 308*235e4399SJeremy L Thompson " sum_2 += uq_array[i]\n", 309edab6123Sjeremylt "\n", 310edab6123Sjeremylt " # Check that (1' * G * u - u' * (G' * 1)) is numerically zero\n", 311*235e4399SJeremy L Thompson " print('1T * G * u - uT * (GT * 1) =', np.abs(sum_1 - sum_2))" 312edab6123Sjeremylt ] 313edab6123Sjeremylt } 314edab6123Sjeremylt ], 315edab6123Sjeremylt "metadata": { 316edab6123Sjeremylt "kernelspec": { 317*235e4399SJeremy L Thompson "display_name": "Python 3 (ipykernel)", 318edab6123Sjeremylt "language": "python", 319edab6123Sjeremylt "name": "python3" 320edab6123Sjeremylt }, 321edab6123Sjeremylt "language_info": { 322edab6123Sjeremylt "codemirror_mode": { 323edab6123Sjeremylt "name": "ipython", 324edab6123Sjeremylt "version": 3 325edab6123Sjeremylt }, 326edab6123Sjeremylt "file_extension": ".py", 327edab6123Sjeremylt "mimetype": "text/x-python", 328edab6123Sjeremylt "name": "python", 329edab6123Sjeremylt "nbconvert_exporter": "python", 330edab6123Sjeremylt "pygments_lexer": "ipython3", 331*235e4399SJeremy L Thompson "version": "3.13.2" 332edab6123Sjeremylt } 333edab6123Sjeremylt }, 334edab6123Sjeremylt "nbformat": 4, 335edab6123Sjeremylt "nbformat_minor": 4 336edab6123Sjeremylt} 337