animate package

Submodules

animate.adapt module

class MetricBasedAdaptor(mesh, metric, name=None, comm=None)[source]

Bases: AdaptorBase

Class for driving metric-based mesh adaptation.

Parameters:
adapted_mesh[source]

Adapt the mesh with respect to the provided metric.

Returns:

the adapted mesh

Return type:

firedrake.mesh.MeshGeometry

interpolate(f)[source]

Interpolate a Function into the corresponding FunctionSpace defined on the adapted mesh.

Parameters:

f (firedrake.function.Function) – a Function on the initial mesh

Returns:

its interpolation onto the adapted mesh

Return type:

firedrake.function.Function

project(f)[source]

Project a Function into the corresponding FunctionSpace defined on the adapted mesh using conservative projection.

Parameters:

f (firedrake.function.Function) – a Function on the initial mesh

Returns:

its projection onto the adapted mesh

Return type:

firedrake.function.Function

adapt(mesh, *metrics, name=None, serialise=None, remove_checkpoints=True)[source]

Adapt a mesh with respect to a metric and some adaptor parameters.

If multiple metrics are provided, then they are intersected.

Parameters:
  • mesh (firedrake.mesh.MeshGeometry) – mesh to be adapted.

  • metrics (list of RiemannianMetrics) – metrics to guide the mesh adaptation

  • name (str) – name for the adapted mesh

  • serialise (bool) – if True, adaptation is done in serial using firedrake.checkpointing.CheckpointFile`s. Defaults to ``True` if the mesh is 2D, and to False if the mesh is 3D or if the code is already run in serial. This is because parallel adaptation is only supported in 3D.

  • remove_checkpoints (bool) – if True, checkpoint files are deleted after use

Returns:

the adapted mesh

Return type:

MeshGeometry

animate.checkpointing module

get_checkpoint_dir()[source]

Make a temporary directory for checkpointing and return its path.

load_checkpoint(filepath, mesh_name, metric_name, comm=<mpi4py.MPI.Intracomm object>)[source]

Load a metric from a CheckpointFile.

Note that the checkpoint will have to be stored within Animate’s .checkpoints subdirectory.

Parameters:
  • filepath (str) – the path to the checkpoint file

  • mesh_name (str) – the name under which the mesh is saved in the checkpoint file

  • metric_name (str) – the name under which the metric is saved in the checkpoint file

  • comm (mpi4py.MPI.Intracom) – MPI communicator for handling the checkpoint file

Returns:

the metric loaded from the checkpoint

Return type:

animate.metric.RiemannianMetric

save_checkpoint(filepath, metric, metric_name=None, comm=<mpi4py.MPI.Intracomm object>)[source]

Write the metric and underlying mesh to a CheckpointFile.

Note that the checkpoint will be stored within Animate’s .checkpoints subdirectory.

Parameters:
  • filepath (str) – the path of the checkpoint file

  • metric (animate.metric.RiemannianMetric) – the metric to save to the checkpoint

  • metric_name (str) – the name under which to save the metric in the checkpoint file

  • comm (mpi4py.MPI.Intracom) – MPI communicator for handling the checkpoint file

animate.interpolation module

Driver functions for mesh-to-mesh data transfer.

clement_interpolant(source, target_space=None, boundary=False)[source]

Compute the Clement interpolant of a \(\mathbb P0\) source field, i.e. take the volume average over neighbouring cells at each vertex. See [].

Parameters:
  • source – the \(\mathbb P0\) source field

  • target_space – the \(\mathbb P1\) space to interpolate into

  • boundary – interpolate over boundary facets or cells?

interpolate(source, target_space, **kwargs)[source]

Overload function firedrake.__future__.interpolate() to account for the case of two mixed function spaces defined on different meshes and for the adjoint interpolation operator when applied to firedrake.cofunction.Cofunctions.

Parameters:
Returns:

the transferred function

Return type:

firedrake.function.Function or firedrake.cofunction.Cofunction

Extra keyword arguments are passed to firedrake.__future__.interpolate()

project(source, target_space, **kwargs)[source]

Overload function firedrake.projection.project() to account for the case of two mixed function spaces defined on different meshes and for the adjoint projection operator when applied to firedrake.cofunction.Cofunctions.

For details on the approach for achieving boundedness through mass lumping and post-processing, see [FPP+09].

Parameters:
Returns:

the transferred function

Return type:

