xref: /libCEED/examples/python/tutorial-4-qfunction.ipynb (revision 7d38c034d40a24f4b961cafe9ae1ebfab6a4cf3f)
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    "## CeedQFunction\n",
37edab6123Sjeremylt    "\n",
3813964f07SJed Brown    "Here we show some basic examples to illustrate the `libceed.QFunction` class. In libCEED, QFunctions represent the spatial terms of the point-wise functions describing the physics at the quadrature points (see [the API documentation](https://libceed.org/en/latest/libCEEDapi.html#api-description)). As shown in the following sketch, QFunctions (such as the one depicted, which defines the Laplacian) are point-wise functions defined at quadrature points. Hence, QFunctions are independent from element shape, resolution and order.\n",
39edab6123Sjeremylt    "\n",
40edab6123Sjeremylt    "![alt text][QFunctionSchematic]\n",
41edab6123Sjeremylt    "\n",
42edab6123Sjeremylt    "[QFunctionSchematic]: ./img/QFunctionSketch.svg \"Schematic of point-wise QFunctions, defined at quadrature points, belonging to elements that can have different shape, resolution and order.\""
43edab6123Sjeremylt   ]
44edab6123Sjeremylt  },
45edab6123Sjeremylt  {
46edab6123Sjeremylt   "cell_type": "markdown",
47edab6123Sjeremylt   "metadata": {},
48edab6123Sjeremylt   "source": [
49edab6123Sjeremylt    "* In the following example, we create and view two QFunctions (for the setup and apply, respectively, of the mass operator in 1D) from the gallery of available built-in QFunctions in libCEED"
50edab6123Sjeremylt   ]
51edab6123Sjeremylt  },
52edab6123Sjeremylt  {
53edab6123Sjeremylt   "cell_type": "code",
54edab6123Sjeremylt   "execution_count": null,
55edab6123Sjeremylt   "metadata": {},
56edab6123Sjeremylt   "outputs": [],
57edab6123Sjeremylt   "source": [
58edab6123Sjeremylt    "import libceed\n",
59edab6123Sjeremylt    "import numpy as np\n",
60edab6123Sjeremylt    "\n",
61edab6123Sjeremylt    "ceed = libceed.Ceed()\n",
62edab6123Sjeremylt    "\n",
63edab6123Sjeremylt    "qf_setup = ceed.QFunctionByName(\"Mass1DBuild\")\n",
64edab6123Sjeremylt    "qf_mass = ceed.QFunctionByName(\"MassApply\")\n",
65edab6123Sjeremylt    "\n",
66edab6123Sjeremylt    "print(qf_setup)\n",
67edab6123Sjeremylt    "print(qf_mass)"
68edab6123Sjeremylt   ]
69edab6123Sjeremylt  },
70edab6123Sjeremylt  {
71edab6123Sjeremylt   "cell_type": "markdown",
72edab6123Sjeremylt   "metadata": {},
73edab6123Sjeremylt   "source": [
74edab6123Sjeremylt    "* In the following example, we create and evaluate a built-in identity QFunction."
75edab6123Sjeremylt   ]
76edab6123Sjeremylt  },
77edab6123Sjeremylt  {
78edab6123Sjeremylt   "cell_type": "code",
79edab6123Sjeremylt   "execution_count": null,
80edab6123Sjeremylt   "metadata": {},
81edab6123Sjeremylt   "outputs": [],
82edab6123Sjeremylt   "source": [
83edab6123Sjeremylt    "qf = ceed.IdentityQFunction(1, libceed.EVAL_INTERP, libceed.EVAL_INTERP)\n",
84edab6123Sjeremylt    "\n",
85edab6123Sjeremylt    "q = 8\n",
86edab6123Sjeremylt    "\n",
87edab6123Sjeremylt    "u_array = np.zeros(q, dtype=\"float64\")\n",
88edab6123Sjeremylt    "for i in range(q):\n",
89edab6123Sjeremylt    "  u_array[i] = i*i\n",
90edab6123Sjeremylt    "\n",
91edab6123Sjeremylt    "u = ceed.Vector(q)\n",
92edab6123Sjeremylt    "u.set_array(u_array, cmode=libceed.USE_POINTER)\n",
93edab6123Sjeremylt    "v = ceed.Vector(q)\n",
94edab6123Sjeremylt    "v.set_value(0)\n",
95edab6123Sjeremylt    "\n",
96edab6123Sjeremylt    "inputs = [ u ]\n",
97edab6123Sjeremylt    "outputs = [ v ]\n",
98edab6123Sjeremylt    "qf.apply(q, inputs, outputs)\n",
99edab6123Sjeremylt    "\n",
100edab6123Sjeremylt    "print('v =', v)"
101edab6123Sjeremylt   ]
102edab6123Sjeremylt  },
103edab6123Sjeremylt  {
104edab6123Sjeremylt   "cell_type": "markdown",
105edab6123Sjeremylt   "metadata": {},
106edab6123Sjeremylt   "source": [
107edab6123Sjeremylt    "* In the following example, we create and evaluate a QFunction (for the mass operator in 1D) from the gallery of available built-in QFunctions in libCEED."
108edab6123Sjeremylt   ]
109edab6123Sjeremylt  },
110edab6123Sjeremylt  {
111edab6123Sjeremylt   "cell_type": "code",
112edab6123Sjeremylt   "execution_count": null,
113edab6123Sjeremylt   "metadata": {},
114edab6123Sjeremylt   "outputs": [],
115edab6123Sjeremylt   "source": [
116edab6123Sjeremylt    "qf_setup = ceed.QFunctionByName(\"Mass1DBuild\")\n",
117edab6123Sjeremylt    "qf_mass = ceed.QFunctionByName(\"MassApply\")\n",
118edab6123Sjeremylt    "\n",
119edab6123Sjeremylt    "q = 8\n",
120edab6123Sjeremylt    "\n",
121edab6123Sjeremylt    "j_array = np.zeros(q, dtype=\"float64\")\n",
122edab6123Sjeremylt    "w_array = np.zeros(q, dtype=\"float64\")\n",
123edab6123Sjeremylt    "u_array = np.zeros(q, dtype=\"float64\")\n",
124edab6123Sjeremylt    "v_true  = np.zeros(q, dtype=\"float64\")\n",
125edab6123Sjeremylt    "for i in range(q):\n",
126edab6123Sjeremylt    "  x = 2.*i/(q-1) - 1\n",
127edab6123Sjeremylt    "  j_array[i] = 1\n",
128edab6123Sjeremylt    "  w_array[i] = 1 - x*x\n",
129edab6123Sjeremylt    "  u_array[i] = 2 + 3*x + 5*x*x\n",
130edab6123Sjeremylt    "  v_true[i]  = w_array[i] * u_array[i]\n",
131edab6123Sjeremylt    "\n",
132edab6123Sjeremylt    "j = ceed.Vector(q)\n",
133edab6123Sjeremylt    "j.set_array(j_array, cmode=libceed.USE_POINTER)\n",
134edab6123Sjeremylt    "w = ceed.Vector(q)\n",
135edab6123Sjeremylt    "w.set_array(w_array, cmode=libceed.USE_POINTER)\n",
136edab6123Sjeremylt    "u = ceed.Vector(q)\n",
137edab6123Sjeremylt    "u.set_array(u_array, cmode=libceed.USE_POINTER)\n",
138edab6123Sjeremylt    "v = ceed.Vector(q)\n",
139edab6123Sjeremylt    "v.set_value(0)\n",
140edab6123Sjeremylt    "qdata = ceed.Vector(q)\n",
141edab6123Sjeremylt    "qdata.set_value(0)\n",
142edab6123Sjeremylt    "\n",
143edab6123Sjeremylt    "inputs = [ j, w ]\n",
144edab6123Sjeremylt    "outputs = [ qdata ]\n",
145edab6123Sjeremylt    "qf_setup.apply(q, inputs, outputs)\n",
146edab6123Sjeremylt    "\n",
147edab6123Sjeremylt    "inputs = [ w, u ]\n",
148edab6123Sjeremylt    "outputs = [ v ]\n",
149edab6123Sjeremylt    "qf_mass.apply(q, inputs, outputs)\n",
150edab6123Sjeremylt    "\n",
151edab6123Sjeremylt    "print('v =', v)"
152edab6123Sjeremylt   ]
153edab6123Sjeremylt  },
154edab6123Sjeremylt  {
155edab6123Sjeremylt   "cell_type": "markdown",
156edab6123Sjeremylt   "metadata": {},
157edab6123Sjeremylt   "source": [
158edab6123Sjeremylt    "* In the following example, we create and evaluate a built-in identity QFunction 3 fields per quadrature point."
159edab6123Sjeremylt   ]
160edab6123Sjeremylt  },
161edab6123Sjeremylt  {
162edab6123Sjeremylt   "cell_type": "code",
163edab6123Sjeremylt   "execution_count": null,
164edab6123Sjeremylt   "metadata": {},
165edab6123Sjeremylt   "outputs": [],
166edab6123Sjeremylt   "source": [
167edab6123Sjeremylt    "fields = 3\n",
168edab6123Sjeremylt    "\n",
169edab6123Sjeremylt    "qf = ceed.IdentityQFunction(fields, libceed.EVAL_INTERP, libceed.EVAL_INTERP)\n",
170edab6123Sjeremylt    "\n",
171edab6123Sjeremylt    "q = 8\n",
172edab6123Sjeremylt    "\n",
173edab6123Sjeremylt    "u_array = np.zeros(q*fields, dtype=\"float64\")\n",
174edab6123Sjeremylt    "for i in range(q*fields):\n",
175edab6123Sjeremylt    "  u_array[i] = i*i\n",
176edab6123Sjeremylt    "\n",
177edab6123Sjeremylt    "u = ceed.Vector(q*fields)\n",
178edab6123Sjeremylt    "u.set_array(u_array, cmode=libceed.USE_POINTER)\n",
179edab6123Sjeremylt    "v = ceed.Vector(q*fields)\n",
180edab6123Sjeremylt    "v.set_value(0)\n",
181edab6123Sjeremylt    "\n",
182edab6123Sjeremylt    "inputs = [ u ]\n",
183edab6123Sjeremylt    "outputs = [ v ]\n",
184edab6123Sjeremylt    "qf.apply(q, inputs, outputs)\n",
185edab6123Sjeremylt    "\n",
186edab6123Sjeremylt    "print('v =', v)"
187edab6123Sjeremylt   ]
188edab6123Sjeremylt  }
189edab6123Sjeremylt ],
190edab6123Sjeremylt "metadata": {
191edab6123Sjeremylt  "kernelspec": {
192*235e4399SJeremy L Thompson   "display_name": "Python 3 (ipykernel)",
193edab6123Sjeremylt   "language": "python",
194edab6123Sjeremylt   "name": "python3"
195edab6123Sjeremylt  },
196edab6123Sjeremylt  "language_info": {
197edab6123Sjeremylt   "codemirror_mode": {
198edab6123Sjeremylt    "name": "ipython",
199edab6123Sjeremylt    "version": 3
200edab6123Sjeremylt   },
201edab6123Sjeremylt   "file_extension": ".py",
202edab6123Sjeremylt   "mimetype": "text/x-python",
203edab6123Sjeremylt   "name": "python",
204edab6123Sjeremylt   "nbconvert_exporter": "python",
205edab6123Sjeremylt   "pygments_lexer": "ipython3",
206*235e4399SJeremy L Thompson   "version": "3.13.2"
207edab6123Sjeremylt  }
208edab6123Sjeremylt },
209edab6123Sjeremylt "nbformat": 4,
210edab6123Sjeremylt "nbformat_minor": 4
211edab6123Sjeremylt}
212