Utility
utility
This module provides several classes and functions to perform a number of pre-, syn-, and post-processing tasks. Users incorporate utility as required in their code, depending on what they would like to achieve.
TimestepAdaptor(dt_const, u, V, target_cfl=1.0, increase_tolerance=1.5, maximum_timestep=None)
Computes timestep based on CFL condition for provided velocity field
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
dt_const
|
Constant whose value will be updated by the timestep adaptor |
required | |
u
|
Velocity to base CFL condition on |
required | |
V
|
FunctionSpace for reference velocity, usually velocity space |
required | |
target_cfl
|
CFL number to target with chosen timestep |
1.0
|
|
increase_tolerance
|
Maximum tolerance timestep is allowed to change by |
1.5
|
|
maximum_timestep
|
Maximum allowable timestep |
None
|
Source code in g-adopt/gadopt/utility.py
67 68 69 70 71 72 73 74 75 76 77 78 | |
CombinedSurfaceMeasure(domain, degree)
Bases: Measure
A surface measure that combines ds_v, the integral over vertical boundary facets, and ds_t and ds_b, the integral over horizontal top and bottom facets. The vertical boundary facets are identified with the same surface ids as ds_v. The top and bottom surfaces are identified via the "top" and "bottom" ids.
Source code in g-adopt/gadopt/utility.py
139 140 141 142 | |
__rmul__(other)
This is to handle terms to be integrated over all surfaces in the form of other*ds. Here the CombinedSurfaceMeasure ds is not called, instead we just split it up as below.
Source code in g-adopt/gadopt/utility.py
152 153 154 155 | |
ExtrudedFunction(*args, mesh_3d=None, **kwargs)
Bases: Function
A 2D Function that provides a 3D view on the extruded domain.
The 3D function can be accessed as ExtrudedFunction.view_3d.
The 3D function resides in V x R function space, where V is the function
space of the source function. The 3D function shares the data of the 2D
function.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
mesh_3d
|
Extruded 3D mesh where the function will be extended to. |
None
|
Source code in g-adopt/gadopt/utility.py
254 255 256 257 258 259 | |
LayerAveraging(mesh, r1d=None, quad_degree=None)
A manager for computing a vertical profile of horizontal layer averages.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
mesh
|
The mesh over which to compute averages |
required | |
r1d
|
An array of either depth coordinates or radii, at which to compute layer averages. If not provided, and mesh is extruded, it uses the same layer heights. If mesh is not extruded, r1d is required. |
None
|
Source code in g-adopt/gadopt/utility.py
313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 | |
get_layer_average(T)
Compute the layer averages of Function T at the predefined depths.
Returns:
| Type | Description |
|---|---|
|
A numpy array containing the averages. |
Source code in g-adopt/gadopt/utility.py
394 395 396 397 398 399 400 401 402 | |
extrapolate_layer_average(u, avg)
Given an array of layer averages avg, extrapolate to Function u
Source code in g-adopt/gadopt/utility.py
404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 | |
InteriorBC
Bases: DirichletBC
DirichletBC applied to anywhere that is not on the specified boundary
log(*args)
Log output to stdout from root processor only
Source code in g-adopt/gadopt/utility.py
33 34 35 | |
depends_on(ufl_expr, terminal)
Does ufl_expr depend on terminal (Function/Constant/...)?
Source code in g-adopt/gadopt/utility.py
182 183 184 | |
tensor_jump(v, n)
Jump term for vector functions based on the tensor product
This is the discrete equivalent of grad(u) as opposed to the
vectorial UFL jump operator ufl.jump which represents div(u).
The equivalent of nabla_grad(u) is given by tensor_jump(n, u).
Source code in g-adopt/gadopt/utility.py
205 206 207 208 209 210 211 212 213 214 215 | |
extend_function_to_3d(func, mesh_extruded)
Returns a 3D view of a 2D Function on the extruded domain.
The 3D function resides in V x R function space, where V is the function
space of the source function. The 3D function shares the data of the 2D
function.
Source code in g-adopt/gadopt/utility.py
218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 | |
absv(u)
Component-wise absolute value of vector for SU stabilisation
Source code in g-adopt/gadopt/utility.py
450 451 452 | |
node_coordinates(function)
Interpolates mesh coordinates at each node of the provided function.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
function
|
Function
|
A Firedrake function |
required |
Returns:
| Type | Description |
|---|---|
Function
|
A Firedrake function for the interpolated mesh coordinates |
Source code in g-adopt/gadopt/utility.py
464 465 466 467 468 469 470 471 472 473 474 475 | |
interpolate_1d_profile(function, one_d_filename)
Assign a one-dimensional profile to a Function function from a file.
The function reads a one-dimensional profile (e.g., viscosity) together with the
radius/height input for the profile, from a file, broadcasts the two arrays to all
processes, and then interpolates this array onto the function space of function.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
function
|
Function
|
The function onto which the 1D profile will be assigned |
required |
one_d_filename
|
str
|
The path to the file containing the 1D radial profile |
required |
Note
- This is designed to read a file with one process and distribute in parallel with MPI.
- The input file should contain an array of radius/height and an array of values, separated by a comma.
Source code in g-adopt/gadopt/utility.py
478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 | |
extruded_layer_heights(DG0_layers, radii)
Calculates layer heights for extruded mesh using rheological boundary
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
DG0_layers
|
int | list[int]
|
Number of layers per rheological layer |
required |
radii
|
list[float]
|
Radii of rheological boundaries to match when extruding the mesh. Assumes sequence starts with outer radii -> inner radii. |
required |
Returns:
| Type | Description |
|---|---|
list[float]
|
list of layer heights |
Source code in g-adopt/gadopt/utility.py
593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 | |
initialise_background_field(f, background_values, X, radii, shift=0.0)
Initialises discontinuous field with sharp jumps at rheological boundaries
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
f
|
Function
|
Function to interpolate background values into. N.b. |
required |
background_values
|
list[float]
|
Discrete values to interpolate into |
required |
X
|
SpatialCoordinate
|
Spatial coordinates associated with function |
required |
radii
|
list[float]
|
Radii of rheological discontinuities. N.b. length of |
required |
shift
|
float
|
Shift vertical radii by a constant, e.g. when mesh is offset from |
0.0
|
Source code in g-adopt/gadopt/utility.py
622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 | |