firedrake.function.Function or firedrake.cofunction.Cofunction

Extra keyword arguments are passed to firedrake.projection.project().

transfer(source, target_space, transfer_method='project', **kwargs)[source]

Overload functions firedrake.__future__.interpolate() and firedrake.projection.project() to account for the case of two mixed function spaces defined on different meshes and for the adjoint interpolation operator when applied to firedrake.cofunction.Cofunctions.

Note that the “interpolate” option works straightforwardly with MPI parallelism, whereas the “project” option can be difficult to set up to make use of this.

Parameters:
Returns:

the transferred function

Return type:

firedrake.function.Function or firedrake.cofunction.Cofunction

Extra keyword arguments are passed to firedrake.__future__.interpolate() or

firedrake.projection.project().

animate.math module

construct_basis(vector, normalise=True)[source]

Construct a basis from a given vector.

Parameters:

vector – the starting vector

Kwargs normalise:

do we want an orthonormal basis?

gram_schmidt(*vectors, normalise=False)[source]

Given some vectors, construct an orthogonal basis using Gram-Schmidt orthogonalisation.

Args vectors:

the vectors to orthogonalise

Kwargs normalise:

do we want an orthonormal basis?

animate.metric module

class RiemannianMetric(*args, **kwargs)[source]

Bases: Function

Class for defining a Riemannian metric over a given mesh.

A metric is a symmetric positive-definite field, which conveys how the mesh is to be adapted. If the mesh is of dimension \(d\) then the metric takes the value of a square \(d\times d\) matrix at each point.

The implementation of metric-based mesh adaptation used in PETSc assumes that the metric is piece-wise linear and continuous, with its degrees of freedom at the mesh vertices.

For details, see the PETSc manual entry:

https://petsc.org/release/docs/manual/dmplex/#metric-based-mesh-adaptation

Parameters:
  • function_space – the tensor FunctionSpace, on which to build this RiemannianMetric. Alternatively, another Function may be passed here and its function space will be used to build this Function. In this case, the function values are copied. If a MeshGeometry is passed here then a tensor \(\mathbb P1\) space is built on top of it.

  • metric_parameters – same as for set_parameters.

assemble_eigendecomposition(evectors, evalues)[source]

Assemble a matrix from its eigenvectors and eigenvalues.

Parameters:
average(*metrics, weights=None)[source]

Average the metric with other metrics.

Args metrics:

the metrics to be averaged with

Parameters:

weights – list of weights to apply to each metric

Returns:

the averaged RiemannianMetric, modified in-place

combine(*metrics, average: bool = True, **kwargs)[source]

Combine metrics using either averaging or intersection.

Parameters:
  • metrics – the list of metrics to combine with

  • average – toggle between averaging and intersection

All other keyword arguments are passed to the relevant method.

complexity(boundary=False)[source]

Compute the metric complexity - the continuous analogue of the (inherently discrete) mesh vertex count.

Parameters:

boundary – should the complexity be computed over the domain boundary?

Returns:

the complexity of the RiemannianMetric

compute_anisotropic_dwr_metric(error_indicator, hessian=None, convergence_rate=1.0, min_eigenvalue=1e-05, interpolant='Clement')[source]

Compute an anisotropic metric from some error indicator, given a Hessian field.

The formulation used is based on that presented in []. Note that normalisation is implicit in the metric construction and involves the convergence_rate parameter, named \(alpha\) in [].

If a Hessian is not provided then an isotropic formulation is used.

Whilst an element-based formulation is used to derive the metric, the result is projected into \(\mathbb P1\) space, by default.

Parameters:
  • error_indicator – the error indicator

  • hessian – the Hessian

  • convergence_rate – normalisation parameter

  • min_eigenvalue – minimum tolerated eigenvalue

  • interpolant – choose from ‘Clement’ or ‘L2’

compute_boundary_hessian(f, method='mixed_L2', **kwargs)[source]

Recover the Hessian of a scalar field on the domain boundary.

Parameters:
  • f – field to recover over the domain boundary

  • method – choose from ‘mixed_L2’ and ‘Clement’

compute_eigendecomposition(reorder=False)[source]

Compute the eigenvectors and eigenvalues of a matrix-valued function.

Parameters:

reorder – should the eigendecomposition be reordered in order of descending eigenvalue magnitude?

Returns:

