xref: /libCEED/examples/python/tutorial-3-basis.ipynb (revision edab612303ce5c583510db4ec9520fcda426d3c1)
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    "## CeedBasis\n",
37*edab6123Sjeremylt    "\n",
38*edab6123Sjeremylt    "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.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    "First we declare some auxiliary functions needed in the following examples"
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    "%matplotlib inline\n",
55*edab6123Sjeremylt    "import numpy as np\n",
56*edab6123Sjeremylt    "import matplotlib.pyplot as plt\n",
57*edab6123Sjeremylt    "plt.style.use('ggplot')\n",
58*edab6123Sjeremylt    "\n",
59*edab6123Sjeremylt    "def eval(dim, x):\n",
60*edab6123Sjeremylt    "  result, center = 1, 0.1\n",
61*edab6123Sjeremylt    "  for d in range(dim):\n",
62*edab6123Sjeremylt    "    result *= np.tanh(x[d] - center)\n",
63*edab6123Sjeremylt    "    center += 0.1\n",
64*edab6123Sjeremylt    "  return result\n",
65*edab6123Sjeremylt    "\n",
66*edab6123Sjeremylt    "def feval(x1, x2):\n",
67*edab6123Sjeremylt    "  return x1*x1 + x2*x2 + x1*x2 + 1\n",
68*edab6123Sjeremylt    "\n",
69*edab6123Sjeremylt    "def dfeval(x1, x2):\n",
70*edab6123Sjeremylt    "  return 2*x1 + x2"
71*edab6123Sjeremylt   ]
72*edab6123Sjeremylt  },
73*edab6123Sjeremylt  {
74*edab6123Sjeremylt   "cell_type": "markdown",
75*edab6123Sjeremylt   "metadata": {},
76*edab6123Sjeremylt   "source": [
77*edab6123Sjeremylt    "## $H^1$ Lagrange bases in 1D\n",
78*edab6123Sjeremylt    "\n",
79*edab6123Sjeremylt    "The Lagrange interpolation nodes are at the Gauss-Lobatto points, so interpolation to Gauss-Lobatto quadrature points is the identity."
80*edab6123Sjeremylt   ]
81*edab6123Sjeremylt  },
82*edab6123Sjeremylt  {
83*edab6123Sjeremylt   "cell_type": "code",
84*edab6123Sjeremylt   "execution_count": null,
85*edab6123Sjeremylt   "metadata": {},
86*edab6123Sjeremylt   "outputs": [],
87*edab6123Sjeremylt   "source": [
88*edab6123Sjeremylt    "import libceed\n",
89*edab6123Sjeremylt    "\n",
90*edab6123Sjeremylt    "ceed = libceed.Ceed()\n",
91*edab6123Sjeremylt    "\n",
92*edab6123Sjeremylt    "b = ceed.BasisTensorH1Lagrange(\n",
93*edab6123Sjeremylt    "    dim=1,   # topological dimension\n",
94*edab6123Sjeremylt    "    ncomp=1, # number of components\n",
95*edab6123Sjeremylt    "    P=4,     # number of basis functions (nodes) per dimension\n",
96*edab6123Sjeremylt    "    Q=4,     # number of quadrature points per dimension\n",
97*edab6123Sjeremylt    "    qmode=libceed.GAUSS_LOBATTO)\n",
98*edab6123Sjeremylt    "print(b)"
99*edab6123Sjeremylt   ]
100*edab6123Sjeremylt  },
101*edab6123Sjeremylt  {
102*edab6123Sjeremylt   "cell_type": "markdown",
103*edab6123Sjeremylt   "metadata": {},
104*edab6123Sjeremylt   "source": [
105*edab6123Sjeremylt    "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."
106*edab6123Sjeremylt   ]
107*edab6123Sjeremylt  },
108*edab6123Sjeremylt  {
109*edab6123Sjeremylt   "cell_type": "code",
110*edab6123Sjeremylt   "execution_count": null,
111*edab6123Sjeremylt   "metadata": {},
112*edab6123Sjeremylt   "outputs": [],
113*edab6123Sjeremylt   "source": [
114*edab6123Sjeremylt    "P = b.get_num_nodes()\n",
115*edab6123Sjeremylt    "nviz = 50\n",
116*edab6123Sjeremylt    "bviz = ceed.BasisTensorH1Lagrange(1, 1, P, nviz, libceed.GAUSS_LOBATTO)\n",
117*edab6123Sjeremylt    "\n",
118*edab6123Sjeremylt    "# Construct P \"elements\" with one node activated\n",
119*edab6123Sjeremylt    "I = ceed.Vector(P * P)\n",
120*edab6123Sjeremylt    "with I.array(P, P) as x:\n",
121*edab6123Sjeremylt    "    x[...] = np.eye(P)\n",
122*edab6123Sjeremylt    "\n",
123*edab6123Sjeremylt    "Bvander = ceed.Vector(P * nviz)\n",
124*edab6123Sjeremylt    "bviz.apply(4, libceed.EVAL_INTERP, I, Bvander)\n",
125*edab6123Sjeremylt    "\n",
126*edab6123Sjeremylt    "qviz, _weight = ceed.lobatto_quadrature(nviz)\n",
127*edab6123Sjeremylt    "with Bvander.array_read(nviz, P) as B:\n",
128*edab6123Sjeremylt    "    plt.plot(qviz, B)\n",
129*edab6123Sjeremylt    "\n",
130*edab6123Sjeremylt    "# Mark tho Lobatto nodes\n",
131*edab6123Sjeremylt    "qb, _weight = ceed.lobatto_quadrature(P)\n",
132*edab6123Sjeremylt    "plt.plot(qb, 0*qb, 'ok');"
133*edab6123Sjeremylt   ]
134*edab6123Sjeremylt  },
135*edab6123Sjeremylt  {
136*edab6123Sjeremylt   "cell_type": "markdown",
137*edab6123Sjeremylt   "metadata": {},
138*edab6123Sjeremylt   "source": [
139*edab6123Sjeremylt    "In contrast, the Gauss quadrature points are not collocated, and thus all basis functions are generally nonzero at every quadrature point."
140*edab6123Sjeremylt   ]
141*edab6123Sjeremylt  },
142*edab6123Sjeremylt  {
143*edab6123Sjeremylt   "cell_type": "code",
144*edab6123Sjeremylt   "execution_count": null,
145*edab6123Sjeremylt   "metadata": {},
146*edab6123Sjeremylt   "outputs": [],
147*edab6123Sjeremylt   "source": [
148*edab6123Sjeremylt    "b = ceed.BasisTensorH1Lagrange(1, 1, 4, 4, libceed.GAUSS)\n",
149*edab6123Sjeremylt    "print(b)\n",
150*edab6123Sjeremylt    "\n",
151*edab6123Sjeremylt    "with Bvander.array_read(nviz, P) as B:\n",
152*edab6123Sjeremylt    "    plt.plot(qviz, B)\n",
153*edab6123Sjeremylt    "# Mark tho Gauss quadrature points\n",
154*edab6123Sjeremylt    "qb, _weight = ceed.gauss_quadrature(P)\n",
155*edab6123Sjeremylt    "plt.plot(qb, 0*qb, 'ok');"
156*edab6123Sjeremylt   ]
157*edab6123Sjeremylt  },
158*edab6123Sjeremylt  {
159*edab6123Sjeremylt   "cell_type": "markdown",
160*edab6123Sjeremylt   "metadata": {},
161*edab6123Sjeremylt   "source": [
162*edab6123Sjeremylt    "Although the underlying functions are not an intrinsic property of a `libceed.Basis` in libCEED, the sizes are.\n",
163*edab6123Sjeremylt    "Here, we create a 3D tensor product element with more quadrature points than Lagrange interpolation nodes."
164*edab6123Sjeremylt   ]
165*edab6123Sjeremylt  },
166*edab6123Sjeremylt  {
167*edab6123Sjeremylt   "cell_type": "code",
168*edab6123Sjeremylt   "execution_count": null,
169*edab6123Sjeremylt   "metadata": {},
170*edab6123Sjeremylt   "outputs": [],
171*edab6123Sjeremylt   "source": [
172*edab6123Sjeremylt    "b = ceed.BasisTensorH1Lagrange(3, 1, 4, 5, libceed.GAUSS_LOBATTO)\n",
173*edab6123Sjeremylt    "\n",
174*edab6123Sjeremylt    "p = b.get_num_nodes()\n",
175*edab6123Sjeremylt    "print('p =', p)\n",
176*edab6123Sjeremylt    "\n",
177*edab6123Sjeremylt    "q = b.get_num_quadrature_points()\n",
178*edab6123Sjeremylt    "print('q =', q)"
179*edab6123Sjeremylt   ]
180*edab6123Sjeremylt  },
181*edab6123Sjeremylt  {
182*edab6123Sjeremylt   "cell_type": "markdown",
183*edab6123Sjeremylt   "metadata": {},
184*edab6123Sjeremylt   "source": [
185*edab6123Sjeremylt    "* In the following example, we demonstrate the application of an interpolatory basis in multiple dimensions"
186*edab6123Sjeremylt   ]
187*edab6123Sjeremylt  },
188*edab6123Sjeremylt  {
189*edab6123Sjeremylt   "cell_type": "code",
190*edab6123Sjeremylt   "execution_count": null,
191*edab6123Sjeremylt   "metadata": {},
192*edab6123Sjeremylt   "outputs": [],
193*edab6123Sjeremylt   "source": [
194*edab6123Sjeremylt    "for dim in range(1, 4):\n",
195*edab6123Sjeremylt    "  Q = 4\n",
196*edab6123Sjeremylt    "  Qdim = Q**dim\n",
197*edab6123Sjeremylt    "  Xdim = 2**dim\n",
198*edab6123Sjeremylt    "  x = np.empty(Xdim*dim, dtype=\"float64\")\n",
199*edab6123Sjeremylt    "  uq = np.empty(Qdim, dtype=\"float64\")\n",
200*edab6123Sjeremylt    "\n",
201*edab6123Sjeremylt    "  for d in range(dim):\n",
202*edab6123Sjeremylt    "    for i in range(Xdim):\n",
203*edab6123Sjeremylt    "      x[d*Xdim + i] = 1 if (i % (2**(dim-d))) // (2**(dim-d-1)) else -1\n",
204*edab6123Sjeremylt    "\n",
205*edab6123Sjeremylt    "  X = ceed.Vector(Xdim*dim)\n",
206*edab6123Sjeremylt    "  X.set_array(x, cmode=libceed.USE_POINTER)\n",
207*edab6123Sjeremylt    "  Xq = ceed.Vector(Qdim*dim)\n",
208*edab6123Sjeremylt    "  Xq.set_value(0)\n",
209*edab6123Sjeremylt    "  U = ceed.Vector(Qdim)\n",
210*edab6123Sjeremylt    "  U.set_value(0)\n",
211*edab6123Sjeremylt    "  Uq = ceed.Vector(Qdim)\n",
212*edab6123Sjeremylt    "\n",
213*edab6123Sjeremylt    "  bxl = ceed.BasisTensorH1Lagrange(dim, dim, 2, Q, libceed.GAUSS_LOBATTO)\n",
214*edab6123Sjeremylt    "  bul = ceed.BasisTensorH1Lagrange(dim, 1, Q, Q, libceed.GAUSS_LOBATTO)\n",
215*edab6123Sjeremylt    "\n",
216*edab6123Sjeremylt    "  bxl.apply(1, libceed.EVAL_INTERP, X, Xq)\n",
217*edab6123Sjeremylt    "\n",
218*edab6123Sjeremylt    "  with Xq.array_read() as xq:\n",
219*edab6123Sjeremylt    "    for i in range(Qdim):\n",
220*edab6123Sjeremylt    "      xx = np.empty(dim, dtype=\"float64\")\n",
221*edab6123Sjeremylt    "      for d in range(dim):\n",
222*edab6123Sjeremylt    "        xx[d] = xq[d*Qdim + i]\n",
223*edab6123Sjeremylt    "      uq[i] = eval(dim, xx)\n",
224*edab6123Sjeremylt    "\n",
225*edab6123Sjeremylt    "  Uq.set_array(uq, cmode=libceed.USE_POINTER)\n",
226*edab6123Sjeremylt    "\n",
227*edab6123Sjeremylt    "  # This operation is the identity because the quadrature is collocated\n",
228*edab6123Sjeremylt    "  bul.T.apply(1, libceed.EVAL_INTERP, Uq, U)\n",
229*edab6123Sjeremylt    "\n",
230*edab6123Sjeremylt    "  bxg = ceed.BasisTensorH1Lagrange(dim, dim, 2, Q, libceed.GAUSS)\n",
231*edab6123Sjeremylt    "  bug = ceed.BasisTensorH1Lagrange(dim, 1, Q, Q, libceed.GAUSS)\n",
232*edab6123Sjeremylt    "\n",
233*edab6123Sjeremylt    "  bxg.apply(1, libceed.EVAL_INTERP, X, Xq)\n",
234*edab6123Sjeremylt    "  bug.apply(1, libceed.EVAL_INTERP, U, Uq)\n",
235*edab6123Sjeremylt    "\n",
236*edab6123Sjeremylt    "  with Xq.array_read() as xq, Uq.array_read() as u:\n",
237*edab6123Sjeremylt    "    #print('xq =', xq)\n",
238*edab6123Sjeremylt    "    #print('u =', u)\n",
239*edab6123Sjeremylt    "    if dim == 2:\n",
240*edab6123Sjeremylt    "        # Default ordering is contiguous in x direction, but\n",
241*edab6123Sjeremylt    "        # pyplot expects meshgrid convention, which is transposed.\n",
242*edab6123Sjeremylt    "        x, y = xq.reshape(2, Q, Q).transpose(0, 2, 1)\n",
243*edab6123Sjeremylt    "        plt.scatter(x, y, c=np.array(u).reshape(Q, Q))\n",
244*edab6123Sjeremylt    "        plt.xlim(-1, 1)\n",
245*edab6123Sjeremylt    "        plt.ylim(-1, 1)\n",
246*edab6123Sjeremylt    "        plt.colorbar(label='u')"
247*edab6123Sjeremylt   ]
248*edab6123Sjeremylt  },
249*edab6123Sjeremylt  {
250*edab6123Sjeremylt   "cell_type": "markdown",
251*edab6123Sjeremylt   "metadata": {},
252*edab6123Sjeremylt   "source": [
253*edab6123Sjeremylt    "* In the following example, we demonstrate the application of the gradient of the shape functions in multiple dimensions"
254*edab6123Sjeremylt   ]
255*edab6123Sjeremylt  },
256*edab6123Sjeremylt  {
257*edab6123Sjeremylt   "cell_type": "code",
258*edab6123Sjeremylt   "execution_count": null,
259*edab6123Sjeremylt   "metadata": {},
260*edab6123Sjeremylt   "outputs": [],
261*edab6123Sjeremylt   "source": [
262*edab6123Sjeremylt    "for dim in range (1, 4):\n",
263*edab6123Sjeremylt    "  P, Q = 8, 10\n",
264*edab6123Sjeremylt    "  Pdim = P**dim\n",
265*edab6123Sjeremylt    "  Qdim = Q**dim\n",
266*edab6123Sjeremylt    "  Xdim = 2**dim\n",
267*edab6123Sjeremylt    "  sum1 = sum2 = 0\n",
268*edab6123Sjeremylt    "  x = np.empty(Xdim*dim, dtype=\"float64\")\n",
269*edab6123Sjeremylt    "  u = np.empty(Pdim, dtype=\"float64\")\n",
270*edab6123Sjeremylt    "\n",
271*edab6123Sjeremylt    "  for d in range(dim):\n",
272*edab6123Sjeremylt    "    for i in range(Xdim):\n",
273*edab6123Sjeremylt    "      x[d*Xdim + i] = 1 if (i % (2**(dim-d))) // (2**(dim-d-1)) else -1\n",
274*edab6123Sjeremylt    "\n",
275*edab6123Sjeremylt    "  X = ceed.Vector(Xdim*dim)\n",
276*edab6123Sjeremylt    "  X.set_array(x, cmode=libceed.USE_POINTER)\n",
277*edab6123Sjeremylt    "  Xq = ceed.Vector(Pdim*dim)\n",
278*edab6123Sjeremylt    "  Xq.set_value(0)\n",
279*edab6123Sjeremylt    "  U = ceed.Vector(Pdim)\n",
280*edab6123Sjeremylt    "  Uq = ceed.Vector(Qdim*dim)\n",
281*edab6123Sjeremylt    "  Uq.set_value(0)\n",
282*edab6123Sjeremylt    "  Ones = ceed.Vector(Qdim*dim)\n",
283*edab6123Sjeremylt    "  Ones.set_value(1)\n",
284*edab6123Sjeremylt    "  Gtposeones = ceed.Vector(Pdim)\n",
285*edab6123Sjeremylt    "  Gtposeones.set_value(0)\n",
286*edab6123Sjeremylt    "\n",
287*edab6123Sjeremylt    "  # Get function values at quadrature points\n",
288*edab6123Sjeremylt    "  bxl = ceed.BasisTensorH1Lagrange(dim, dim, 2, P, libceed.GAUSS_LOBATTO)\n",
289*edab6123Sjeremylt    "  bxl.apply(1, libceed.EVAL_INTERP, X, Xq)\n",
290*edab6123Sjeremylt    "\n",
291*edab6123Sjeremylt    "  with Xq.array_read() as xq:\n",
292*edab6123Sjeremylt    "    for i in range(Pdim):\n",
293*edab6123Sjeremylt    "      xx = np.empty(dim, dtype=\"float64\")\n",
294*edab6123Sjeremylt    "      for d in range(dim):\n",
295*edab6123Sjeremylt    "        xx[d] = xq[d*Pdim + i]\n",
296*edab6123Sjeremylt    "      u[i] = eval(dim, xx)\n",
297*edab6123Sjeremylt    "\n",
298*edab6123Sjeremylt    "  U.set_array(u, cmode=libceed.USE_POINTER)\n",
299*edab6123Sjeremylt    "\n",
300*edab6123Sjeremylt    "  # Calculate G u at quadrature points, G' * 1 at dofs\n",
301*edab6123Sjeremylt    "  bug = ceed.BasisTensorH1Lagrange(dim, 1, P, Q, libceed.GAUSS)\n",
302*edab6123Sjeremylt    "  bug.apply(1, libceed.EVAL_GRAD, U, Uq)\n",
303*edab6123Sjeremylt    "  bug.T.apply(1, libceed.EVAL_GRAD, Ones, Gtposeones)\n",
304*edab6123Sjeremylt    "\n",
305*edab6123Sjeremylt    "  # Check if 1' * G * u = u' * (G' * 1)\n",
306*edab6123Sjeremylt    "  with Gtposeones.array_read() as gtposeones, Uq.array_read() as uq:\n",
307*edab6123Sjeremylt    "    for i in range(Pdim):\n",
308*edab6123Sjeremylt    "      sum1 += gtposeones[i]*u[i]\n",
309*edab6123Sjeremylt    "    for i in range(dim*Qdim):\n",
310*edab6123Sjeremylt    "      sum2 += uq[i]\n",
311*edab6123Sjeremylt    "\n",
312*edab6123Sjeremylt    "  # Check that (1' * G * u - u' * (G' * 1)) is numerically zero\n",
313*edab6123Sjeremylt    "  print('1T * G * u - uT * (GT * 1) =', np.abs(sum1 - sum2))"
314*edab6123Sjeremylt   ]
315*edab6123Sjeremylt  },
316*edab6123Sjeremylt  {
317*edab6123Sjeremylt   "cell_type": "markdown",
318*edab6123Sjeremylt   "metadata": {},
319*edab6123Sjeremylt   "source": [
320*edab6123Sjeremylt    "### Advanced topics"
321*edab6123Sjeremylt   ]
322*edab6123Sjeremylt  },
323*edab6123Sjeremylt  {
324*edab6123Sjeremylt   "cell_type": "markdown",
325*edab6123Sjeremylt   "metadata": {},
326*edab6123Sjeremylt   "source": [
327*edab6123Sjeremylt    "* In the following example, we demonstrate the QR factorization of a basis matrix.\n",
328*edab6123Sjeremylt    "The representation is similar to LAPACK's [`dgeqrf`](https://www.netlib.org/lapack/explore-html/dd/d9a/group__double_g_ecomputational_ga3766ea903391b5cf9008132f7440ec7b.html#ga3766ea903391b5cf9008132f7440ec7b), with elementary reflectors in the lower triangular block, scaled by `tau`."
329*edab6123Sjeremylt   ]
330*edab6123Sjeremylt  },
331*edab6123Sjeremylt  {
332*edab6123Sjeremylt   "cell_type": "code",
333*edab6123Sjeremylt   "execution_count": null,
334*edab6123Sjeremylt   "metadata": {},
335*edab6123Sjeremylt   "outputs": [],
336*edab6123Sjeremylt   "source": [
337*edab6123Sjeremylt    "qr = np.array([1, -1, 4, 1, 4, -2, 1, 4, 2, 1, -1, 0], dtype=\"float64\")\n",
338*edab6123Sjeremylt    "tau = np.empty(3, dtype=\"float64\")\n",
339*edab6123Sjeremylt    "\n",
340*edab6123Sjeremylt    "qr, tau = ceed.qr_factorization(qr, tau, 4, 3)\n",
341*edab6123Sjeremylt    "\n",
342*edab6123Sjeremylt    "print('qr =')\n",
343*edab6123Sjeremylt    "print(qr.reshape(4, 3))\n",
344*edab6123Sjeremylt    "\n",
345*edab6123Sjeremylt    "print('tau =')\n",
346*edab6123Sjeremylt    "print(tau)"
347*edab6123Sjeremylt   ]
348*edab6123Sjeremylt  },
349*edab6123Sjeremylt  {
350*edab6123Sjeremylt   "cell_type": "markdown",
351*edab6123Sjeremylt   "metadata": {},
352*edab6123Sjeremylt   "source": [
353*edab6123Sjeremylt    "* In the following example, we demonstrate the symmetric Schur decomposition of a basis matrix"
354*edab6123Sjeremylt   ]
355*edab6123Sjeremylt  },
356*edab6123Sjeremylt  {
357*edab6123Sjeremylt   "cell_type": "code",
358*edab6123Sjeremylt   "execution_count": null,
359*edab6123Sjeremylt   "metadata": {},
360*edab6123Sjeremylt   "outputs": [],
361*edab6123Sjeremylt   "source": [
362*edab6123Sjeremylt    "A = np.array([0.19996678, 0.0745459, -0.07448852, 0.0332866,\n",
363*edab6123Sjeremylt    "              0.0745459, 1., 0.16666509, -0.07448852,\n",
364*edab6123Sjeremylt    "              -0.07448852, 0.16666509, 1., 0.0745459,\n",
365*edab6123Sjeremylt    "              0.0332866, -0.07448852, 0.0745459, 0.19996678], dtype=\"float64\")\n",
366*edab6123Sjeremylt    "\n",
367*edab6123Sjeremylt    "lam = ceed.symmetric_schur_decomposition(A, 4)\n",
368*edab6123Sjeremylt    "\n",
369*edab6123Sjeremylt    "print(\"Q =\")\n",
370*edab6123Sjeremylt    "for i in range(4):\n",
371*edab6123Sjeremylt    "  for j in range(4):\n",
372*edab6123Sjeremylt    "    if A[j+4*i] <= 1E-14 and A[j+4*i] >= -1E-14:\n",
373*edab6123Sjeremylt    "       A[j+4*i] = 0\n",
374*edab6123Sjeremylt    "    print(\"%12.8f\"%A[j+4*i])\n",
375*edab6123Sjeremylt    "\n",
376*edab6123Sjeremylt    "print(\"lambda =\")\n",
377*edab6123Sjeremylt    "for i in range(4):\n",
378*edab6123Sjeremylt    "  if lam[i] <= 1E-14 and lam[i] >= -1E-14:\n",
379*edab6123Sjeremylt    "    lam[i] = 0\n",
380*edab6123Sjeremylt    "  print(\"%12.8f\"%lam[i])"
381*edab6123Sjeremylt   ]
382*edab6123Sjeremylt  },
383*edab6123Sjeremylt  {
384*edab6123Sjeremylt   "cell_type": "markdown",
385*edab6123Sjeremylt   "metadata": {},
386*edab6123Sjeremylt   "source": [
387*edab6123Sjeremylt    "* In the following example, we demonstrate the simultaneous diagonalization of a basis matrix"
388*edab6123Sjeremylt   ]
389*edab6123Sjeremylt  },
390*edab6123Sjeremylt  {
391*edab6123Sjeremylt   "cell_type": "code",
392*edab6123Sjeremylt   "execution_count": null,
393*edab6123Sjeremylt   "metadata": {},
394*edab6123Sjeremylt   "outputs": [],
395*edab6123Sjeremylt   "source": [
396*edab6123Sjeremylt    "M = np.array([0.19996678, 0.0745459, -0.07448852, 0.0332866,\n",
397*edab6123Sjeremylt    "              0.0745459, 1., 0.16666509, -0.07448852,\n",
398*edab6123Sjeremylt    "              -0.07448852, 0.16666509, 1., 0.0745459,\n",
399*edab6123Sjeremylt    "              0.0332866, -0.07448852, 0.0745459, 0.19996678], dtype=\"float64\")\n",
400*edab6123Sjeremylt    "K = np.array([3.03344425, -3.41501767, 0.49824435, -0.11667092,\n",
401*edab6123Sjeremylt    "              -3.41501767, 5.83354662, -2.9167733, 0.49824435,\n",
402*edab6123Sjeremylt    "              0.49824435, -2.9167733, 5.83354662, -3.41501767,\n",
403*edab6123Sjeremylt    "              -0.11667092, 0.49824435, -3.41501767, 3.03344425], dtype=\"float64\")\n",
404*edab6123Sjeremylt    "\n",
405*edab6123Sjeremylt    "x, lam = ceed.simultaneous_diagonalization(K, M, 4)\n",
406*edab6123Sjeremylt    "\n",
407*edab6123Sjeremylt    "print(\"x =\")\n",
408*edab6123Sjeremylt    "for i in range(4):\n",
409*edab6123Sjeremylt    "  for j in range(4):\n",
410*edab6123Sjeremylt    "    if x[j+4*i] <= 1E-14 and x[j+4*i] >= -1E-14:\n",
411*edab6123Sjeremylt    "      x[j+4*i] = 0\n",
412*edab6123Sjeremylt    "    print(\"%12.8f\"%x[j+4*i])\n",
413*edab6123Sjeremylt    "\n",
414*edab6123Sjeremylt    "print(\"lambda =\")\n",
415*edab6123Sjeremylt    "for i in range(4):\n",
416*edab6123Sjeremylt    "  if lam[i] <= 1E-14 and lam[i] >= -1E-14:\n",
417*edab6123Sjeremylt    "    lam[i] = 0\n",
418*edab6123Sjeremylt    "  print(\"%12.8f\"%lam[i])"
419*edab6123Sjeremylt   ]
420*edab6123Sjeremylt  }
421*edab6123Sjeremylt ],
422*edab6123Sjeremylt "metadata": {
423*edab6123Sjeremylt  "kernelspec": {
424*edab6123Sjeremylt   "display_name": "Python 3",
425*edab6123Sjeremylt   "language": "python",
426*edab6123Sjeremylt   "name": "python3"
427*edab6123Sjeremylt  },
428*edab6123Sjeremylt  "language_info": {
429*edab6123Sjeremylt   "codemirror_mode": {
430*edab6123Sjeremylt    "name": "ipython",
431*edab6123Sjeremylt    "version": 3
432*edab6123Sjeremylt   },
433*edab6123Sjeremylt   "file_extension": ".py",
434*edab6123Sjeremylt   "mimetype": "text/x-python",
435*edab6123Sjeremylt   "name": "python",
436*edab6123Sjeremylt   "nbconvert_exporter": "python",
437*edab6123Sjeremylt   "pygments_lexer": "ipython3",
438*edab6123Sjeremylt   "version": "3.8.5"
439*edab6123Sjeremylt  }
440*edab6123Sjeremylt },
441*edab6123Sjeremylt "nbformat": 4,
442*edab6123Sjeremylt "nbformat_minor": 4
443*edab6123Sjeremylt}
444