Meshes in 1DΒΆ

We will start with the simplest case: 1D domain, $(a,b) \subset \mathbb{R}$, decomposed into $N$ line segments.

InΒ [1]:
from import_hack import *
from methodsnm.mesh_1d import *
from methodsnm.visualize import DrawMesh1D
m = Mesh1D([0,1],10)
DrawMesh1D(m)
source module for methodsNM imported.
No description has been provided for this image
InΒ [2]:
m = Mesh1D([0,0.9,1],10)
DrawMesh1D(m)
No description has been provided for this image
InΒ [3]:
m = m.uniform_refine()
DrawMesh1D(m)
No description has been provided for this image

Mesh functionsΒΆ

We use the term MeshFunction for functions that can be evaluated on a mesh using a combination of

  • a mesh element and
  • a (local) integration point.

To connect the local point of view (integration point) and the global/mesh view we use a transformation from reference element to physical element

In 1D: $\Phi_T: (0,1) \to (a,b)$ where $(a,b)$ is the domain of an element $T$ of the mesh.

There are essentially two types of functions that we want to unify with a MeshFunction:

  • Functions given in world coordinate (not directly related to the mesh)
  • Finite element based functions that are naturally described by integration point and element index

Let's take a look at the implementations of MeshFunctions:

Example of a GlobalFunctionΒΆ

InΒ [4]:
from methodsnm.meshfct import *
from methodsnm.visualize import DrawFunction1D
from math import sin, pi
m = Mesh1D([0,1],3)
uex = GlobalFunction(function=lambda x: sin(4*x*pi), mesh=m)
DrawFunction1D(uex, sampling=25, show_mesh=True)
No description has been provided for this image

Next, we make the connection between local and global for finite elements.

Finite element spaces and the global to local mapΒΆ

  • Discrete objects so far were only defined on the reference element.
  • With the help of the transformation we can define these on any physical element.
  • Next, we need to sum up the contributions and manage the related dofs

Dof handlingΒΆ

The finite element space (FESpace) basically translates local to global dofs.

We again start by a simple example: piecewise linear elements on a 1D mesh:

InΒ [5]:
from methodsnm.fes import *
m = Mesh1D([0,1],10)
fes = P1_Segments_Space(m)
for elnr,el in enumerate(m.elements()):
    print("element number:", elnr, end="\t")
    print("element vertices (indices):", el, end="\t")
    print("global dofs of element:", fes.element_dofs(elnr))
element number: 0	element vertices (indices): [0 1]	global dofs of element: [0 1]
element number: 1	element vertices (indices): [1 2]	global dofs of element: [1 2]
element number: 2	element vertices (indices): [2 3]	global dofs of element: [2 3]
element number: 3	element vertices (indices): [3 4]	global dofs of element: [3 4]
element number: 4	element vertices (indices): [4 5]	global dofs of element: [4 5]
element number: 5	element vertices (indices): [5 6]	global dofs of element: [5 6]
element number: 6	element vertices (indices): [6 7]	global dofs of element: [6 7]
element number: 7	element vertices (indices): [7 8]	global dofs of element: [7 8]
element number: 8	element vertices (indices): [8 9]	global dofs of element: [8 9]
element number: 9	element vertices (indices): [ 9 10]	global dofs of element: [ 9 10]

TasksΒΆ

Task FES1D-1ΒΆ

Add a parameter to the P1_Segments_Space class that identifies the most left and the most right dof, rending the functions periodic. Adjust ndof and the global to local map accordingly.

  • What are the ndofs now?
  • What is the global to local map now?
  • Implement the changes directly in P1_Segments_Space and test your implementation with the following code:
InΒ [6]:
m = Mesh1D([0,1],10)
fes = P1_Segments_Space(m, periodic = True)
for elnr,el in enumerate(m.elements()):
    print("element number:", elnr, end="\t")
    print("element vertices (indices):", el, end="\t")
    print("global dofs of element:", fes.element_dofs(elnr))
element number: 0	element vertices (indices): [0 1]	global dofs of element: [0 1]
element number: 1	element vertices (indices): [1 2]	global dofs of element: [1 2]
element number: 2	element vertices (indices): [2 3]	global dofs of element: [2 3]
element number: 3	element vertices (indices): [3 4]	global dofs of element: [3 4]
element number: 4	element vertices (indices): [4 5]	global dofs of element: [4 5]
element number: 5	element vertices (indices): [5 6]	global dofs of element: [5 6]
element number: 6	element vertices (indices): [6 7]	global dofs of element: [6 7]
element number: 7	element vertices (indices): [7 8]	global dofs of element: [7 8]
element number: 8	element vertices (indices): [8 9]	global dofs of element: [8 9]
element number: 9	element vertices (indices): [ 9 10]	global dofs of element: [np.int64(9), 0]

Global basis functions:ΒΆ

To test the global to local map, let's draw all shape functions of an FESpace:

InΒ [7]:
from methodsnm.meshfct import *
from methodsnm.visualize import DrawFunction1D, DrawShapes
InΒ [8]:
m = Mesh1D([0,1],5)
DrawShapes(P1_Segments_Space(m))
DrawShapes(P1_Segments_Space(m,periodic=True))
No description has been provided for this image
No description has been provided for this image

Task FES1D-2ΒΆ

Implement a finite element space of discontinuous functions corresponding to either (or all):

  • the P1 finite element or
  • the P2 finite element or
  • a new P0 element of element-wise constant functions

Task:

  • Think about the global number of unknowns
  • Think about the global to local dof map
  • Implement the FESpace and draw the basis functions
InΒ [9]:
m = Mesh1D([0,1],3)
DrawShapes(P1disc_Segments_Space(m),sampling=3)
No description has been provided for this image

Task FES1D-3ΒΆ

Implement a higher order Lagrange finite element space and draw the basis functions.

InΒ [10]:
m = Mesh1D([0,1],3)
DrawShapes(Lagrange_Segments_Space(m,3),sampling=25)
No description has been provided for this image
InΒ [11]:
m = Mesh1D([0,1],3)
DrawShapes(Pk_IntLeg_Segments_Space(m,5),sampling=25)
No description has been provided for this image

Task FES1D-4 (optional)ΒΆ

Implement an additional function interpolate that takes another MeshFunction and interpolates the values based on the Lagrange interpolation:

$$ u_h(x) = I_h(v) = \sum_{i=1}^{\texttt{ndof}} \Psi_i(v) \phi_i(x) $$ where $\Psi_i$ are the functionals that are defined by the degrees of freedom. Note that global and local functionals have a correspondence through the finite element space. In the case of the Lagrange space, the functionals are point evaluations.