eigenvector firedrake.function.Function and eigenvalue firedrake.function.Function from the firedrake.functionspace.TensorFunctionSpace() underpinning the metric

compute_hessian(field, method='mixed_L2', **kwargs)[source]

Recover the Hessian of a scalar field.

Parameters:
  • f – the scalar field whose Hessian we seek to recover

  • method – recovery method

All other keyword arguments are passed to the chosen recovery routine.

In the case of the ‘L2’ method, the target_space keyword argument is used for the gradient recovery. The target space for the Hessian recovery is inherited from the metric itself.

compute_isotropic_dwr_metric(error_indicator, convergence_rate=1.0, min_eigenvalue=1e-05, interpolant='Clement')[source]

Compute an isotropic metric from some error indicator using an element-based formulation.

The formulation is based on that presented in []. Note that normalisation is implicit in the metric construction and involves the convergence_rate parameter, named \(alpha\) in [].

Whilst an element-based formulation is used to derive the metric, the result is projected into \(\mathbb P1\) space, by default.

Parameters:
  • error_indicator – the error indicator

  • convergence_rate – normalisation parameter

  • min_eigenvalue – minimum tolerated eigenvalue

  • interpolant – choose from ‘Clement’ or ‘L2’

compute_isotropic_metric(error_indicator, interpolant='Clement', **kwargs)[source]

Compute an isotropic metric from some error indicator.

The result is a \(\mathbb P1\) diagonal tensor field whose entries are projections of the error indicator in modulus.

Parameters:
  • error_indicator – the error indicator

  • interpolant – choose from ‘Clement’ or ‘L2’

compute_weighted_hessian_metric(error_indicators, hessians, average=False, interpolant='Clement')[source]

Compute a vertex-wise anisotropic metric from a list of error indicators, given a list of corresponding Hessian fields.

The formulation used is based on that presented in []. It is assumed that the error indicators have been constructed in the appropriate way.

Parameters:
  • error_indicators – list of error indicators

  • hessians – list of Hessians

  • average – should metric components be averaged or intersected?

  • interpolant – choose from ‘Clement’ or ‘L2’

copy(deepcopy=False)[source]

Copy the metric and any associated parameters.

Parameters:

deepcopy (bool) – If True, the new metric will allocate new space and copy values. If False (default) then the new metric will share the DoF values.

Returns:

a copy of the metric with the same parameters set

Return type:

RiemannianMetric

density_and_quotients(reorder=False)[source]

Extract the density and anisotropy quotients from a metric.

By symmetry, Riemannian metrics admit an orthogonal eigendecomposition,

\[\underline{\mathbf M}(\mathbf x) = \underline{\mathbf V}(\mathbf x)\: \underline{\boldsymbol\Lambda}(\mathbf x)\: \underline{\mathbf V}(\mathbf x)^T,\]

at each point \(\mathbf x\in\Omega\), where \(\underline{\mathbf V}\) and \(\underline{\boldsymbol\Sigma}\) are matrices holding the eigenvectors and eigenvalues, respectively. By positive-definiteness, entries of \(\underline{\boldsymbol\Lambda}\) are all positive.

An alternative decomposition,

\[\underline{\mathbf M}(\mathbf x) = d(\mathbf x)^\frac2n \underline{\mathbf V}(\mathbf x)\: \underline{\mathbf R}(\mathbf x)^{-\frac2n} \underline{\mathbf V}(\mathbf x)^T\]

can also be deduced, in terms of the metric density and anisotropy quotients,

\[d = \prod_{i=1}^n h_i,\qquad r_i = h_i^n d,\qquad \forall i=1:n,\]

where \(h_i := \frac1{\sqrt{\lambda_i}}\).

Parameters:

reorder – should the eigendecomposition be reordered?

Returns:

metric density, anisotropy quotients and eigenvector matrix

enforce_spd(restrict_sizes=False, restrict_anisotropy=False)[source]

Enforce that the metric is symmetric positive-definite.

Parameters:
  • restrict_sizes – should minimum and maximum metric magnitudes be enforced?

  • restrict_anisotropy – should maximum anisotropy be enforced?

Returns:

the RiemannianMetric, modified in-place.

intersect(*metrics)[source]

Intersect the metric with other metrics.

Metric intersection means taking the minimal ellipsoid in the direction of each eigenvector at each point in the domain.

Parameters:

metrics – the metrics to be intersected with

Returns:

the intersected RiemannianMetric, modified in-place

