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:
mesh (
firedrake.mesh.MeshGeometry
) – the physical mesh to be movedtimestep (
float
) – the timestep length used
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:
mesh (
firedrake.mesh.MeshGeometry
) – the physical meshmonitor_function (
Callable
) – a Python function which takes a mesh as inputmethod (
str
) – choose from ‘relaxation’ and ‘quasi_newton’phi_init (
firedrake.function.Function
) – initial guess for the scalar potentialsigma_init (
firedrake.function.Function
) – initial guess for the Hessian
- Returns:
the Monge-Ampere Mover object
- Return type:
- 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 meshmonitor_function (
Callable
) – a Python function which takes a mesh as inputphi_init (
firedrake.function.Function
) – initial guess for the scalar potentialsigma_init (
firedrake.function.Function
) – initial guess for the Hessianmaxiter (
int
) – maximum number of iterations for the Quasi-Newton solverrtol (
float
) – relative tolerance for the residualdtol (
float
) – divergence tolerance for the residual
- property equidistributor¶
Setup the equidistributor for the quasi-newton method.
- Returns:
the equidistributor
- Return type:
NonlinearVariationalSolver
- 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 meshmonitor_function (
Callable
) – a Python function which takes a mesh as inputphi_init (
firedrake.function.Function
) – initial guess for the scalar potentialsigma_init (
firedrake.function.Function
) – initial guess for the Hessianpseudo_timestep (
float
) – pseudo-timestep to use for the relaxationmaxiter (
int
) – maximum number of iterations for the relaxationrtol (
float
) – relative tolerance for the residualdtol (
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:
- 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 meshmonitor_function (
Callable
) – a Python function which takes a mesh as inputraise_convergence_errors – convergence error handling behaviour: if True then
ConvergenceError
s are raised, else warnings are raised and the program is allowed to continueraise_convergence_errors –
bool
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:
mesh (
firedrake.mesh.MeshGeometry
) – the physical mesh to be movedtimestep (
float
) – the timestep length usedmethod (
str
) – flavour of spring-based method to use
- 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:
mesh (
firedrake.mesh.MeshGeometry
) – the physical mesh to be movedtimestep (
float
) – the timestep length used
- 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:
mesh (
firedrake.mesh.MeshGeometry
) – the physical mesh to be movedtimestep (
float
) – the timestep length used
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 tangledraise_error (
bool
) – ifTrue
, 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¶