movement package

Submodules

movement.laplacian_smoothing module

class LaplacianSmoother(mesh, timestep, **kwargs)[source]

Bases: PrimeMover

Movement of a mesh is driven by a mesh velocity \(\mathbf{v}\), which is determined by solving a vector Laplace problem

\[\nabla^2_{\boldsymbol{\xi}}\mathbf{v} = \boldsymbol{0}, \quad \boldsymbol{\xi}\in\Omega,\]

under non-zero Dirichlet boundary conditions on a forced boundary section \(\partial\Omega_f\) and zero Dirichlet boundary conditions elsewhere:

\[\begin{split}\mathbf{v} = \left\{\begin{array}{rl} \mathbf{v}_D, & \boldsymbol{\xi}\in\partial\Omega_f\\ \boldsymbol{0}, & \boldsymbol{\xi}\in\partial\Omega\backslash\partial\Omega_f \end{array}\right.\end{split}\]

where the computational coordinates \(\boldsymbol{\xi} := \mathbf{x}(t_0)\) are the physical coordinates at the beginning of the simulation.

Parameters:
move(time, update_boundary_velocity=None, boundary_conditions=None)[source]

Assemble and solve the Laplacian system and update the coordinates.

Parameters:
  • time (float) – the current time

  • update_boundary_velocity (Callable with a single argument of float type) – function that updates the boundary conditions at the current time

  • boundary_conditions (DirichletBC or list thereof) – Dirichlet boundary conditions to be enforced

movement.monge_ampere module

MongeAmpereMover(mesh, monitor_function, method='relaxation', **kwargs)[source]

Movement of a mesh is determined by a monitor_function \(m\) and the Monge-Ampère type equation

\[m(x)\det(I + H(\phi)) = \theta,\]

for a scalar potential \(\phi\), where \(I\) is the identity matrix, \(\theta\) is a normalisation coefficient and \(H(\phi)\) denotes the Hessian of \(\phi\) with respect to the coordinates \(\xi\) of the computational mesh.

The physical mesh coordinates \(x\) are updated according to

\[x = \xi + \nabla\phi.\]
Parameters:
Returns:

the Monge-Ampere Mover object

Return type:

MongeAmpereMover_Relaxation or MongeAmpereMover_QuasiNewton

class MongeAmpereMover_QuasiNewton(mesh, monitor_function, phi_init=None, sigma_init=None, **kwargs)[source]

Bases: MongeAmpereMover_Base

The elliptic Monge-Ampere equation is solved using a quasi-Newton method (see [MCB18] for details).

Parameters:
  • mesh (firedrake.mesh.MeshGeometry) – the physical mesh

  • monitor_function (Callable) – a Python function which takes a mesh as input

  • phi_init (firedrake.function.Function) – initial guess for the scalar potential

  • sigma_init (firedrake.function.Function) – initial guess for the Hessian

  • maxiter (int) – maximum number of iterations for the Quasi-Newton solver

  • rtol (float) – relative tolerance for the residual

  • dtol (float) – divergence tolerance for the residual

property equidistributor

Setup the equidistributor for the quasi-newton method.

Returns:

the equidistributor

Return type:

NonlinearVariationalSolver

move()[source]

Run the quasi-Newton method to convergence and update the mesh.

Returns:

the iteration count

Return type:

int

class MongeAmpereMover_Relaxation(mesh, monitor_function, phi_init=None, sigma_init=None, **kwargs)[source]

Bases: MongeAmpereMover_Base

The elliptic Monge-Ampere equation is solved in a parabolised form using a pseudo-time relaxation,

\[-\frac\partial{\partial\tau}\Delta\phi = m(x)\det(I + H(\phi)) - \theta,\]

where \(\tau\) is the pseudo-time variable. Forward Euler is used for the pseudo-time integration (see [MCB18] for details).

Parameters:
  • mesh (firedrake.mesh.MeshGeometry) – the physical mesh

  • monitor_function (Callable) – a Python function which takes a mesh as input

  • phi_init (firedrake.function.Function) – initial guess for the scalar potential

  • sigma_init (firedrake.function.Function) – initial guess for the Hessian

  • pseudo_timestep (float) – pseudo-timestep to use for the relaxation

  • maxiter (int) – maximum number of iterations for the relaxation

  • rtol (float) – relative tolerance for the residual

  • dtol (float) – divergence tolerance for the residual

property equidistributor

Setup the equidistributor for the relaxation method.

Returns:

the equidistributor

Return type:

LinearVariationalSolver

move()[source]

Run the relaxation method to convergence and update the mesh.

Returns:

the iteration count

Return type:

int

property pseudotimestepper

Setup the pseudo-timestepper for the relaxation method.

Returns:

the pseudo-timestepper

Return type:

LinearVariationalSolver

movement.mover module

class PrimeMover(mesh, monitor_function=None, raise_convergence_errors=True, **kwargs)[source]

Bases: object

Base class for all mesh movers.

Parameters:
  • mesh (firedrake.mesh.MeshGeometry) – the physical mesh

  • monitor_function (Callable) – a Python function which takes a mesh as input

  • raise_convergence_errors – convergence error handling behaviour: if True then ConvergenceErrors are raised, else warnings are raised and the program is allowed to continue

  • raise_convergence_errorsbool

adapt()[source]

Alias of move.

coordinate(index)[source]

Get the mesh coordinate associated with a given index.

coordinate_offset(index)[source]

Get the DMPlex coordinate section offset for a given index.

edge_vector_offset(index)[source]

Get the DMPlex edge vector section offset for a given index.

move()[source]

Move the mesh according to the method of choice.

movement.solver_parameters module

movement.spring module

SpringMover(*args, method='lineal', **kwargs)[source]

Movement of a mesh is determined by reinterpreting it as a structure of stiff beams and solving an associated discrete linear elasticity problem. (See [] for details.)

Parameters:
class SpringMover_Lineal(mesh, timestep, **kwargs)[source]

Bases: SpringMover_Base

Movement of a mesh is determined by reinterpreting it as a structure of stiff beams and solving an associated discrete linear elasticity problem.

We consider the ‘lineal’ case, as described in [].

Parameters:
move(time, update_boundary_displacement=None, boundary_conditions=None)[source]

Assemble and solve the lineal spring system and update the coordinates.

Parameters:
  • time (float) – the current time

  • update_boundary_displacement (Callable with a single argument of float type) – function that updates the boundary conditions at the current time

  • boundary_conditions (DirichletBC or list thereof) – Dirichlet boundary conditions to be enforced

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

Bases: SpringMover_Lineal

Movement of a mesh is determined by reinterpreting it as a structure of stiff beams and solving an associated discrete linear elasticity problem.

We consider the ‘torsional’ case, as described in [].

Parameters:

movement.tangling module

class MeshTanglingChecker(mesh, raise_error=True)[source]

Bases: object

A class for tracking whether a mesh has tangled, i.e. whether any of its elements have become inverted.

Parameters:
  • mesh (firedrake.mesh.MeshGeometry) – the mesh to track if tangled

  • raise_error (bool) – if True, an error is raised if any element is tangled, otherwise a warning is raised

check()[source]

Check whether any element orientations have changed since the tangling checker was created.

property scaled_jacobian

Module contents