property metric_parameters
normalise(global_factor=None, boundary=False, **kwargs)[source]

Apply \(L^p\) normalisation to the metric.

Parameters:
  • global_factor – pre-computed global normalisation factor

  • boundary – is the normalisation to be done over the boundary?

  • restrict_sizes – should minimum and maximum metric magnitudes be enforced?

  • restrict_anisotropy – should maximum anisotropy be enforced?

Returns:

the normalised RiemannianMetric, modified in-place

set_parameters(metric_parameters=None)[source]

Set metric parameter values internally. All options have the prefix dm_plex_metric_ and are listed below (with prefix dropped for brevity). Note that any parameter which supports Function values must be in a \(\mathbb{P}1\) space defined on the same mesh as the metric, i.e., from FunctionSpace(mesh, "CG", 1).

  • target_complexity: Strictly positive target metric complexity value. No

    default - must be set.

  • h_min: Minimum tolerated metric magnitude, which allows approximate control

    of minimum element size in the adapted mesh. Supports Constant and Function input, as well as float. Default: 1.0e-30.

  • h_max: Maximum tolerated metric magnitude, which allows approximate control

    of maximum element size in the adapted mesh. Supports Constant and Function input, as well as float. Default: 1.0e+30.

  • a_max: Maximum tolerated metric anisotropy, which allows approximate control

    of maximum element anisotropy in the adapted mesh. Supports Constant and Function input, as well as float. Default: 1.0e+05.

  • p: \(L^p\) normalisation order. Supports np.inf as well as

    float values from \([0,\infty)\). Default: 1.0.

  • gradation_factor: Maximum ratio by which adjacent edges in the adapted mesh

    may differ. Default: 1.3. For more detail, see https://www.mmgtools.org/mmg-remesher-try-mmg/mmg-remesher-options/mmg-remesher-option-hgrad.

  • hausdorff_number: Spatial scale factor for the problem. The default value

    0.01 corresponds to an \(\mathcal{O}(1)\) length scale. A rule of thumb is to scale this value appropriately to the length scale of your problem. For more detail, see https://www.mmgtools.org/mmg-remesher-try-mmg/mmg-remesher-options/mmg-remesher-option-hausd.

  • boundary_tag: Mesh boundary tag to restrict attention to during

    boundary-specific metric manipulations. Unset by default, which implies all boundaries are considered. (Note that this parameter does not currently exist in the underlying PETSc implementation.)

  • no_insert: Boolean flag for turning off node insertion and deletion during

    adaptation. Default: False.

  • no_swap: Boolean flag for turning off edge and face swapping during

    adaptation. Default: False.

  • no_move: Boolean flag for turning off node movement during adaptation.

    Default: False.

  • no_surf: Boolean flag for turning off surface modification during adaptation.

    Default: False.

  • num_iterations: Number of adaptation-repartitioning iterations in the

    parallel case. Default: 3. For details on the parallel algorithm, see https://inria.hal.science/hal-02386837.

  • verbosity: Verbosity of the mesh adaptation package (-1 = silent,

    10 = maximum). Default: -1. For more detail, see https://www.mmgtools.org/mmg-remesher-try-mmg/mmg-remesher-options/mmg-remesher-option-v.

  • isotropic: Optimisation for isotropic metrics. (Currently unsupported.)

  • uniform: Optimisation for uniform metrics. (Currently unsupported.)

  • restrict_anisotropy_first: Specify that anisotropy should be restricted

    before normalisation? (Currently unsupported.)

Parameters:

metric_parameters (dict with str keys and value which may take various types) – parameters as above

determine_metric_complexity(H_interior, H_boundary, target, p, **kwargs)[source]

Solve an algebraic problem to obtain coefficients for the interior and boundary metrics to obtain a given metric complexity.

See [] for details. Note that we use a slightly different formulation here.

Parameters:
  • H_interior – Hessian component from domain interior

  • H_boundary – Hessian component from domain boundary

  • target – target metric complexity

  • p – normalisation order

  • H_interior_scaling – optional scaling for interior component

  • H_boundary_scaling – optional scaling for boundary component

intersect_on_boundary(*metrics, boundary_tag='on_boundary')[source]

Combine a list of metrics by intersection.

Parameters:
  • metrics – the metrics to be combined

  • boundary_tag – optional boundary segment physical ID for boundary intersection. Otherwise, the intersection is over the whole boundary.

