1. Solving the Poisson Equation#

The finite element method is a numerical method for solving partial differential equations approximately. A typical example is the Poisson equation:

\[ -\Delta u(x) = f(x) \quad \forall \, x \in \Omega \]

The right hand side \(f\) is a given function, and we search for the solution \(u\). The domain \(\Omega\) is a subset of \({\mathbb R}^d\). The Poisson equation is a model for many physical phenomena:

  • f can be a heat source distribution, and u is the temperature

  • f can be an electric charge distribution, and u is the electrostatic potential

To select a unique solution \(u\) we have to specify boundary conditions, for example homogeneous Dirichlet boundary conditions

\[ u(x) = 0 \quad \forall \, x \in \partial \Omega \]

1.1. Weak formulation#

We derive the weak formulation (also called variational formulation) of the Poisson equation. The formulation above is called the strong form. The weak form is the starting point for the finite element discretization method.

First, we multiply the Poisson equation by a so called test function. It is an arbitrary function, some restriction will be given later as needed. We multiply the strong form by the function v:

\[ - \Delta u(x) v(x) = f(x) v(x) \qquad \forall x \in \Omega \]

We integrate over the domain \(\Omega\):

\[ - \int_\Omega \Delta u(x) v(x) dx = \int_\Omega f(x) v(x) dx \]

From Gauss’ Theorem \(\int_\Omega \operatorname{div} b \, dx = \int_{\partial \Omega} n\cdot b \, ds\) applied to the vector field \(b = \nabla u v\) we obtain

\[ \int_{\partial \Omega} n \nabla u \, v \, ds = \int_\Omega \operatorname{div} (\nabla u \, v) \, dx = \int_{\Omega} \Delta u v + \nabla u \nabla v \, dx. \]

This allows us to rewrite the left hand side such that

\[ \int_\Omega \nabla u \nabla v \, dx - \int_{\partial \Omega} \frac{\partial u}{\partial n} v \, ds = \int_\Omega f v \, dx \]

In the case of Dirichlet boundary conditions we allow only test-functions \(v\) such that \(v(x) = 0\) on the boundary \(\partial \Omega\).

We have derived the weak form: find \(u\) such that \(u = 0\) on \(\partial \Omega\) and

\[ \int_\Omega \nabla u \nabla v \, dx = \int_\Omega f v \, dx \]

holds true for all test-functions \(v\) with \(v = 0\) on \(\partial \Omega\). Note that the weak formulation needs only first order derivatives of \(u\) and \(v\), in contrast to the strong form which requires second order derivatives of \(u\).

1.2. The Sobolev space \(H^1\), linear and bilinear forms#

The proper space to search for the solution is the so called Sobolev space

\[ H^1(\Omega) := \{ u \in L_2(\Omega) : \nabla u \in L_2(\Omega)^d \} \]

The super-script \(1\) indicates that we want to have first order derivatives in \(L_2\). We just note that the derivative is understood in weak sense, which is well defined for functions with kinks. The vector space \(H^1\) comes with the norm

\[ \| u \|_{H^1}^2 := \| u \|_{L_2}^2 + \| \nabla u \|_{L_2}^2 \]

and the inner product

\[ (u,v)_{H^1} = (u,v)_{L_2} + (\nabla u, \nabla v)_{L_2}. \]

It is a complete space with an inner product what is called a Hilbert space.

It does not make sense to take boundary values of \(L_2\)-functions. The so called trace theorem tells us that boundary values of \(H^1\) functions are well defined:

\[ u_{|\partial \Omega} \in L_2(\partial \Omega) \]

Thus it makes sense to define the sub-space satisfying homogeneous Dirichlet boundary conditions

\[ H_0^1(\Omega) = \{ u \in H^1(\Omega) : u_{|\partial \Omega} = 0 \}  \]

Let us consider the term on the left hand side of the variational formulation:

\[ A(u,v) := \int_{\Omega} \nabla u \nabla v \, dx \]

For given functions \(u\) and \(v\) from the Sobolev space, we compute the number \(\int \nabla u \nabla v \, dx\). Thus, \(A(.,.)\) is a function mapping from two elements from \(H^1\) into \({\mathbb R}\):

\[ A(.,.) : H^1 \times H^1 \rightarrow {\mathbb R} \]

The function \(A(.,.)\) is linear in both arguments, and thus we call it a bilinear-form.

Similarly, the right hand side

\[ f(v) := \int_{\Omega} f v \, dx \]

is a linear function

\[ f(.) : H^1 \rightarrow {\mathbb R}, \]

which we call a linear form. We use the same symbol for the function \(f : \Omega \rightarrow {\mathbb R}\) and the linear-form \(f : H^1(\Omega) \rightarrow {\mathbb R}\).

Having these objects defined, the weak formulation reads now

\[ \text{find} \, u \in H_0^1 : \quad A(u,v) = f(v) \quad \forall \, v \in H_0^1 \]

This abstract formalism in Hilbert spaces, bilinear and linear forms apply for a large class of (elliptic) partial differential equations.

1.3. The Finite Element Method#

The weak formulation is the starting point for the finite element method. We cannot compute the solution in an infinite dimensional Hilbert space. But, we can define a finite dimensional sub-space

\[ V_h \subset H^1_0 \]

and restrict the weak formulation to \(V_h\):

\[ \text{find} \, u_h \in V_h : \quad A(u_h,v_h) = f(v_h) \quad \forall \, v_h \in V_h \]

This is called the Galerkin method. The finite element solution \(u_h\) is some approximation to the true solution \(u\). We will analyze the discretization error \(\| u - u_h \|_{H^1}\).

For computing the discrete solution \(u_h\) we have to choose a basis for the function space \(V_h\), where \(N = \operatorname{dim} V_h\)

\[ V_h = \operatorname{span} \{ p_1(x), \ldots, p_N(x) \} \]

By means of this basis we can expand the solution \(u_h\) as

\[ u_h(x) = \sum_{i=1}^N u_i p_i(x) \]

The coefficients \(u_i\) are combined to the coefficient vector \(u = (u_1, \ldots, u_N) \in {\mathbb R}^N\).

Instead of testing with all test-functions from \(V_h\), by linearity of \(A(.,.)\) and \(f(.)\), it is enough to test only with the basis functions \(p_j(x), j = 1, \ldots, N\)

Thus, the finite element problem can be rewritten as

\[ \text{find } u \in {\mathbb R}^N : \quad A(\sum_i u_i p_i, p_j) = f(p_j) \qquad \forall \, j = 1, \ldots N \]

By linearity of \(A(.,.)\) in the first argument we can write

\[ \text{find } u \in {\mathbb R}^N : \quad \sum_{i=1}^N A(p_i, p_j) \, u_i = f(p_j) \qquad \forall \, j = 1, \ldots N \]

Since the basis functions are known, we can compute the matrix \(A \in {\mathbb R}^{N\times N}\) with entries

\[ A_{j,i} = A(p_i,p_j) = \int_\Omega \nabla p_i(x) \nabla p_j(x) \, dx \]

and the vector \(f \in {\mathbb R}^N\) as

\[ f_j = f(p_j) = \int_\Omega f(x) p_j(x) \, dx \]

Solving the finite element problem results in the linear system of equations for the coefficient vector \(u = (u_1, \ldots, u_N)\):

