Diagnostics
diagnostics
This module provides a class to simplify computing of diagnostics typically encountered in geodynamical simulations. Users instantiate the class by providing relevant parameters and call individual class methods to compute associated diagnostics.
FunctionContext(quad_degree, func)
Hold objects that can be derived from a Firedrake Function.
This class gathers references to objects that can be pulled from a Firedrake
Function object and calculates quantities based on those objects that will remain
constant for the duration of a simulation. The set of objects/quantities stored
are: mesh, function_space, dx and ds measures, the FacetNormal of the mesh
(as the .normal attribute) and the volume of the domain.
Typical usage example:
function_contexts[F] = FunctionContext(quad_degree, F)
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
quad_degree
|
int
|
Quadrature degree to use when approximating integrands involving |
required |
func
|
Function
|
Function |
required |
Source code in g-adopt/gadopt/diagnostics.py
66 67 68 | |
function
cached
property
The function associated with the instance
mesh
cached
property
The mesh on which the function has been defined
function_space
cached
property
The function space on which the function has been defined
dx
cached
property
The volume integration measure defined by the mesh and
quad_degree passed when creating this instance
ds
cached
property
The surface integration measure defined by the mesh and
quad_degree passed when creating this instance
normal
cached
property
The facet normal of the mesh belonging to this instance
volume
cached
property
The volume of the mesh belonging to this instance
boundary_ids
cached
property
The boundary IDs of the mesh associated with this instance
surface_area
cached
property
The surface area of all surfaces of the mesh belonging to this instance
get_boundary_nodes(boundary_id)
cached
Return the list of nodes on the boundary owned by this process
Creates a DirichletBC object, then uses the .nodes attribute for that
object to provide a list of indices that reside on the boundary of the domain
of the function associated with this FunctionContext instance. The
dof_dset.size parameter of the FunctionSpace is used to exclude nodes in
the halo region of the domain.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
boundary_id
|
int
|
Integer ID of the domain boundary |
required |
Returns:
| Type | Description |
|---|---|
list[int]
|
List of integers corresponding to nodes on the boundary identified by |
list[int]
|
|
Source code in g-adopt/gadopt/diagnostics.py
121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 | |
BaseDiagnostics(quad_degree, **funcs)
A base class containing useful operations for diagnostics
For each Firedrake function passed as a keyword argument in the funcs parameter,
store that function as an attribute of the class accessible by its keyword, e.g.:
diag = BaseDiagnostics(quad_degree, z=z)
sets the Firedrake function z to the diag.z parameter.
If the function is a MixedFunction, the subfunctions will be accessible by an
index, e.g.:
diag = BaseDiagnostics(quad_degree, z=z)
sets the subfunctions of z to diag.z_0, diag.z_1, etc. A FunctionContext
is created for each function. These attributes are accessed by the
diag._function_contexts dict.
This class is intended to be subclassed by domain-specific diagnostic classes
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
quad_degree
|
int
|
Quadrature degree to use when approximating integrands managed by |
required |
**funcs
|
Function | None
|
Firedrake functions to associate with this instance |
{}
|
Initialise a BaseDiagnostics object.
Sets the quad_degree for measures used by this object and passes the
remaining keyword arguments through to register_functions.
Source code in g-adopt/gadopt/diagnostics.py
168 169 170 171 172 173 174 175 176 177 | |
register_functions(*, quad_degree=None, **funcs)
Register a function with this BaseDiagnostics object.
Creates a FunctionContext object for each function passed in as a keyword
argument. Also creates an attribute on the instance to access the input function
named for the key of the keyword argument. i.e:
> diag.register_functions(self, F=F)
> type(diag.F)
<class 'firedrake.function.Function'>
None, the attribute will still be created
but set to 0.0. If a mixed function is entered, each subfunction will have
a FunctionContext object associated with it, and the attribute will be named
with an additional number to denote the index of the subfunction i.e.:
```
diag.register_functions(self, F) type(diag.F) AttributeError: 'Requested 'F', which lives on a mixed space. Instead, access subfunctions via F_0, F_1, ..." type(diag.F_0)
type(diag.F_1)
Args:
quad_degree (optional): The quadrature degree for the measures to be used
by this function. If None, the quad_degree passed at object
instantiation time is used. Defaults to None.
**funcs: key-value pairs of Firedrake functions to associate with this
instance
Source code in g-adopt/gadopt/diagnostics.py
179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 | |
get_upward_component(f)
cached
Get the upward (against gravity) component of a function.
Returns a UFL expression for the upward component of a function. Uses the
G-ADOPT vertical_component function and caches the result such that the
UFL expression only needs to be constructed once per run.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
f
|
Function
|
Function |
required |
Returns:
| Type | Description |
|---|---|
Operator
|
UFL expression for the vertical component of |
Source code in g-adopt/gadopt/diagnostics.py
420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 | |
min(func_or_op, boundary_id=None, dim=None)
Calculate the minimum value of a function. See _minmax
docstring for more information.
Source code in g-adopt/gadopt/diagnostics.py
495 496 497 498 499 500 501 502 503 504 505 | |
max(func_or_op, boundary_id=None, dim=None)
Calculate the maximum value of a function See _minmax
docstring for more information.
Source code in g-adopt/gadopt/diagnostics.py
507 508 509 510 511 512 513 514 515 516 517 | |
integral(f, boundary_id=None)
Calculate the integral of a function over the domain associated with it
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
f
|
Function
|
Function. |
required |
boundary_id
|
optional
|
Boundary ID. If not provided or set to |
None
|
Returns:
| Type | Description |
|---|---|
float
|
Result of integration |
Source code in g-adopt/gadopt/diagnostics.py
519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 | |
l1norm(f, boundary_id=None)
Calculate the L1norm of a function over the domain associated with it
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
f
|
Function
|
Function. |
required |
boundary_id
|
optional
|
Boundary ID .If not provided or set to |
None
|
Returns:
| Name | Type | Description |
|---|---|---|
float |
float
|
L1 norm |
Source code in g-adopt/gadopt/diagnostics.py
535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 | |
l2norm(f, boundary_id=None)
Calculate the L2norm of a function over the domain associated with it
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
f
|
Function
|
Function. |
required |
boundary_id
|
optional
|
Boundary ID. If not provided or set to |
None
|
Returns:
| Name | Type | Description |
|---|---|---|
float |
float
|
L2 norm |
Source code in g-adopt/gadopt/diagnostics.py
551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 | |
rms(f)
Calculate the RMS of a function over the domain associated with it
For the purposes of this function, RMS is defined as L2norm/volume
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
f
|
Function
|
Function. |
required |
boundary_id
|
optional
|
Boundary ID. If not provided or set to |
required |
Returns:
| Name | Type | Description |
|---|---|---|
float |
float
|
RMS |
Source code in g-adopt/gadopt/diagnostics.py
567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 | |
GeodynamicalDiagnostics(z, T=None, /, bottom_id=None, top_id=None, *, quad_degree=4)
Bases: BaseDiagnostics
Typical simulation diagnostics used in geodynamical simulations.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
z
|
Function
|
Firedrake function for mixed Stokes function space (velocity, pressure) |
required |
T
|
Function | None
|
Firedrake function for temperature |
None
|
bottom_id
|
int | str | None
|
Bottom boundary identifier |
None
|
top_id
|
int | str | None
|
Top boundary identifier |
None
|
quad_degree
|
int
|
Degree of polynomial quadrature approximation |
4
|
Note
All diagnostics are returned as floats.
Methods:
| Name | Description |
|---|---|
u_rms |
Root-mean-square velocity |
u_rms_top |
Root-mean-square velocity along the top boundary |
Nu_top |
Nusselt number at the top boundary |
Nu_bottom |
Nusselt number at the bottom boundary |
T_avg |
Average temperature in the domain |
T_min |
Minimum temperature in domain |
T_max |
Maximum temperature in domain |
ux_max |
Maximum velocity (first component, optionally over a given boundary) |
Source code in g-adopt/gadopt/diagnostics.py
609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 | |
GIADiagnostics(u, /, bottom_id=None, top_id=None, *, quad_degree=4)
Bases: BaseDiagnostics
Typical simulation diagnostics used in glacial isostatic adjustment simulations.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
d
|
Firedrake function for displacement |
required | |
bottom_id
|
int | str | None
|
Bottom boundary identifier |
None
|
top_id
|
int | str | None
|
Top boundary identifier |
None
|
quad_degree
|
int
|
Degree of polynomial quadrature approximation |
4
|
Note
All diagnostics are returned as floats.
Methods:
| Name | Description |
|---|---|
u_rms |
Root-mean-square displacement |
u_rms_top |
Root-mean-square displacement along the top boundary |
ux_max |
Maximum displacement (first component, optionally over a given boundary) |
uv_min |
Minimum vertical displacement, optionally over a given boundary |
uv_max |
Maximum vertical displacement, optionally over a given boundary |
l2_norm_top |
L2 norm of displacement on top surface |
l1_norm_top |
L1 norm of displacement on top surface |
integrated_displacement |
integral of displacement on top surface |
Source code in g-adopt/gadopt/diagnostics.py
691 692 693 694 695 696 697 698 699 700 701 702 703 704 705 706 | |
uv_min(boundary_id=None)
Minimum value of vertical component of velocity/displacement
Source code in g-adopt/gadopt/diagnostics.py
717 718 719 | |
uv_max(boundary_id=None)
Maximum value of vertical component of velocity/displacement
Source code in g-adopt/gadopt/diagnostics.py
721 722 723 | |