animate.quality module

Functions for computing mesh quality measures.

class QualityMeasure(mesh, metric=None, python=False)[source]

Bases: object

Class for computing quality measures associated with a given mesh.

Choices of quality measure:
  • min_angle: the minimum angle of each cell

  • area: the area of each cell in a 2D triangular mesh

  • volume: the volume of each cell in a 3D tetrahedral mesh

  • facet_area: the area of each facet.

  • aspect_ratio: the aspect ratio of each cell

  • eskew: the equiangle skew of each cell

  • skewness: the skewness of each cell in a 2D triangular mesh

  • scaled_jacobian: the scaled Jacobian of each cell

  • metric: given a Riemannian metric, this function outputs the value of the quality measure Q_M based on the transformation encoded by the metric.

Parameters:
  • mesh – the input mesh to do computations on

  • metric – the tensor field representing the metric space transformation

  • python – compute the measure using Python?

animate.recovery module

Driver functions for derivative recovery.

recover_boundary_hessian(f, method='Clement', target_space=None, **kwargs)[source]

Recover the Hessian of a scalar field on the domain boundary.

Parameters:
recover_gradient_l2(f, target_space=None)[source]

Recover the gradient of a scalar or vector field using \(L^2\) projection.

Parameters:
recover_hessian_clement(f)[source]

Recover the gradient and Hessian of a scalar field using two applications of Clement interpolation.

Note that if the field is of degree 2 then projection will be used to obtain the gradient. If the field is of degree 3 or greater then projection will be used for the Hessian recovery, too.

Parameters:

f – the scalar field whose derivatives we seek to recover

animate.utility module

Utility functions and classes for metric-based mesh adaptation.

Mesh(arg, **kwargs)[source]

Overload firedrake.mesh.Mesh() to endow the output mesh with useful quantities.

The following quantities are computed by default:
  • cell size;

  • facet area.

The argument and keyword arguments are passed to Firedrake’s firedrake.mesh.Mesh() constructor, modified so that the argument could also be a mesh.

Arguments and keyword arguments are the same as for firedrake.mesh.Mesh().

Returns:

the constructed mesh

Return type:

firedrake.mesh.MeshGeometry

class VTKFile(*args, **kwargs)[source]

Bases: VTKFile

Overload firedrake.output.VTKFile so that it uses adaptive mode by default.

Whilst this means that the mesh topology is recomputed at every export, it removes any need for the user to reset it manually.

Create an object for outputting data for visualisation.

This produces output in VTU format, suitable for visualisation with Paraview or other VTK-capable visualisation packages.

Parameters:
  • filename – The name of the output file (must end in .pvd).

  • project_output – Should the output be projected to a computed output space? Default is to use interpolation.

  • comm – The MPI communicator to use.

  • mode – “w” to overwrite any existing file, “a” to append to an existing file.

  • target_degree – override the degree of the output space.

  • target_continuity – override the continuity of the output space; A UFL ufl.sobolevspace.SobolevSpace object: H1 for a continuous output and L2 for a discontinuous output.

  • adaptive – allow different meshes at different exports if True.

Note

Visualisation is only possible for Lagrange fields (either continuous or discontinuous). All other fields are first either projected or interpolated to Lagrange elements before storing for visualisation purposes.

errornorm(u, uh, norm_type='L2', boundary=False, **kwargs)[source]

Overload firedrake.norms.errornorm() to allow for \(\ell^p\) norms.

Currently supported norm_type options: * 'l1' * 'l2' * 'linf' * 'L2' * 'Linf' * 'H1' * 'Hdiv' * 'Hcurl' * or any 'Lp' with \(p >= 1\).

Note that this version is case sensitive, i.e. 'l2' and 'L2' will give different results in general.

Parameters:
Returns:

the error norm value

Return type:

float

Any other keyword arguments are passed to firedrake.norms.errornorm().

norm(v, norm_type='L2', condition=None, boundary=False)[source]

Overload firedrake.norms.norm() to allow for \(\ell^p\) norms.

Currently supported norm_type options: * 'l1' * 'l2' * 'linf' * 'L2' * 'Linf' * 'H1' * 'Hdiv' * 'Hcurl' * or any 'Lp' with \(p >= 1\).

Note that this version is case sensitive, i.e. 'l2' and 'L2' will give different results in general.

Parameters:
Returns:

the norm value

Return type:

float

Module contents