\[ \text{find } u \in {\mathbb R}^N : \quad A u = f \]

By means of the coefficient vector, we have a representation of the finite element solution

\[ u_h(x) = \sum_{i=1}^N u_i p_i(x) \]

1.4. Poisson equation in NGSolve:#

The Python interface to NGSolve allows us to enter the equation very close to its mathematical formulation.

# load Netgen/NGSolve 
from ngsolve import *
from ngsolve.webgui import Draw

The unit-square \(\Omega = (0,1)^2\) is a predefined domain in NGSolve:

Draw(unit_square.shape);

We generate a mesh (also called a triangulation) for \(\Omega\) of mesh-size \(h=0.2\):

mesh = Mesh(unit_square.GenerateMesh(maxh=0.2))
Draw (mesh);

Number of vertices and elements:

mesh.nv, mesh.ne
(39, 56)

define the finite element space (dof is short for degree of freedom):

fes = H1(mesh, order=3, dirichlet=".*")
print ("number of dofs =", fes.ndof)
number of dofs = 283

Define the bilinear-form. We define it by means of trial- and test-functions of the space:

u = fes.TrialFunction()
v = fes.TestFunction()
a = BilinearForm(grad(u)*grad(v)*dx)

Similarly, we define the linear-form. The function funcf is defined by means of x and y, which are pre-defined symbols for Cartesian coordinates:

funcf = 50*x*y
f = LinearForm(funcf*v*dx)

compute the matrix and right hand side vector:

a.Assemble()
f.Assemble();

solve the linear system. Restrict the set of basis functions to non-Dirichlet basis functions (freedofs):

gfu = GridFunction(fes)
gfu.vec.data = a.mat.Inverse(freedofs=fes.FreeDofs()) * f.vec
Draw (gfu);
Draw (grad(gfu), mesh, order=3, vectors = True);

We can inspect the matrix entries. Either directly printing the matrix, or convert the sparse matrix to scipy:

# (i,j,val) = a.mat.COO()
# print (list(i),list(j),list(val))
# print (a.mat)
from scipy.sparse import csr_matrix
print (csr_matrix(a.mat.CSR()))
  (0, 0)	0.9999999999999964
  (0, 4)	-0.499999999999999
  (0, 19)	-0.4999999999999974
  (0, 39)	-0.08333333333333291
  (0, 40)	-2.662800535624399e-16
  (0, 41)	-0.08333333333333307
  (0, 42)	-3.1181654480683108e-16
  (0, 61)	0.16666666666666605
  (0, 62)	-4.85722573273506e-17
  (0, 227)	-5.204170427930421e-17
  (1, 1)	0.8287315681024927
  (1, 7)	-0.1957737541722982
  (1, 8)	-0.19794827628847217
  (1, 23)	-0.4350095376417224
  (1, 43)	-0.036081925033399936
  (1, 44)	-1.2706849461530112e-16
  (1, 45)	-0.036419664573553744
  (1, 46)	-1.0793232627093197e-16
  (1, 47)	-0.06562033841012835
  (1, 48)	-2.478486166301863e-16
  (1, 79)	0.06871088406211634
  (1, 80)	1.734723475976807e-17
  (1, 83)	0.06941104395496575
  (1, 84)	-2.949029909160572e-17
  (1, 233)	-1.9081958235744878e-17
  :	:
  (280, 221)	-0.004044215842409425
  (280, 222)	-0.00010949506857410921
  (280, 225)	-0.005055833317030188
  (280, 226)	-0.0004467008934476962
  (280, 280)	0.009656245121461408
  (281, 34)	-2.2551405187698492e-17
  (281, 36)	5.724587470723463e-17
  (281, 38)	-3.8163916471489756e-17
  (281, 215)	-0.003295602986137172
  (281, 216)	-0.00037798057948934126
  (281, 217)	-0.006250616150810586
  (281, 218)	0.0006070238087351294
  (281, 223)	-0.0051166744123425645
  (281, 224)	0.0009850043882244734
  (281, 281)	0.009775262366193539
  (282, 28)	-2.0816681711721685e-17
  (282, 35)	7.632783294297951e-17
  (282, 38)	-5.204170427930421e-17
  (282, 187)	-0.004200374856019856
  (282, 188)	0.0017754830163756964
  (282, 191)	-0.0028410571431563554
  (282, 192)	0.0013223771120878588
  (282, 221)	-0.008167506192283442
  (282, 222)	-0.00045310590428783134
  (282, 282)	0.010139292127639758

the right hand side vector:

