Skip to content

Equations

equations

This module contains abstract classes to define the structure of mathematical terms and equations within the G-ADOPT library. Users should not interact with these classes; instead, please use the solvers provided in other modules.

BaseEquation(test_space, trial_space, quad_degree=None, **kwargs)

Produces the UFL for the registered terms constituting an equation.

The equation instance is initialised given test and trial function spaces. Test and trial spaces are only used to determine the employed discretisation (i.e. UFL elements); test and trial functions are provided separately in the residual method.

Keyword arguments provided here are passed on to each collected equation term.

Attributes:

Name Type Description
terms list[BaseTerm]

List of equation terms defined through inheritance from BaseTerm.

Parameters:

Name Type Description Default
test_space WithGeometry

Firedrake function space of the test function

required
trial_space WithGeometry

Firedrake function space of the trial function

required
quad_degree Optional[int]

Quadrature degree. Default value is 2p + 1, where p the polynomial degree of the trial space

None
Source code in g-adopt/gadopt/equations.py
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
def __init__(
    self,
    test_space: firedrake.functionspaceimpl.WithGeometry,
    trial_space: firedrake.functionspaceimpl.WithGeometry,
    quad_degree: Optional[int] = None,
    **kwargs
):
    self.test_space = test_space
    self.trial_space = trial_space
    self.mesh = trial_space.mesh()

    p = trial_space.ufl_element().degree()
    if isinstance(p, int):  # isotropic element
        if quad_degree is None:
            quad_degree = 2*p + 1
    else:  # tensorproduct element
        p_h, p_v = p
        if quad_degree is None:
            quad_degree = 2*max(p_h, p_v) + 1

    if trial_space.extruded:
        # Create surface measures that treat the bottom and top boundaries similarly
        # to lateral boundaries. This way, integration using the ds and dS measures
        # occurs over both horizontal and vertical boundaries, and we can also use
        # "bottom" and "top" as surface identifiers, for example, ds("top").
        self.ds = CombinedSurfaceMeasure(self.mesh, quad_degree)
        self.dS = (
            firedrake.dS_v(domain=self.mesh, degree=quad_degree) +
            firedrake.dS_h(domain=self.mesh, degree=quad_degree)
        )
    else:
        self.ds = firedrake.ds(domain=self.mesh, degree=quad_degree)
        self.dS = firedrake.dS(domain=self.mesh, degree=quad_degree)

    self.dx = firedrake.dx(domain=self.mesh, degree=quad_degree)

    # self._terms stores the actual instances of the BaseTerm-classes in self.terms
    self._terms = []
    for TermClass in self.terms:
        self._terms.append(TermClass(test_space, trial_space, self.dx, self.ds, self.dS, **kwargs))

mass_term(test, trial)

Typical mass term used in time discretisations.

Parameters:

Name Type Description Default
test Argument

Firedrake test function

required
trial Argument | Function

Firedrake trial function

required

Returns:

Type Description
Expr

The UFL expression associated with the mass term of the equation.

Source code in g-adopt/gadopt/equations.py
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
def mass_term(
    self,
    test: firedrake.ufl_expr.Argument,
    trial: firedrake.ufl_expr.Argument | firedrake.Function,
) -> firedrake.ufl.core.expr.Expr:
    """Typical mass term used in time discretisations.

    Arguments:
      test:  Firedrake test function
      trial: Firedrake trial function

    Returns:
      The UFL expression associated with the mass term of the equation.

    """
    return firedrake.inner(test, trial) * self.dx

residual(test, trial, trial_lagged=None, fields=None, bcs=None)

Finite element residual.

The final residual is calculated as a sum of all individual term residuals.

Parameters:

Name Type Description Default
test Argument

Firedrake test function

required
trial Argument | Function

Firedrake trial function

required
trial_lagged Optional[Argument | Function]

Firedrake trial function from the previous time step

None
fields Optional[dict[str, Constant | Function]]

Dictionary of physical fields from the simulation's state

None
bcs Optional[dict[int, dict[str, int | float]]]

Dictionary of identifier-value pairs specifying boundary conditions

None

Returns:

Type Description
Expr

