File: documentation.rst

package info (click to toggle)
dolfin 2019.2.0~legacy20240219.1c52e83-18
  • links: PTS, VCS
  • area: main
  • in suites: trixie
  • size: 31,700 kB
  • sloc: xml: 104,040; cpp: 102,227; python: 24,356; sh: 460; makefile: 330; javascript: 226
file content (155 lines) | stat: -rw-r--r-- 5,270 bytes parent folder | download | duplicates (5)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
.. Documentation for the Stokes mini demo from DOLFIN.

.. _demo_pde_stokes-mini_python_documentation:

Stokes equations with Mini elements
===================================

This demo is implemented in a single Python file,
:download:`demo_stokes-mini.py`, which contains both the variational
form and the solver.

.. include:: ../common.txt

Implementation
--------------

This description goes through the implementation (in
:download:`demo_stokes-mini.py`) of a solver for the Stokes equation
step-by-step.

First, the :py:mod:`dolfin` module is imported:

.. code-block:: python

    from dolfin import *

In this example, different boundary conditions are prescribed on
different parts of the boundaries. This information must be made
available to the solver. One way of doing this, is to tag the
different sub-regions with different (integer) labels. DOLFIN provides
a class :py:class:`MeshFunction <dolfin.cpp.mesh.MeshFunction>` which
is useful for these types of operations: instances of this class
represent a functions over mesh entities (such as over cells or over
facets). Mesh and mesh functions can be read from file in the
following way:

.. code-block:: python

    # Load mesh and subdomains
    mesh = Mesh("../dolfin_fine.xml.gz")
    sub_domains = MeshFunction("size_t", mesh, "../dolfin_fine_subdomains.xml.gz")

Next, we define a :py:class:`FunctionSpace
<dolfin.functions.functionspace.FunctionSpace>` on Mini element
``V * Q``. UFL object ``V = VectorElement(NodalEnrichedElement(P1, B))`` stands for
the vectorial Lagrange element of degree 1 enriched with the cubic
Bubble. ``V * Q`` defines the mixed element for velocity and pressure.

.. code-block:: python

    # Build function spaces on Mini element
    P1 = FiniteElement("Lagrange", mesh.ufl_cell(), 1)
    B = FiniteElement("Bubble",   mesh.ufl_cell(), mesh.topology().dim() + 1)
    V = VectorElement(NodalEnrichedElement(P1, B))
    Q = P1
    W = FunctionSpace(mesh, V * Q)

Now that we have our mixed function space and marked subdomains
defining the boundaries, we define boundary conditions:

.. code-block:: python

    # No-slip boundary condition for velocity
    bc0 = DirichletBC(W.sub(0), (0, 0), sub_domains, 0)

    # Inflow boundary condition for velocity
    inflow = Expression(("-sin(x[1]*pi)", "0.0"), degree=2)
    bc1 = DirichletBC(W.sub(0), inflow, sub_domains, 1)

    # Collect boundary conditions
    bcs = [bc0, bc1]

Here, we have given four arguments in the call of
:py:class:`DirichletBC <dolfin.cpp.fem.DirichletBC>`. The first
specifies the :py:class:`FunctionSpace
<dolfin.cpp.function.FunctionSpace>`. Since we have a
mixed function space, we write
``W.sub(0)`` for the velocity componenet of the space,
and ``W.sub(1)`` for the pressure component of the space. The second
argument specifies the value on the Dirichlet boundary. The two last
arguments specify the marking of the subdomains; ``sub_domains`` contains
the subdomain markers and the number given as the last argument is the
subdomain index.

The bilinear and linear forms corresponding to the weak mixed
formulation of the Stokes equations are defined as follows:

.. code-block:: python

    # Define variational problem
    u, p = TrialFunctions(W)
    v, q = TestFunctions(W)
    f = Constant((0, 0))
    a = (inner(grad(u), grad(v)) - div(v)*p - q*div(u))*dx
    L = inner(f, v)*dx

To compute the solution we use the bilinear and linear forms, and the
boundary condition, but we also need to create a :py:class:`Function
<dolfin.cpp.function.Function>` to store the solution(s). The (full)
solution will be stored in ``w``, which we initialize using the mixed
function space ``W``. The actual
computation is performed by calling solve with the arguments ``a``,
``L``, ``w`` and ``bcs``. The separate components ``u`` and ``p`` of
the solution can be extracted by calling the :py:meth:`split
<dolfin.functions.function.Function.split>` function. Here we use an
optional argument True in the split function to specify that we want a
deep copy. If no argument is given we will get a shallow copy. We want
a deep copy for further computations on the coefficient vectors.

.. code-block:: python

    # Compute solution
    w = Function(W)
    solve(a == L, w, bcs)

    # Split the mixed solution using deepcopy
    # (needed for further computation on coefficient vector)
    u, p = w.split(True)

We may be interested in the :math:`l^2` norms of u and p, they can be
calculated and printed by writing

.. code-block:: python

    print("Norm of velocity coefficient vector: %.15g" % u.vector().norm("l2"))
    print("Norm of pressure coefficient vector: %.15g" % p.vector().norm("l2"))

One can also split functions using shallow copies (which is enough
when we just plot the result) by writing

.. code-block:: python

    # Split the mixed solution using a shallow copy
    u, p = w.split()

Finally, we can store to file and plot the solutions.

.. code-block:: python

    # Save solution in VTK format
    ufile_pvd = File("velocity.pvd")
    ufile_pvd << u
    pfile_pvd = File("pressure.pvd")
    pfile_pvd << p

    # Plot solution
    plot(u)
    plot(p)
    interactive()

Complete code
-------------

.. literalinclude:: demo_stokes-mini.py
   :start-after: # Begin demo