import numpy as np
print (np.array(f.vec))
[ 6.66666667e-04  2.62349341e-02  3.00666667e-01  2.41938614e-02
  1.58768986e-02  3.27043380e-02  3.18198359e-02  2.75906396e-02
  1.61373488e-01  3.37816023e-01  4.59571120e-01  7.65080144e-01
  9.90491892e-01  5.15375584e-01  3.72201759e-01  1.38794309e-01
  2.46500883e-02  3.01059709e-02  3.22396274e-02  1.56496750e-02
  1.66524538e-01  1.92210142e-01  1.90956415e-01  1.39617577e-01
  5.06041231e-01  6.86702751e-01  1.06967160e+00  7.97947213e-01
  8.91277642e-01  4.11175628e-01  1.15852596e-01  1.72752866e-01
  1.95764202e-01  3.52020762e-01  3.33003422e-01  6.68347587e-01
  4.23554864e-01  3.78842972e-01  5.04632471e-01 -1.11111111e-04
 -1.58730159e-05 -1.11111111e-04 -1.58730159e-05 -7.12517048e-04
  7.65287341e-06 -3.04844263e-03 -3.80767289e-04 -4.05785815e-03
 -4.57756990e-04 -3.67777778e-02  5.39682540e-04 -3.67777778e-02
  5.39682540e-04 -2.86914888e-03 -3.63421433e-04 -6.35674978e-04
  6.81090441e-06 -3.70825255e-03 -4.12198974e-04 -6.99252178e-04
 -2.37848562e-05 -1.15456810e-03  5.48026033e-07 -2.69786564e-03
 -4.22507913e-04 -7.64996884e-04 -1.46795278e-05 -4.01918863e-03
 -4.25169512e-04 -4.24131829e-03 -5.88970269e-04 -8.02538052e-04
 -1.09561861e-05 -3.95165835e-03 -4.89298369e-04 -4.00668231e-03
 -5.56080516e-04 -3.50759366e-03 -4.62624400e-04 -3.23170524e-03
 -4.38776040e-04 -1.14247326e-02 -5.05147539e-04 -9.20161982e-03
  3.17104988e-04 -1.68410527e-02 -2.96251115e-04 -1.63058946e-02
 -4.33316884e-04 -2.54836020e-02  6.99673786e-04 -3.04333400e-02
 -2.28614435e-04 -2.28786796e-02 -4.36373064e-04 -3.42984684e-02
  6.67975607e-04 -4.05017114e-02 -5.21422613e-05 -7.22749714e-02
 -1.01231228e-05 -5.80040475e-02  1.03235807e-03 -2.36233754e-02
  4.56415614e-04 -6.19613792e-02  1.11052790e-03 -4.97669399e-02
  1.13946661e-03 -2.16227787e-02  5.83208976e-04 -4.16902662e-02
  1.47682560e-04 -3.98236278e-02  1.10835203e-03 -1.00587075e-02
  4.54285092e-04 -3.56260820e-02  2.85855089e-05 -2.46723283e-02
  8.19414012e-04 -1.39435204e-02 -1.63404990e-04 -7.76492876e-03
  2.99724319e-04 -7.30871063e-04  9.91997734e-06 -2.86300492e-03
 -3.89843315e-04 -3.14967093e-03 -4.20416263e-04 -7.35262210e-04
  1.40034946e-05 -3.75317153e-03 -5.20003175e-04 -3.80161548e-03
 -4.77881017e-04 -6.75704791e-04  2.28519735e-05 -3.98506472e-03
 -4.24637796e-04 -4.20246277e-03 -5.86533680e-04 -2.65603819e-03
 -4.17798587e-04 -8.39661832e-03 -2.42512937e-04 -8.38610419e-03
 -2.41806640e-04 -1.08513052e-02 -7.15531330e-04 -7.84374591e-03
 -1.03591803e-04 -1.14717576e-02 -5.97567491e-04 -1.18746770e-02
 -4.11825842e-04 -7.02642613e-03 -1.82433803e-05 -1.28400582e-02
 -6.89738810e-04 -1.26215991e-02 -4.78832306e-04 -1.13978020e-02
 -5.94066240e-04 -2.46507838e-02 -8.45688875e-04 -1.63499185e-02
  2.47657413e-04 -1.89614036e-02 -2.81073187e-04 -3.41063837e-02
 -7.04901124e-04 -2.63998469e-02  1.67903582e-05 -2.15873168e-02
  4.26095800e-04 -4.15866512e-02  2.03394806e-04 -3.20360315e-02
  6.84372631e-04 -3.44292495e-02  8.22649519e-04 -3.22214879e-02
  5.30138177e-04 -2.69926933e-02  8.68249759e-04 -3.13317477e-02
 -2.38371554e-04 -2.56081511e-02  1.01554149e-03 -2.86452589e-02
  5.65359773e-04 -8.76577828e-03  4.32474227e-04 -1.06852126e-02
  5.29948886e-04 -1.83915107e-02  1.29215552e-04 -5.75494085e-03
  6.93996993e-06 -7.81598785e-03  8.70080819e-05 -1.20710259e-02
 -5.14225686e-04 -1.23298683e-02 -4.26887948e-04 -1.22447412e-02
 -6.65982942e-04 -1.46644119e-02  2.43566373e-04 -1.62363235e-02
 -4.92491430e-04 -1.63463979e-02 -3.55921383e-04 -1.58003557e-02
 -2.97584791e-04 -1.72386981e-02 -6.69680035e-04 -2.04304343e-02
  3.88629681e-04 -2.45816144e-02  6.74252409e-04 -2.19168273e-02
 -4.23588458e-04 -1.73491052e-02  1.46177993e-04  6.34920635e-05
  4.56649872e-04  5.08916336e-04  1.18856044e-03  5.34503153e-04
  1.07109762e-03  4.72977177e-04  8.24123503e-04  1.25322175e-03
  4.50016337e-03  2.37266314e-03  6.43407504e-03  5.68927196e-03
  8.99761394e-03  7.24523562e-03  1.45079365e-02  9.20278558e-03
  8.39634548e-03  7.47646191e-03  3.91395496e-03  5.94214740e-03
  1.17034605e-03  1.83371241e-03  4.22267988e-04  4.87581228e-04
  7.27688891e-04  4.90189805e-04  1.02504840e-03  4.42063694e-04
  5.08754364e-04  1.43069641e-02  1.04322115e-02  1.19502289e-03
  2.05886485e-03  2.17266415e-03  2.52444753e-03  6.37948302e-03
  6.38792594e-03  2.06372227e-03  3.03307384e-03  1.60426442e-03
  2.68795971e-03  4.79711823e-03  6.46811512e-03  2.12007716e-03
  2.15479975e-03  2.75456096e-03  4.23770460e-03  3.41373102e-03
  5.32578507e-03  4.28331106e-03  3.12604564e-03  3.34032616e-03
  3.83432663e-03  3.54285512e-03  6.03732654e-03]

and the solution vector:

print (np.array(gfu.vec))
[ 0.00000000e+00  0.00000000e+00  0.00000000e+00  0.00000000e+00
  0.00000000e+00  0.00000000e+00  0.00000000e+00  0.00000000e+00
  0.00000000e+00  0.00000000e+00  0.00000000e+00  0.00000000e+00
  0.00000000e+00  0.00000000e+00  0.00000000e+00  0.00000000e+00
  0.00000000e+00  0.00000000e+00  0.00000000e+00  0.00000000e+00
  3.02651352e-01  4.17444909e-01  3.71729962e-01  2.04295345e-01
  5.53986083e-01  7.31340616e-01  7.73825402e-01  8.10673815e-01
  8.42524628e-01  4.68063272e-01  1.86388258e-01  3.49944425e-01
  4.12612009e-01  7.63003331e-01  7.05615795e-01  1.04452565e+00
  7.52636441e-01  9.34678613e-01  9.89120579e-01  0.00000000e+00
  0.00000000e+00  0.00000000e+00  0.00000000e+00  0.00000000e+00
  0.00000000e+00  0.00000000e+00  0.00000000e+00  3.45788548e-01
 -6.14103409e-02  0.00000000e+00  0.00000000e+00  0.00000000e+00
  0.00000000e+00  0.00000000e+00  0.00000000e+00  0.00000000e+00
  0.00000000e+00  3.18773436e-01 -5.28179168e-02  0.00000000e+00
  0.00000000e+00 -4.42879081e-01  1.06815699e-05  5.31053506e-02
 -1.77984470e-02  0.00000000e+00  0.00000000e+00 -3.97926710e-01
 -3.30288929e-02  5.50442051e-02 -2.21448524e-02  0.00000000e+00
  0.00000000e+00 -1.05261663e-01 -3.13658800e-02 -6.36909436e-02
 -2.30269301e-02  8.45569001e-02 -3.52031126e-02 -1.22889271e-01
 -1.06173136e-02  0.00000000e+00  0.00000000e+00 -3.22902254e-01
  3.12502359e-02  3.45261752e-03 -8.05072273e-02  0.00000000e+00
  0.00000000e+00 -7.47899861e-01  7.95574753e-02 -2.65063196e-01
 -6.78173444e-02  0.00000000e+00  0.00000000e+00 -6.91358750e-01
  5.93682062e-02 -7.50358631e-01 -5.90154230e-02 -2.92176436e+00
  3.42272101e-04 -6.33243573e-01  5.39347185e-02  0.00000000e+00
  0.00000000e+00 -1.82286681e+00  3.99272270e-01 -6.88735676e-01
  3.33327795e-02  0.00000000e+00  0.00000000e+00 -7.23052075e-01
 -1.59289558e-03 -1.17484926e+00  1.38832225e-01  0.00000000e+00
  0.00000000e+00 -6.11274286e-01 -4.02914557e-02 -7.24856511e-01
  7.46255587e-02 -2.46730851e-02 -4.39444956e-02 -3.23493810e-01
  3.09680027e-02  0.00000000e+00  0.00000000e+00 -1.30408705e-01
 -1.18844850e-02  7.93070738e-02 -2.73128488e-02  0.00000000e+00
  0.00000000e+00 -7.24729233e-02 -2.64119171e-02 -8.88135271e-02
 -2.67198008e-02  0.00000000e+00  0.00000000e+00 -3.76652503e-01
 -3.14929094e-02  6.10455003e-02 -2.42626166e-02  6.40022829e-02
 -1.75053646e-02 -2.32202735e-01 -5.64428432e-03 -2.40405870e-01
 -5.83425690e-03  3.70940954e-02 -3.89899776e-02 -1.58213525e-01
 -3.13724910e-03 -8.02840051e-02 -2.56969874e-02 -2.64893941e-01
 -2.17444576e-02 -9.90065562e-02  3.76870051e-03 -2.78117344e-01
 -3.60372397e-02 -1.11519878e-01 -3.09574110e-02  4.41019727e-02
 -1.91093951e-02 -1.83757108e-01 -1.19954481e-02 -3.89402445e-01
  2.30823297e-02 -1.72535447e-01 -4.07867885e-02 -2.55036412e-01
 -3.63605128e-02 -4.76999840e-01 -9.47830128e-03 -3.41001347e-01
  2.89867542e-02 -6.32947036e-01  4.43052514e-02 -3.84225488e-01
  2.27470834e-02 -2.50682327e-01  1.04445533e-02 -4.82746652e-01
  5.23222011e-02 -9.08658007e-02  3.93362996e-02 -5.27068308e-01
 -4.82734566e-02 -4.04378948e-01  4.76404167e-02 -6.14677104e-01
  6.25805960e-02  3.71933788e-02  9.70398951e-03 -2.72312767e-01
  3.55951239e-02 -3.23468220e-01  8.17363736e-04 -8.43405785e-02
 -4.43606881e-03 -1.57030863e-01  1.14401498e-03 -8.27373870e-02
 -2.67211217e-02 -2.46279346e-01 -2.12424681e-02 -1.08686469e-01
 -3.31923341e-02 -2.75245664e-01  1.18378034e-02 -1.54910123e-01
 -1.48356595e-02 -2.63609902e-01 -2.71076796e-02 -3.21921127e-01
 -1.70420964e-02 -1.49706348e-01 -2.91715363e-02 -2.16192264e-01
  2.31462394e-02 -2.66339889e-01  2.07213824e-02 -2.66713703e-01
 -2.90411155e-02 -2.71013089e-01  1.06479309e-02  5.71428571e-03
  1.11326140e-02  3.52915266e-02 -3.09862199e-02  6.21240144e-02
 -4.81353155e-02  7.81662709e-02 -4.85774805e-02 -5.99280637e-03
  4.06924526e-02 -5.34826800e-02  7.93598132e-02 -3.53171058e-02
  2.10459024e-01 -1.29849936e-01  1.30571429e+00  2.25121637e-01
  1.58097206e-01 -1.09384602e-01  5.10999455e-02 -1.34576114e-01
 -2.34808772e-03 -3.81873156e-02  7.58369941e-02  6.07623135e-02
 -4.48304492e-02  3.58390741e-02 -4.99116470e-02  1.17018031e-02
 -1.05473882e-02 -5.94354799e-01 -9.68751578e-02 -3.52761847e-02
  4.34701629e-02  2.15674039e-02 -1.91494610e-02  9.67167662e-02
 -7.43369229e-02  5.55632805e-02 -6.84206913e-02  4.43752689e-02
 -6.62395745e-02  8.34223115e-02  6.40049428e-02  4.88721028e-02
  2.65192015e-02 -3.09425754e-02  2.30844242e-02  4.25440808e-02
 -5.64691905e-02 -5.23352169e-02 -1.19756385e-02  1.34050016e-02
  2.30074832e-02  2.26447785e-02 -8.47010385e-03]

1.5. Visualizing the basis functions#

For a finite element space of order=1, the basis-functions are associated with mesh vertices. For higher order, we also have edge and face based basis-functions.

Let’s have a look at the basis functions. Well know are the so-called hat-functions, which have value 1 in one vertex, 0 in all other, and a linear interpolation in between.

Open controls, and use the multidim slider to select basis function

gf = GridFunction(fes, multidim=mesh.nv)
for i in range (mesh.nv):
    gf.vecs[i][:] = 0
    gf.vecs[i][i] = 1
scene = Draw (gf, mesh, deformation=True, animate=True);

For higher order spaces there are also basis functions on edges:

edgedofs = fes.GetDofNrs(NodeId(EDGE,40))
print ("dofs on edge #40: ", edgedofs)

gf = GridFunction(fes)
for edof in edgedofs:
    gf.vec[:] = 0
    gf.vec[edof] = 1
    Draw (-3*gf, mesh, deformation=True, order=3)
dofs on edge #40:  (119, 120)

and also basis functions associated with faces:

facedofs = fes.GetDofNrs(NodeId(FACE,0))
print ("dofs on face #0: ", facedofs)

for fdof in facedofs:
    gf.vec[:] = 0
    gf.vec[fdof] = 1  
    Draw (5*gf, mesh, deformation=True, order=3)
dofs on face #0:  (227,)

To get help on types or objects use the help function:

help (Mesh)
Help on class Mesh in module ngsolve.comp:

class Mesh(pybind11_builtins.pybind11_object)
 |  NGSolve interface to the Netgen mesh. Provides access and functionality
 |  to use the mesh for finite element calculations.
 |
 |  Parameters:
 |
 |  mesh (netgen.Mesh): a mesh generated from Netgen
 |
 |  Method resolution order:
 |      Mesh
 |      pybind11_builtins.pybind11_object
 |      builtins.object
 |
 |  Methods defined here:
 |
 |  BBBoundaries(...)
 |      BBBoundaries(self: ngsolve.comp.Mesh, pattern: str) -> ngsolve.comp.Region
 |
 |      Return co dim 3 boundary mesh-region matching the given regex pattern
 |
 |  BBoundaries(...)
 |      BBoundaries(self: ngsolve.comp.Mesh, pattern: str) -> ngsolve.comp.Region
 |
 |      Return co dim 2 boundary mesh-region matching the given regex pattern
 |
 |  Boundaries(...)
 |      Boundaries(*args, **kwargs)
 |      Overloaded function.
 |
 |      1. Boundaries(self: ngsolve.comp.Mesh, pattern: str) -> ngsolve.comp.Region
 |
 |      Return boundary mesh-region matching the given regex pattern
 |
 |      2. Boundaries(self: ngsolve.comp.Mesh, bnds: List[int]) -> ngsolve.comp.Region
 |
 |      Generate boundary mesh-region by boundary condition numbers
 |
 |  BoundaryCF(...)
 |      BoundaryCF(self: ngsolve.comp.Mesh, values: dict, default: ngsolve.fem.CoefficientFunction = None) -> ngsolve.fem.CoefficientFunction
 |
 |      Boundary wise CoefficientFunction.
 |      First argument is a dict from either boundary names or Region objects to
 |      CoefficientFunction (-values). Later given names/regions override earlier
 |      values. Optional last argument (default) is the value for not given boundaries.
 |      >>> penalty = mesh.BoundaryCF({ "top" : 1e6 }, default=0)
 |      will create a CF being 1e6 on the top boundary and 0. elsewhere.
 |
 |  BuildRefinementTree(...)
 |      BuildRefinementTree(self: ngsolve.comp.Mesh) -> pyngcore.pyngcore.Array_y_S
 |
 |  Contains(...)
 |      Contains(self: ngsolve.comp.Mesh, x: float = 0.0, y: float = 0.0, z: float = 0.0) -> bool
 |
 |      Check if the point (x,y,z) is in the meshed domain (is inside a volume element)
 |
 |  Curve(...)
 |      Curve(self: ngsolve.comp.Mesh, order: int) -> ngsolve.comp.Mesh
 |
 |      Curve the mesh elements for geometry approximation of given order
 |
 |  Elements(...)
 |      Elements(self: ngsolve.comp.Mesh, VOL_or_BND: ngsolve.comp.VorB = <VorB.VOL: 0>) -> ngcomp::ElementRange
 |
 |      Return an iterator over elements on VOL/BND
 |
 |  GeoParamCF(...)
 |      GeoParamCF(self: ngsolve.comp.Mesh) -> ngsolve.fem.CoefficientFunction
 |
 |  GetBBBoundaries(...)
 |      GetBBBoundaries(self: ngsolve.comp.Mesh) -> tuple
 |
 |      Return list of boundary conditions for co dimension 3
 |
 |  GetBBoundaries(...)
 |      GetBBoundaries(self: ngsolve.comp.Mesh) -> tuple
 |
 |      Return list of boundary conditions for co dimension 2
 |
 |  GetBoundaries(...)
 |      GetBoundaries(self: ngsolve.comp.Mesh) -> tuple
 |
 |      Return list of boundary condition names
 |
 |  GetCurveOrder(...)
 |      GetCurveOrder(self: ngsolve.comp.Mesh) -> int
 |
 |  GetHPElementLevel(...)
 |      GetHPElementLevel(self: ngsolve.comp.Mesh, ei: ngsolve.comp.ElementId) -> int
 |
 |      THIS FUNCTION IS WIP!
 |       Return HP-refinement level of element
 |
 |  GetMaterials(...)
 |      GetMaterials(self: ngsolve.comp.Mesh) -> tuple
 |
 |      Return list of material names
 |
 |  GetNE(...)
 |      GetNE(self: ngsolve.comp.Mesh, arg0: ngsolve.comp.VorB) -> int
 |
 |      Return number of elements of codimension VorB.
 |
 |  GetPMLTrafo(...)
 |      GetPMLTrafo(self: ngsolve.comp.Mesh, dom: int = 1) -> ngsolve.comp.pml.PML
 |
 |      Return pml transformation on domain dom
 |
 |  GetPMLTrafos(...)
 |      GetPMLTrafos(self: ngsolve.comp.Mesh) -> list
 |
 |      Return list of pml transformations
 |
 |  GetParentElement(...)
 |      GetParentElement(self: ngsolve.comp.Mesh, ei: ngsolve.comp.ElementId) -> ngsolve.comp.ElementId
 |
 |      Return parent element id on refined mesh
 |
 |  GetParentFaces(...)
 |      GetParentFaces(self: ngsolve.comp.Mesh, fnum: int) -> tuple
 |
 |      Return parent faces
 |
 |  GetParentVertices(...)
 |      GetParentVertices(self: ngsolve.comp.Mesh, vnum: int) -> tuple
 |
 |      Return parent vertex numbers on refined mesh
 |
 |  GetPeriodicNodePairs(...)
 |      GetPeriodicNodePairs(self: ngsolve.comp.Mesh, arg0: ngsolve.fem.NODE_TYPE) -> list
 |
 |      returns list of periodic nodes with their identification number as [((master_nr, minion_nr),idnr),...]
 |
 |  GetTrafo(...)
 |      GetTrafo(self: ngsolve.comp.Mesh, eid: ngsolve.comp.ElementId) -> ngsolve.fem.ElementTransformation
 |
 |      returns element transformation of given element id
 |
 |  LocalHCF(...)
 |      LocalHCF(self: ngsolve.comp.Mesh) -> ngsolve.fem.CoefficientFunction
 |
 |  MapToAllElements(...)
 |      MapToAllElements(*args, **kwargs)
 |      Overloaded function.
 |
 |      1. MapToAllElements(self: ngsolve.comp.Mesh, arg0: ngsolve.fem.IntegrationRule, arg1: Union[ngsolve.comp.VorB, ngsolve.comp.Region]) -> numpy.ndarray[ngsolve.fem.MeshPoint]
 |
 |      2. MapToAllElements(self: ngsolve.comp.Mesh, arg0: Dict[ngsolve.fem.ET, ngsolve.fem.IntegrationRule], arg1: Union[ngsolve.comp.VorB, ngsolve.comp.Region]) -> numpy.ndarray[ngsolve.fem.MeshPoint]
 |
 |  MaterialCF(...)
 |      MaterialCF(self: ngsolve.comp.Mesh, values: dict, default: ngsolve.fem.CoefficientFunction = None) -> ngsolve.fem.CoefficientFunction
 |
 |      Domain wise CoefficientFunction.
 |      First argument is a dict from either material names or Region objects to
 |      CoefficientFunction (-values). Later given names/regions override earlier
 |      values. Optional last argument (default) is the value for not given materials.
 |      >>> sigma = mesh.MaterialCF({ "steel_.*" : 2e6 }, default=0)
 |      will create a CF being 2e6 on all domains starting with 'steel_' and 0 elsewhere.
 |
 |  Materials(...)
 |      Materials(*args, **kwargs)
 |      Overloaded function.
 |
 |      1. Materials(self: ngsolve.comp.Mesh, pattern: str) -> ngsolve.comp.Region
 |
 |      Return mesh-region matching the given regex pattern
 |
 |      2. Materials(self: ngsolve.comp.Mesh, domains: List[int]) -> ngsolve.comp.Region
 |
 |      Generate mesh-region by domain numbers
 |
 |  Refine(...)
 |      Refine(self: ngsolve.comp.Mesh, mark_surface_elements: bool = False, onlyonce: bool = False) -> None
 |
 |      Local mesh refinement based on marked elements, uses element-bisection algorithm
 |
 |  RefineFromTree(...)
 |      RefineFromTree(self: ngsolve.comp.Mesh, arg0: pyngcore.pyngcore.Array_y_S) -> None
 |
 |  RefineHP(...)
 |      RefineHP(self: ngsolve.comp.Mesh, levels: int, factor: float = 0.125) -> None
 |
 |      Geometric mesh refinement towards marked vertices and edges, uses factor for placement of new points
 |
 |  Region(...)
 |      Region(self: ngsolve.comp.Mesh, vb: ngsolve.comp.VorB, pattern: Optional[str] = '.*') -> ngsolve.comp.Region
 |
 |      Return boundary mesh-region matching the given regex pattern
 |
 |  RegionCF(...)
 |      RegionCF(self: ngsolve.comp.Mesh, VorB: ngsolve.comp.VorB, value: dict, default: ngsolve.fem.CoefficientFunction = None) -> ngsolve.fem.CoefficientFunction
 |
 |      Region wise CoefficientFunction.
 |      First argument is VorB, defining the co-dimension,
 |      second argument is a dict from either region names or Region objects to
 |      CoefficientFunction (-values). Later given names/regions override earlier
 |      values. Optional last argument (default) is the value for not given regions.
 |      >>> sigma = mesh.RegionCF(VOL, { "steel_.*" : 2e6 }, default=0)
 |      will create a CF being 2e6 on all domains starting with 'steel_' and 0 elsewhere.
 |
 |  SetDeformation(...)
 |      SetDeformation(self: ngsolve.comp.Mesh, gf: ngcomp::GridFunction) -> None
 |
 |      Deform the mesh with the given GridFunction
 |
 |  SetElementOrder(...)
 |      SetElementOrder(self: ngsolve.comp.Mesh, eid: ngsolve.comp.ElementId, order: int) -> None
 |
 |      For backward compatibility, not recommended to use
 |
 |  SetPML(...)
 |      SetPML(self: ngsolve.comp.Mesh, pmltrafo: ngsolve.comp.pml.PML, definedon: object) -> None
 |
 |      Set PML transformation on domain
 |
 |  SetRefinementFlag(...)
 |      SetRefinementFlag(self: ngsolve.comp.Mesh, ei: ngsolve.comp.ElementId, refine: bool) -> None
 |
 |      Set refinementflag for mesh-refinement
 |
 |  SetRefinementFlags(...)
 |      SetRefinementFlags(self: ngsolve.comp.Mesh, refine: List[bool]) -> None
 |
 |      Set refinementflags for mesh-refinement
 |
 |  SplitElements_Alfeld(...)
 |      SplitElements_Alfeld(self: ngsolve.comp.Mesh) -> None
 |
 |  UnSetPML(...)
 |      UnSetPML(self: ngsolve.comp.Mesh, definedon: object) -> None
 |
 |      Unset PML transformation on domain
 |
 |  UnsetDeformation(...)
 |      UnsetDeformation(self: ngsolve.comp.Mesh) -> None
 |
 |      Unset the deformation
 |
 |  __call__(...)
 |      __call__(self: ngsolve.comp.Mesh, x: numpy.ndarray[numpy.float64] = 0.0, y: numpy.ndarray[numpy.float64] = 0.0, z: numpy.ndarray[numpy.float64] = 0.0, VOL_or_BND: ngsolve.comp.VorB = <VorB.VOL: 0>) -> object
 |
 |      Get a MappedIntegrationPoint in the point (x,y,z) on the matching volume (VorB=VOL, default) or surface (VorB=BND) element. BBND elements aren't supported
 |
 |  __eq__(...)
 |      __eq__(self: ngsolve.comp.Mesh, mesh: ngsolve.comp.Mesh) -> bool
 |
 |  __getitem__(...)
 |      __getitem__(*args, **kwargs)
 |      Overloaded function.
 |
 |      1. __getitem__(self: ngsolve.comp.Mesh, arg0: ngsolve.comp.ElementId) -> ngsolve.comp.Ngs_Element
 |
 |      Return Ngs_Element from given ElementId
 |
 |      2. __getitem__(self: ngsolve.comp.Mesh, arg0: ngsolve.comp.NodeId) -> ngsolve.comp.MeshNode
 |
 |      Return MeshNode from given NodeId
 |
 |  __getstate__(...)
 |      __getstate__(self: ngsolve.comp.Mesh) -> tuple
 |
 |  __init__(...)
 |      __init__(*args, **kwargs)
 |      Overloaded function.
 |
 |      1. __init__(self: ngsolve.comp.Mesh, ngmesh: netgen.libngpy._meshing.Mesh) -> None
 |
 |      Make an NGSolve-mesh from a Netgen-mesh
 |
 |      2. __init__(self: ngsolve.comp.Mesh, filename: str, comm: netgen.libngpy._meshing.MPI_Comm = <netgen.libngpy._meshing.MPI_Comm object at 0x10627a8b0>) -> None
 |
 |      Load a mesh from file.
 |      In MPI-parallel mode the mesh is distributed over the MPI-group given by the communicator (WIP!)
 |
 |  __setstate__(...)
 |      __setstate__(self: ngsolve.comp.Mesh, arg0: tuple) -> None
 |
 |  nnodes(...)
 |      nnodes(self: ngsolve.comp.Mesh, arg0: ngsolve.fem.NODE_TYPE) -> int
 |
 |      number of nodes given type
 |
 |  nodes(...)
 |      nodes(self: ngsolve.comp.Mesh, node_type: ngsolve.fem.NODE_TYPE) -> ngsolve.comp.MeshNodeRange
 |
 |      iterable of mesh nodes of type node_type
 |
 |  ----------------------------------------------------------------------
 |  Readonly properties defined here:
 |
 |  comm
 |      MPI-communicator the Mesh lives in
 |
 |  dim
 |      mesh dimension
 |
 |  edges
 |      iterable of mesh edges
 |
 |  faces
 |      iterable of mesh faces
 |
 |  facets
 |      iterable of mesh facets
 |
 |  levels
 |      multigrid levels
 |
 |  ne
 |      number of volume elements
 |
 |  nedge
 |      number of edges
 |
 |  nface
 |      number of faces
 |
 |  nfacet
 |      number of facets
 |
 |  ngmesh
 |      the Netgen mesh
 |
 |  nv
 |      number of vertices
 |
 |  vertices
 |      iterable of mesh vertices
 |
 |  ----------------------------------------------------------------------
 |  Data descriptors defined here:
 |
 |  __dict__
 |
 |  deformation
 |      mesh deformation
 |
 |  ----------------------------------------------------------------------
 |  Data and other attributes defined here:
 |
 |  __hash__ = None
 |
 |  ----------------------------------------------------------------------
 |  Static methods inherited from pybind11_builtins.pybind11_object:
 |
 |  __new__(*args, **kwargs) class method of pybind11_builtins.pybind11_object
 |      Create and return a new object.  See help(type) for accurate signature.