The UFL expression associated with all equation terms except the mass term.

Source code in g-adopt/gadopt/equations.py
 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
def residual(
    self,
    test: firedrake.ufl_expr.Argument,
    trial: firedrake.ufl_expr.Argument | firedrake.Function,
    trial_lagged: Optional[firedrake.ufl_expr.Argument | firedrake.Function] = None,
    fields: Optional[dict[str, firedrake.Constant | firedrake.Function]] = None,
    bcs: Optional[dict[int, dict[str, int | float]]] = None,
) -> firedrake.ufl.core.expr.Expr:
    """Finite element residual.

    The final residual is calculated as a sum of all individual term residuals.

    Arguments:
      test:         Firedrake test function
      trial:        Firedrake trial function
      trial_lagged: Firedrake trial function from the previous time step
      fields:       Dictionary of physical fields from the simulation's state
      bcs:          Dictionary of identifier-value pairs specifying boundary conditions

    Returns:
      The UFL expression associated with all equation terms except the mass term.

    """
    if trial_lagged is None:
        trial_lagged = trial
    if fields is None:
        fields = {}
    if bcs is None:
        bcs = {}

    F = 0
    for term in self._terms:
        F += term.residual(test, trial, trial_lagged, fields, bcs)

    return F

BaseTerm(test_space, trial_space, dx, ds, dS, **kwargs)

Bases: ABC

Defines an equation's term using an UFL expression.

The implemented expression describes the term's contribution to the residual in the finite element discretisation.

Parameters:

Name Type Description Default
test_space WithGeometry

Firedrake function space of the test function

required
trial_space WithGeometry

Firedrake function space of the trial function

required
dx Measure

UFL measure for the domain, boundaries excluded

required
ds Measure

UFL measure for the domain's outer boundaries

required
dS Measure

UFL measure for the domain's inner boundaries when using a discontinuous function space

required
Source code in g-adopt/gadopt/equations.py
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
def __init__(
    self,
    test_space: firedrake.functionspaceimpl.WithGeometry,
    trial_space: firedrake.functionspaceimpl.WithGeometry,
    dx: firedrake.Measure,
    ds: firedrake.Measure,
    dS: firedrake.Measure,
    **kwargs,
):
    self.test_space = test_space
    self.trial_space = trial_space

    self.dx = dx
    self.ds = ds
    self.dS = dS

    self.mesh = test_space.mesh()
    self.dim = self.mesh.geometric_dimension()
    self.n = firedrake.FacetNormal(self.mesh)

    self.term_kwargs = kwargs

residual(test, trial, trial_lagged=None, fields=None, bcs=None) abstractmethod

Residual associated with the equation's term.

Parameters:

Name Type Description Default
test Argument

Firedrake test function

required
trial Argument | Function

Firedrake trial function

required
trial_lagged Optional[Argument | Function]

Firedrake trial function from the previous time step

None
fields Optional[dict[str, Constant | Function]]

Dictionary of physical fields from the simulation's state

None
bcs Optional[dict[int, dict[str, int | float]]]

Dictionary of identifier-value pairs specifying boundary conditions

None

Returns:

Type Description
Expr

A UFL expression for the term's contribution to the finite element residual.

Source code in g-adopt/gadopt/equations.py
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
@abstractmethod
def residual(
    self,
    test: firedrake.ufl_expr.Argument,
    trial: firedrake.ufl_expr.Argument | firedrake.Function,
    trial_lagged: Optional[firedrake.ufl_expr.Argument | firedrake.Function] = None,
    fields: Optional[dict[str, firedrake.Constant | firedrake.Function]] = None,
    bcs: Optional[dict[int, dict[str, int | float]]] = None,
) -> firedrake.ufl.core.expr.Expr:
    """Residual associated with the equation's term.

    Arguments:
      test:         Firedrake test function
      trial:        Firedrake trial function
      trial_lagged: Firedrake trial function from the previous time step
      fields:       Dictionary of physical fields from the simulation's state
      bcs:          Dictionary of identifier-value pairs specifying boundary conditions

    Returns:
      A UFL expression for the term's contribution to the finite element residual.

    """
    pass