1*edab6123Sjeremylt{ 2*edab6123Sjeremylt "cells": [ 3*edab6123Sjeremylt { 4*edab6123Sjeremylt "cell_type": "markdown", 5*edab6123Sjeremylt "metadata": {}, 6*edab6123Sjeremylt "source": [ 7*edab6123Sjeremylt "# libCEED for Python examples\n", 8*edab6123Sjeremylt "\n", 9*edab6123Sjeremylt "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", 10*edab6123Sjeremylt "\n", 11*edab6123Sjeremylt "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/)." 12*edab6123Sjeremylt ] 13*edab6123Sjeremylt }, 14*edab6123Sjeremylt { 15*edab6123Sjeremylt "cell_type": "markdown", 16*edab6123Sjeremylt "metadata": {}, 17*edab6123Sjeremylt "source": [ 18*edab6123Sjeremylt "## Setting up libCEED for Python\n", 19*edab6123Sjeremylt "\n", 20*edab6123Sjeremylt "Install libCEED for Python by running" 21*edab6123Sjeremylt ] 22*edab6123Sjeremylt }, 23*edab6123Sjeremylt { 24*edab6123Sjeremylt "cell_type": "code", 25*edab6123Sjeremylt "execution_count": null, 26*edab6123Sjeremylt "metadata": {}, 27*edab6123Sjeremylt "outputs": [], 28*edab6123Sjeremylt "source": [ 29*edab6123Sjeremylt "! python -m pip install libceed" 30*edab6123Sjeremylt ] 31*edab6123Sjeremylt }, 32*edab6123Sjeremylt { 33*edab6123Sjeremylt "cell_type": "markdown", 34*edab6123Sjeremylt "metadata": {}, 35*edab6123Sjeremylt "source": [ 36*edab6123Sjeremylt "## CeedOperator\n", 37*edab6123Sjeremylt "\n", 38*edab6123Sjeremylt "Here we show some basic examples to illustrate the `libceed.Operator` class. In libCEED, a `libceed.Operator` defines the finite/spectral element operator associated to a `libceed.QFunction` (see [the API documentation](https://libceed.readthedocs.io/en/latest/libCEEDapi.html#finite-element-operator-decomposition))." 39*edab6123Sjeremylt ] 40*edab6123Sjeremylt }, 41*edab6123Sjeremylt { 42*edab6123Sjeremylt "cell_type": "markdown", 43*edab6123Sjeremylt "metadata": {}, 44*edab6123Sjeremylt "source": [ 45*edab6123Sjeremylt "* In the following example, we create and apply a CeedOperator for the mass matrix in 1D. By applying this operator to a vector of 1's, we compute the length of this 1D domain, similar to Ex1-Volume in the [tutorial-6-shell tutorial](./tutorial-6-shell.ipynb)" 46*edab6123Sjeremylt ] 47*edab6123Sjeremylt }, 48*edab6123Sjeremylt { 49*edab6123Sjeremylt "cell_type": "code", 50*edab6123Sjeremylt "execution_count": null, 51*edab6123Sjeremylt "metadata": {}, 52*edab6123Sjeremylt "outputs": [], 53*edab6123Sjeremylt "source": [ 54*edab6123Sjeremylt "import libceed\n", 55*edab6123Sjeremylt "import numpy as np\n", 56*edab6123Sjeremylt "\n", 57*edab6123Sjeremylt "ceed = libceed.Ceed()\n", 58*edab6123Sjeremylt "\n", 59*edab6123Sjeremylt "nelem = 15\n", 60*edab6123Sjeremylt "p = 5\n", 61*edab6123Sjeremylt "q = 8\n", 62*edab6123Sjeremylt "nx = nelem + 1\n", 63*edab6123Sjeremylt "nu = nelem*(p-1) + 1\n", 64*edab6123Sjeremylt "\n", 65*edab6123Sjeremylt "# Vectors\n", 66*edab6123Sjeremylt "x = ceed.Vector(nx)\n", 67*edab6123Sjeremylt "x_array = np.zeros(nx)\n", 68*edab6123Sjeremylt "for i in range(nx):\n", 69*edab6123Sjeremylt " x_array[i] = i / (nx - 1.0)\n", 70*edab6123Sjeremylt "x.set_array(x_array, cmode=libceed.USE_POINTER)\n", 71*edab6123Sjeremylt "\n", 72*edab6123Sjeremylt "qdata = ceed.Vector(nelem*q)\n", 73*edab6123Sjeremylt "u = ceed.Vector(nu)\n", 74*edab6123Sjeremylt "v = ceed.Vector(nu)\n", 75*edab6123Sjeremylt "\n", 76*edab6123Sjeremylt "# Restrictions\n", 77*edab6123Sjeremylt "indx = np.zeros(nx*2, dtype=\"int32\")\n", 78*edab6123Sjeremylt "for i in range(nx):\n", 79*edab6123Sjeremylt " indx[2*i+0] = i\n", 80*edab6123Sjeremylt " indx[2*i+1] = i+1\n", 81*edab6123Sjeremylt "rx = ceed.ElemRestriction(nelem, 2, 1, 1, nx, indx, cmode=libceed.USE_POINTER)\n", 82*edab6123Sjeremylt "\n", 83*edab6123Sjeremylt "indu = np.zeros(nelem*p, dtype=\"int32\")\n", 84*edab6123Sjeremylt "for i in range(nelem):\n", 85*edab6123Sjeremylt " for j in range(p):\n", 86*edab6123Sjeremylt " indu[p*i+j] = i*(p-1) + j\n", 87*edab6123Sjeremylt "ru = ceed.ElemRestriction(nelem, p, 1, 1, nu, indu, cmode=libceed.USE_POINTER)\n", 88*edab6123Sjeremylt "strides = np.array([1, q, q], dtype=\"int32\")\n", 89*edab6123Sjeremylt "rui = ceed.StridedElemRestriction(nelem, q, 1, q*nelem, strides)\n", 90*edab6123Sjeremylt "\n", 91*edab6123Sjeremylt "# Bases\n", 92*edab6123Sjeremylt "bx = ceed.BasisTensorH1Lagrange(1, 1, 2, q, libceed.GAUSS)\n", 93*edab6123Sjeremylt "bu = ceed.BasisTensorH1Lagrange(1, 1, p, q, libceed.GAUSS)\n", 94*edab6123Sjeremylt "\n", 95*edab6123Sjeremylt "# QFunctions\n", 96*edab6123Sjeremylt "qf_setup = ceed.QFunctionByName(\"Mass1DBuild\")\n", 97*edab6123Sjeremylt "qf_mass = ceed.QFunctionByName(\"MassApply\")\n", 98*edab6123Sjeremylt "\n", 99*edab6123Sjeremylt "# Setup operator\n", 100*edab6123Sjeremylt "op_setup = ceed.Operator(qf_setup)\n", 101*edab6123Sjeremylt "op_setup.set_field(\"dx\", rx, bx, libceed.VECTOR_ACTIVE)\n", 102*edab6123Sjeremylt "op_setup.set_field(\"weights\", libceed.ELEMRESTRICTION_NONE, bx,\n", 103*edab6123Sjeremylt " libceed.VECTOR_NONE)\n", 104*edab6123Sjeremylt "op_setup.set_field(\"qdata\", rui, libceed.BASIS_COLLOCATED,\n", 105*edab6123Sjeremylt " libceed.VECTOR_ACTIVE)\n", 106*edab6123Sjeremylt "print('Setup operator: ', op_setup)\n", 107*edab6123Sjeremylt "\n", 108*edab6123Sjeremylt "# Mass operator\n", 109*edab6123Sjeremylt "op_mass = ceed.Operator(qf_mass)\n", 110*edab6123Sjeremylt "op_mass.set_field(\"u\", ru, bu, libceed.VECTOR_ACTIVE)\n", 111*edab6123Sjeremylt "op_mass.set_field(\"qdata\", rui, libceed.BASIS_COLLOCATED, qdata)\n", 112*edab6123Sjeremylt "op_mass.set_field(\"v\", ru, bu, libceed.VECTOR_ACTIVE)\n", 113*edab6123Sjeremylt "print('Mass operator: ', op_mass)\n", 114*edab6123Sjeremylt "\n", 115*edab6123Sjeremylt "# Setup\n", 116*edab6123Sjeremylt "op_setup.apply(x, qdata)\n", 117*edab6123Sjeremylt "\n", 118*edab6123Sjeremylt "# Apply mass matrix\n", 119*edab6123Sjeremylt "u.set_value(1)\n", 120*edab6123Sjeremylt "op_mass.apply(u, v)\n", 121*edab6123Sjeremylt "\n", 122*edab6123Sjeremylt "# Check\n", 123*edab6123Sjeremylt "with v.array_read() as v_array:\n", 124*edab6123Sjeremylt " print('The length of the domain is l = %4.2f'%np.sum(v_array))" 125*edab6123Sjeremylt ] 126*edab6123Sjeremylt } 127*edab6123Sjeremylt ], 128*edab6123Sjeremylt "metadata": { 129*edab6123Sjeremylt "kernelspec": { 130*edab6123Sjeremylt "display_name": "Python 3", 131*edab6123Sjeremylt "language": "python", 132*edab6123Sjeremylt "name": "python3" 133*edab6123Sjeremylt }, 134*edab6123Sjeremylt "language_info": { 135*edab6123Sjeremylt "codemirror_mode": { 136*edab6123Sjeremylt "name": "ipython", 137*edab6123Sjeremylt "version": 3 138*edab6123Sjeremylt }, 139*edab6123Sjeremylt "file_extension": ".py", 140*edab6123Sjeremylt "mimetype": "text/x-python", 141*edab6123Sjeremylt "name": "python", 142*edab6123Sjeremylt "nbconvert_exporter": "python", 143*edab6123Sjeremylt "pygments_lexer": "ipython3", 144*edab6123Sjeremylt "version": "3.8.5" 145*edab6123Sjeremylt } 146*edab6123Sjeremylt }, 147*edab6123Sjeremylt "nbformat": 4, 148*edab6123Sjeremylt "nbformat_minor": 4 149*edab6123Sjeremylt} 150