help (fes)
Help on H1 in module ngsolve.comp object:

class H1(FESpace)
 |  An H1-conforming finite element space.
 |
 |  The H1 finite element space consists of continuous and
 |  element-wise polynomial functions. It uses a hierarchical (=modal)
 |  basis built from integrated Legendre polynomials on tensor-product elements,
 |  and Jaboci polynomials on simplicial elements.
 |
 |  Boundary values are well defined. The function can be used directly on the
 |  boundary, using the trace operator is optional.
 |
 |  The H1 space supports variable order, which can be set individually for edges,
 |  faces and cells.
 |
 |  Internal degrees of freedom are declared as local dofs and are eliminated
 |  if static condensation is on.
 |
 |  The wirebasket consists of all vertex dofs. Optionally, one can include the
 |  first (the quadratic bubble) edge basis function, or all edge basis functions
 |  into the wirebasket.
 |
 |  Keyword arguments can be:
 |
 |  order: int = 1
 |    order of finite element space
 |  complex: bool = False
 |    Set if FESpace should be complex
 |  dirichlet: regexpr
 |    Regular expression string defining the dirichlet boundary.
 |    More than one boundary can be combined by the | operator,
 |    i.e.: dirichlet = 'top|right'
 |  dirichlet_bbnd: regexpr
 |    Regular expression string defining the dirichlet bboundary,
 |    i.e. points in 2D and edges in 3D.
 |    More than one boundary can be combined by the | operator,
 |    i.e.: dirichlet_bbnd = 'top|right'
 |  dirichlet_bbbnd: regexpr
 |    Regular expression string defining the dirichlet bbboundary,
 |    i.e. points in 3D.
 |    More than one boundary can be combined by the | operator,
 |    i.e.: dirichlet_bbbnd = 'top|right'
 |  definedon: Region or regexpr
 |    FESpace is only defined on specific Region, created with mesh.Materials('regexpr')
 |    or mesh.Boundaries('regexpr'). If given a regexpr, the region is assumed to be
 |    mesh.Materials('regexpr').
 |  dim: int = 1
 |    Create multi dimensional FESpace (i.e. [H1]^3)
 |  dgjumps: bool = False
 |    Enable discontinuous space for DG methods, this flag is needed for DG methods,
 |    since the dofs have a different coupling then and this changes the sparsity
 |    pattern of matrices.
 |  autoupdate: bool = False
 |    Automatically update on a change to the mesh.
 |  low_order_space: bool = True
 |    Generate a lowest order space together with the high-order space,
 |    needed for some preconditioners.
 |  order_policy: ORDER_POLICY = ORDER_POLICY.OLDSTYLE
 |    CONSTANT .. use the same fixed order for all elements,
 |    NODAL ..... use the same order for nodes of same shape,
 |    VARIABLE ... use an individual order for each edge, face and cell,
 |    OLDSTYLE .. as it used to be for the last decade
 |  wb_withedges: bool = true(3D) / false(2D)
 |    use lowest-order edge dofs for BDDC wirebasket
 |  wb_fulledges: bool = false
 |    use all edge dofs for BDDC wirebasket
 |
 |  Method resolution order:
 |      H1
 |      FESpace
 |      NGS_Object
 |      pybind11_builtins.pybind11_object
 |      builtins.object
 |
 |  Methods defined here:
 |
 |  __getstate__(...)
 |      __getstate__(self: ngsolve.comp.FESpace) -> tuple
 |
 |  __init__(...)
 |      __init__(self: ngsolve.comp.H1, mesh: ngsolve.comp.Mesh, **kwargs) -> None
 |
 |  __setstate__(...)
 |      __setstate__(self: ngsolve.comp.H1, arg0: tuple) -> None
 |
 |  ----------------------------------------------------------------------
 |  Static methods defined here:
 |
 |  __flags_doc__(...) method of builtins.PyCapsule instance
 |      __flags_doc__() -> dict
 |
 |  ----------------------------------------------------------------------
 |  Data descriptors defined here:
 |
 |  __dict__
 |
 |  ----------------------------------------------------------------------
 |  Methods inherited from FESpace:
 |
 |  ApplyM(...)
 |      ApplyM(self: ngsolve.comp.FESpace, vec: ngsolve.la.BaseVector, rho: ngsolve.fem.CoefficientFunction = None, definedon: ngsolve.comp.Region = None) -> None
 |
 |      Apply mass-matrix. Available only for L2-like spaces
 |
 |  ConvertL2Operator(...)
 |      ConvertL2Operator(self: ngsolve.comp.FESpace, l2space: ngsolve.comp.FESpace) -> BaseMatrix
 |
 |  CouplingType(...)
 |      CouplingType(self: ngsolve.comp.FESpace, dofnr: int) -> ngsolve.comp.COUPLING_TYPE
 |
 |
 |               Get coupling type of a degree of freedom.
 |
 |      Parameters:
 |
 |      dofnr : int
 |        input dof number
 |
 |  CreateDirectSolverCluster(...)
 |      CreateDirectSolverCluster(self: ngsolve.comp.FESpace, **kwargs) -> list
 |
 |  CreateSmoothingBlocks(...)
 |      CreateSmoothingBlocks(self: ngsolve.comp.FESpace, **kwargs) -> pyngcore.pyngcore.Table_I
 |
 |  Elements(...)
 |      Elements(*args, **kwargs)
 |      Overloaded function.
 |
 |      1. Elements(self: ngsolve.comp.FESpace, VOL_or_BND: ngsolve.comp.VorB = <VorB.VOL: 0>) -> ngsolve.comp.FESpaceElementRange
 |
 |
 |      Returns an iterable range of elements.
 |
 |      Parameters:
 |
 |      VOL_or_BND : ngsolve.comp.VorB
 |        input VOL, BND, BBND,...
 |
 |
 |
 |      2. Elements(self: ngsolve.comp.FESpace, arg0: ngsolve.comp.Region) -> Iterator
 |
 |  FinalizeUpdate(...)
 |      FinalizeUpdate(self: ngsolve.comp.FESpace) -> None
 |
 |      finalize update
 |
 |  FreeDofs(...)
 |      FreeDofs(self: ngsolve.comp.FESpace, coupling: bool = False) -> pyngcore.pyngcore.BitArray
 |
 |
 |
 |      Return BitArray of free (non-Dirichlet) dofs\n
 |      coupling=False ... all free dofs including local dofs\n
 |      coupling=True .... only element-boundary free dofs
 |
 |      Parameters:
 |
 |      coupling : bool
 |        input coupling
 |
 |  GetDofNrs(...)
 |      GetDofNrs(*args, **kwargs)
 |      Overloaded function.
 |
 |      1. GetDofNrs(self: ngsolve.comp.FESpace, ei: ngsolve.comp.ElementId) -> tuple
 |
 |
 |
 |      Parameters:
 |
 |      ei : ngsolve.comp.ElementId
 |        input element id
 |
 |
 |
 |      2. GetDofNrs(self: ngsolve.comp.FESpace, ni: ngsolve.comp.NodeId) -> tuple
 |
 |
 |
 |      Parameters:
 |
 |      ni : ngsolve.comp.NodeId
 |        input node id
 |
 |  GetDofs(...)
 |      GetDofs(self: ngsolve.comp.FESpace, region: ngsolve.comp.Region) -> pyngcore.pyngcore.BitArray
 |
 |
 |      Returns all degrees of freedom in given region.
 |
 |      Parameters:
 |
 |      region : ngsolve.comp.Region
 |        input region
 |
 |  GetFE(...)
 |      GetFE(self: ngsolve.comp.FESpace, ei: ngsolve.comp.ElementId) -> object
 |
 |
 |      Get the finite element to corresponding element id.
 |
 |      Parameters:
 |
 |      ei : ngsolve.comp.ElementId
 |         input element id
 |
 |  GetOrder(...)
 |      GetOrder(self: ngsolve.comp.FESpace, nodeid: ngsolve.comp.NodeId) -> int
 |
 |      return order of node.
 |      by now, only isotropic order is supported here
 |
 |  GetTrace(...)
 |      GetTrace(self: ngsolve.comp.FESpace, arg0: ngsolve.comp.FESpace, arg1: ngsolve.la.BaseVector, arg2: ngsolve.la.BaseVector, arg3: bool) -> None
 |
 |  GetTraceTrans(...)
 |      GetTraceTrans(self: ngsolve.comp.FESpace, arg0: ngsolve.comp.FESpace, arg1: ngsolve.la.BaseVector, arg2: ngsolve.la.BaseVector, arg3: bool) -> None
 |
 |  HideAllDofs(...)
 |      HideAllDofs(self: ngsolve.comp.FESpace, component: object = <ngsolve.ngstd.DummyArgument>) -> None
 |
 |      set all visible coupling types to HIDDEN_DOFs (will be overwritten by any Update())
 |
 |  InvM(...)
 |      InvM(self: ngsolve.comp.FESpace, rho: ngsolve.fem.CoefficientFunction = None) -> BaseMatrix
 |
 |  Mass(...)
 |      Mass(self: ngsolve.comp.FESpace, rho: ngsolve.fem.CoefficientFunction = None, definedon: Optional[ngsolve.comp.Region] = None) -> BaseMatrix
 |
 |  ParallelDofs(...)
 |      ParallelDofs(self: ngsolve.comp.FESpace) -> ngsolve.la.ParallelDofs
 |
 |      Return dof-identification for MPI-distributed meshes
 |
 |  Prolongation(...)
 |      Prolongation(self: ngsolve.comp.FESpace) -> ngmg::Prolongation
 |
 |      Return prolongation operator for use in multi-grid
 |
 |  Range(...)
 |      Range(self: ngsolve.comp.FESpace, arg0: int) -> ngsolve.la.DofRange
 |
 |      deprecated, will be only available for ProductSpace
 |
 |  SetCouplingType(...)
 |      SetCouplingType(*args, **kwargs)
 |      Overloaded function.
 |
 |      1. SetCouplingType(self: ngsolve.comp.FESpace, dofnr: int, coupling_type: ngsolve.comp.COUPLING_TYPE) -> None
 |
 |
 |               Set coupling type of a degree of freedom.
 |
 |      Parameters:
 |
 |      dofnr : int
 |        input dof number
 |
 |      coupling_type : ngsolve.comp.COUPLING_TYPE
 |        input coupling type
 |
 |
 |
 |      2. SetCouplingType(self: ngsolve.comp.FESpace, dofnrs: ngsolve.ngstd.IntRange, coupling_type: ngsolve.comp.COUPLING_TYPE) -> None
 |
 |
 |               Set coupling type for interval of dofs.
 |
 |      Parameters:
 |
 |      dofnrs : Range
 |        range of dofs
 |
 |      coupling_type : ngsolve.comp.COUPLING_TYPE
 |        input coupling type
 |
 |  SetDefinedOn(...)
 |      SetDefinedOn(self: ngsolve.comp.FESpace, region: ngsolve.comp.Region) -> None
 |
 |
 |      Set the regions on which the FESpace is defined.
 |
 |      Parameters:
 |
 |      region : ngsolve.comp.Region
 |        input region
 |
 |  SetOrder(...)
 |      SetOrder(*args, **kwargs)
 |      Overloaded function.
 |
 |      1. SetOrder(self: ngsolve.comp.FESpace, element_type: ngsolve.fem.ET, order: int) -> None
 |
 |
 |
 |      Parameters:
 |
 |      element_type : ngsolve.fem.ET
 |        input element type
 |
 |      order : object
 |        input polynomial order
 |
 |
 |      2. SetOrder(self: ngsolve.comp.FESpace, nodeid: ngsolve.comp.NodeId, order: int) -> None
 |
 |
 |
 |      Parameters:
 |
 |      nodeid : ngsolve.comp.NodeId
 |        input node id
 |
 |      order : int
 |        input polynomial order
 |
 |  SolveM(...)
 |      SolveM(self: ngsolve.comp.FESpace, vec: ngsolve.la.BaseVector, rho: ngsolve.fem.CoefficientFunction = None, definedon: ngsolve.comp.Region = None) -> None
 |
 |
 |               Solve with the mass-matrix. Available only for L2-like spaces.
 |
 |      Parameters:
 |
 |      vec : ngsolve.la.BaseVector
 |        input right hand side vector
 |
 |      rho : ngsolve.fem.CoefficientFunction
 |        input CF
 |
 |  TestFunction(...)
 |      TestFunction(self: ngsolve.comp.FESpace) -> object
 |
 |      Return a proxy to be used as a testfunction for Symbolic Integrators
 |
 |  TnT(...)
 |      TnT(self: ngsolve.comp.FESpace) -> Tuple[object, object]
 |
 |      Return a tuple of trial and testfunction
 |
 |  TraceOperator(...)
 |      TraceOperator(self: ngsolve.comp.FESpace, tracespace: ngsolve.comp.FESpace, average: bool) -> BaseMatrix
 |
 |  TrialFunction(...)
 |      TrialFunction(self: ngsolve.comp.FESpace) -> object
 |
 |      Return a proxy to be used as a trialfunction in Symbolic Integrators
 |
 |  Update(...)
 |      Update(self: ngsolve.comp.FESpace) -> None
 |
 |      update space after mesh-refinement
 |
 |  UpdateDofTables(...)
 |      UpdateDofTables(self: ngsolve.comp.FESpace) -> None
 |
 |      update dof-tables after changing polynomial order distribution
 |
 |  __eq__(...)
 |      __eq__(self: ngsolve.comp.FESpace, space: ngsolve.comp.FESpace) -> bool
 |
 |  __mul__(...)
 |      __mul__(self: ngsolve.comp.FESpace, arg0: ngsolve.comp.FESpace) -> ngcomp::CompoundFESpace
 |
 |  __pow__(...)
 |      __pow__(self: ngsolve.comp.FESpace, arg0: int) -> ngcomp::CompoundFESpaceAllSame
 |
 |  __str__(...)
 |      __str__(self: ngsolve.comp.FESpace) -> str
 |
 |  __timing__(...)
 |      __timing__(self: ngsolve.comp.FESpace) -> object
 |
 |  ----------------------------------------------------------------------
 |  Static methods inherited from FESpace:
 |
 |  __special_treated_flags__(...) method of builtins.PyCapsule instance
 |      __special_treated_flags__() -> dict
 |
 |  ----------------------------------------------------------------------
 |  Readonly properties inherited from FESpace:
 |
 |  autoupdate
 |
 |  components
 |      deprecated, will be only available for ProductSpace
 |
 |  couplingtype
 |
 |  dim
 |      multi-dim of FESpace
 |
 |  globalorder
 |      query global order of space
 |
 |  is_complex
 |
 |  loembedding
 |
 |  lospace
 |
 |  mesh
 |      mesh on which the FESpace is created
 |
 |  ndof
 |      number of degrees of freedom
 |
 |  ndofglobal
 |      global number of dofs on MPI-distributed mesh
 |
 |  type
 |      type of finite element space
 |
 |  ----------------------------------------------------------------------
 |  Data and other attributes inherited from FESpace:
 |
 |  __hash__ = None
 |
 |  ----------------------------------------------------------------------
 |  Readonly properties inherited from NGS_Object:
 |
 |  __memory__
 |
 |  flags
 |
 |  ----------------------------------------------------------------------
 |  Data descriptors inherited from NGS_Object:
 |
 |  name
 |
 |  ----------------------------------------------------------------------
 |  Static methods inherited from pybind11_builtins.pybind11_object:
 |
 |  __new__(*args, **kwargs) class method of pybind11_builtins.pybind11_object
 |      Create and return a new object.  See help(type) for accurate signature.