underworld.systems module¶
This module contains routines relating to differential system.
Module Summary¶
submodules:¶
functions:¶
underworld.systems.Solver |
This method simply returns a necessary solver for the provided system. |
classes:¶
underworld.systems.Stokes |
This class provides functionality for a discrete representation of the Stokes flow equations. |
underworld.systems.SteadyStateHeat |
This class provides functionality for a discrete representation of the steady state heat equation. |
underworld.systems.SteadyStateDarcyFlow |
This class provides functionality for a discrete representation of the steady state darcy flow equation. |
underworld.systems.TimeIntegration |
Abstract class for integrating numerical objects (fields, swarms, etc.) in time. |
underworld.systems.AdvectionDiffusion |
This class provides functionality for a discrete representation of an advection-diffusion equation. |
underworld.systems.SwarmAdvector |
Objects of this class advect a swarm through time using the provided velocity field. |
Module Details¶
functions:¶
classes:¶
-
class
underworld.systems.
Stokes
(velocityField, pressureField, fn_viscosity, fn_bodyforce=None, fn_one_on_lambda=None, fn_lambda=None, fn_source=None, voronoi_swarm=None, conditions=[], _removeBCs=True, _fn_viscosity2=None, _fn_director=None, fn_stresshistory=None, _fn_stresshistory=None, _fn_v0=None, _fn_p0=None, _callback_post_solve=None, **kwargs)[source]¶ Bases:
underworld._stgermain.StgCompoundComponent
This class provides functionality for a discrete representation of the Stokes flow equations.
Specifically, the class uses a mixed finite element method to construct a system of linear equations which may then be solved using an object of the underworld.system.Solver class.
The underlying element types are determined by the supporting mesh used for the ‘velocityField’ and ‘pressureField’ parameters.
The strong form of the given boundary value problem, for \(f\), \(g\) and \(h\) given, is
\[\begin{split}\begin{align} \sigma_{ij,j} + f_i =& \: 0 & \text{ in } \Omega \\ u_{k,k} + \frac{p}{\lambda} =& \: H & \text{ in } \Omega \\ u_i =& \: g_i & \text{ on } \Gamma_{g_i} \\ \sigma_{ij}n_j =& \: h_i & \text{ on } \Gamma_{h_i} \\ \end{align}\end{split}\]where,
- \(\sigma_{i,j}\) is the stress tensor
- \(u_i\) is the velocity,
- \(p\) is the pressure,
- \(f_i\) is a body force,
- \(\lambda\) is a bulk viscosity,
- \(H\) is the compressible equation source term,
- \(g_i\) are the velocity boundary conditions (DirichletCondition)
- \(h_i\) are the traction boundary conditions (NeumannCondition).
The problem boundary, \(\Gamma\), admits the decompositions \(\Gamma=\Gamma_{g_i}\cup\Gamma_{h_i}\) where \(\emptyset=\Gamma_{g_i}\cap\Gamma_{h_i}\). The equivalent weak form is:
\[\int_{\Omega} w_{(i,j)} \sigma_{ij} \, d \Omega = \int_{\Omega} w_i \, f_i \, d\Omega + \sum_{j=1}^{n_{sd}} \int_{\Gamma_{h_j}} w_i \, h_i \, d \Gamma\]where we must find \(u\) which satisfies the above for all \(w\) in some variational space.
Parameters: - velocityField (underworld.mesh.MeshVariable) – Variable used to record system velocity.
- pressureField (underworld.mesh.MeshVariable) – Variable used to record system pressure.
- fn_viscosity (underworld.function.Function) – Function which reports a viscosity value. Function must return scalar float values.
- fn_bodyforce (underworld.function.Function, Default = None) – Function which reports a body force for the system. Function must return float values of identical dimensionality to the provided velocity variable.
- fn_minus_one_on_lambda (underworld.function.Function, Default = None) – Function which defines a non solenoidal velocity field via the relationship div(velocityField) = -fn_minus_one_on_lambda * pressurefield + fn_source When this is left as None a incompressible formulation of the stokes equation is formed, ie, div(velocityField) = 0. fn_minus_one_on_lambda is incompatible with the ‘penalty’ stokes solver, ensure a ‘penalty’ equal to 0 is used when fn_minus_one_on_lambda is used. By default this is the case.
- fn_source (underworld.function.Function, Default = None) – Function which defines a non solenoidal velocity field via the relationship div(velocityField) = -fn_minus_one_on_lambda * pressurefield + fn_source. fn_minus_one_on_lambda is incompatible with the ‘penalty’ stokes solver, ensure the ‘penalty’ of 0, is used when fn_lambda is used. By default this is the case.
- fn_stresshistory (underworld.function.Function, Default = None) – Function which defines the stress history term used for viscoelasticity. Function is a vector of size 3 (dim=2) or 6 (dim=3) representing a symetric tensor.
- voronoi_swarm (underworld.swarm.Swarm) – If a voronoi_swarm is provided, voronoi type numerical integration is utilised. The provided swarm is used as the basis for the voronoi integration. If no voronoi_swarm is provided, Gauss integration is used.
- conditions (underworld.conditions.SystemCondition) – Numerical conditions to impose on the system. This should be supplied as the condition itself, or a list object containing the conditions.
Notes
Constructor must be called by collectively all processes.
-
eqResiduals
¶ Returns the stokes flow equations’ residuals from the latest solve. Residual calculations use the matrices and vectors of the discretised problem. The residuals correspond to the momentum equation and the continuity equation.
Returns: r1 is the momentum equation residual r2 is the continuity equation residual Return type: (r1, r2) - 2 tuple of doubles Notes
This method must be called collectively by all processes.
-
fn_bodyforce
¶ The body force function. You may change this function directly via this property.
-
fn_one_on_lambda
¶ A bulk viscosity parameter
-
fn_source
¶ The volumetric source term function. You may change this function directly via this property.
-
fn_viscosity
¶ The viscosity function. You may change this function directly via this property.
-
stokes_callback
¶ Return the callback function used by this system
-
class
underworld.systems.
SteadyStateHeat
(temperatureField, fn_diffusivity, fn_heating=0.0, voronoi_swarm=None, conditions=[], _removeBCs=True, **kwargs)[source]¶ Bases:
underworld._stgermain.StgCompoundComponent
This class provides functionality for a discrete representation of the steady state heat equation.
The class uses a standard Galerkin finite element method to construct a system of linear equations which may then be solved using an object of the underworld.system.Solver class.
The underlying element types are determined by the supporting mesh used for the ‘temperatureField’.
The strong form of the given boundary value problem, for \(f\), \(h\) and \(h\) given, is
\[\begin{split}\begin{align} q_i =& - \alpha \, u_{,i} & \\ q_{i,i} =& \: f & \text{ in } \Omega \\ u =& \: g & \text{ on } \Gamma_g \\ -q_i n_i =& \: h & \text{ on } \Gamma_h \\ \end{align}\end{split}\]where, \(\alpha\) is the diffusivity, \(u\) is the temperature, \(f\) is a source term, \(g\) is the Dirichlet condition, and \(h\) is a Neumann condition. The problem boundary, \(\Gamma\), admits the decomposition \(\Gamma=\Gamma_g\cup\Gamma_h\) where \(\emptyset=\Gamma_g\cap\Gamma_h\). The equivalent weak form is:
\[-\int_{\Omega} w_{,i} \, q_i \, d \Omega = \int_{\Omega} w \, f \, d\Omega + \int_{\Gamma_h} w \, h \, d \Gamma\]where we must find \(u\) which satisfies the above for all \(w\) in some variational space.
Parameters: - temperatureField (underworld.mesh.MeshVariable) – The solution field for temperature.
- fn_diffusivity (underworld.function.Function) – The function that defines the diffusivity across the domain.
- fn_heating (underworld.function.Function) – A function that defines the heating across the domain. Optional.
- voronoi_swarm (underworld.swarm.Swarm) – If a voronoi_swarm is provided, voronoi type numerical integration is utilised. The provided swarm is used as the basis for the voronoi integration. If no voronoi_swarm is provided, Gauss integration is used.
- conditions (underworld.conditions.SystemCondition) – Numerical conditions to impose on the system. This should be supplied as the condition itself, or a list object containing the conditions.
Notes
Constructor must be called collectively by all processes.
Example
Setup a basic thermal system:
>>> linearMesh = uw.mesh.FeMesh_Cartesian( elementType='Q1/dQ0', elementRes=(4,4), minCoord=(0.,0.), maxCoord=(1.,1.) ) >>> tField = uw.mesh.MeshVariable( linearMesh, 1 ) >>> topNodes = linearMesh.specialSets["MaxJ_VertexSet"] >>> bottomNodes = linearMesh.specialSets["MinJ_VertexSet"] >>> tbcs = uw.conditions.DirichletCondition(tField, topNodes + bottomNodes) >>> tField.data[topNodes.data] = 0.0 >>> tField.data[bottomNodes.data] = 1.0 >>> tSystem = uw.systems.SteadyStateHeat(temperatureField=tField, fn_diffusivity=1.0, conditions=[tbcs])
Example with non diffusivity:
>>> k = tField + 1.0 >>> tSystem = uw.systems.SteadyStateHeat(temperatureField=tField, fn_diffusivity=k, conditions=[tbcs]) >>> solver = uw.systems.Solver(tSystem) >>> solver.solve() Traceback (most recent call last): ... RuntimeError: Nonlinearity detected. Diffusivity function depends on the temperature field provided to the system. Please set the 'nonLinearIterate' solve parameter to 'True' or 'False' to continue. >>> solver.solve(nonLinearIterate=True)
-
fn_diffusivity
¶ The diffusivity function. You may change this function directly via this property.
-
fn_heating
¶ The heating function. You may change this function directly via this property.
-
class
underworld.systems.
SteadyStateDarcyFlow
(pressureField, fn_diffusivity, fn_bodyforce=None, voronoi_swarm=None, conditions=[], velocityField=None, swarmVarVelocity=None, _removeBCs=True, **kwargs)[source]¶ Bases:
underworld._stgermain.StgCompoundComponent
This class provides functionality for a discrete representation of the steady state darcy flow equation.
The class uses a standard Galerkin finite element method to construct a system of linear equations which may then be solved using an object of the underworld.system.Solver class.
The underlying element types are determined by the supporting mesh used for the ‘pressureField’.
The strong form of the given boundary value problem, for \(f\), \(q\) and \(h\) given, is
\[\begin{split}\begin{align} q_i =& \kappa \, ( -u_{,i} + S_i ) & \\ q_{i,i} =& \: f & \text{ in } \Omega \\ u =& \: q & \text{ on } \Gamma_q \\ -q_i n_i =& \: h & \text{ on } \Gamma_h \\ \end{align}\end{split}\]where,
- \(\kappa\) is the diffusivity, \(u\) is the pressure,
- \(S\) is a flow body-source, for example due to gravity,
- \(f\) is a source term, \(q\) is the Dirichlet condition, and
- \(h\) is a Neumann condition.
The problem boundary, \(\Gamma\), admits the decomposition \(\Gamma=\Gamma_q\cup\Gamma_h\) where \(\emptyset=\Gamma_q\cap\Gamma_h\). The equivalent weak form is:
\[-\int_{\Omega} w_{,i} \, q_i \, d \Omega = \int_{\Omega} w \, f \, d\Omega + \int_{\Gamma_h} w \, h \, d \Gamma\]where we must find \(u\) which satisfies the above for all \(w\) in some variational space.
Parameters: - pressureField (underworld.mesh.MeshVariable) – The solution field for pressure.
- fn_diffusivity (underworld.function.Function) – The function that defines the diffusivity across the domain.
- fn_bodyforce (underworld.function.Function) – A function that defines the flow body-force across the domain, for example gravity. Must be a vector. Optional.
- voronoi_swarm (underworld.swarm.Swarm) – A swarm with just one particle within each cell should be provided. This avoids the evaluation of the velocity on nodes and inaccuracies arising from diffusivity changes within cells. If a swarm is provided, voronoi type numerical integration is utilised. The provided swarm is used as the basis for the voronoi integration. If no voronoi_swarm is provided, Gauss integration is used.
- conditions (underworld.conditions.SystemCondition) – Numerical conditions to impose on the system. This should be supplied as the condition itself, or a list object containing the conditions.
- velocityField (underworld.mesh.MeshVariable) – Solution field for darcy flow velocity. Optional.
- swarmVarVelocity (undeworld.swarm.SwarmVariable) – If a swarm variable is provided, the velocity calculated on the swarm will be stored. This is the most representative velocity data object, as the velocity calculation occurs on the swarm, away from mesh nodes. Optional.
Notes
Constructor must be called collectively by all processes.
-
fn_bodyforce
¶ The heating function. You may change this function directly via this property.
-
fn_diffusivity
¶ The diffusivity function. You may change this function directly via this property.
-
class
underworld.systems.
TimeIntegration
(order, **kwargs)[source]¶ Bases:
underworld._stgermain.StgCompoundComponent
Abstract class for integrating numerical objects (fields, swarms, etc.) in time.
The integration algorithm is a modified Runge Kutta method that only evaluates midpoint information varying in space - using only the present timestep solution. The order of the integration used can be 1,2,4
Parameters: order (int {1,2,4}) – Defines the numerical order ‘in space’ of the Runge Kutta like integration scheme. -
dt
¶ Time integrator timestep size.
-
time
¶ Time integrator time value.
-
-
class
underworld.systems.
AdvectionDiffusion
(phiField, phiDotField, velocityField, fn_diffusivity, fn_sourceTerm=None, conditions=[], _allow_non_q1=False, **kwargs)[source]¶ Bases:
underworld._stgermain.StgCompoundComponent
This class provides functionality for a discrete representation of an advection-diffusion equation.
The class uses the Streamline Upwind Petrov Galerkin SUPG method to integrate through time.
\[\frac{\partial\phi}{\partial t} + {\bf u } \cdot \nabla \phi= \nabla { ( k \nabla \phi ) } + H\]Parameters: - phiField (underworld.mesh.MeshVariable) – The concentration field, typically the temperature field
- phiDotField (underworld.mesh.MeshVariable) – A MeshVariable that defines the initial time derivative of the phiField. Typically 0 at the beginning of a model, e.g. phiDotField.data[:]=0 When using a phiField loaded from disk one should also load the phiDotField to ensure the solving method has the time derivative information for a smooth restart. No dirichlet conditions are required for this field as the phiField degrees of freedom map exactly to this field’s dirichlet conditions, the value of which ought to be 0 for constant values of phi.
- velocityField (underworld.mesh.MeshVariable) – The velocity field.
- fn_diffusivity (underworld.function.Function) – A function that defines the diffusivity within the domain.
- fn_sourceTerm (underworld.function.Function) – A function that defines the heating within the domain. Optional.
- conditions (underworld.conditions.SystemCondition) – Numerical conditions to impose on the system. This should be supplied as the condition itself, or a list object containing the conditions.
Notes
Constructor must be called by collectively all processes.
-
class
underworld.systems.
SwarmAdvector
(velocityField, swarm, order=2, **kwargs)[source]¶ Bases:
underworld.systems._timeintegration.TimeIntegration
Objects of this class advect a swarm through time using the provided velocity field.
Parameters: - velocityField (underworld.mesh.MeshVariable) – The MeshVariable field used for evaluating the velocity field that advects the swarm particles
- swarm (underworld.swarm.Swarm) – Particle swarm that will be advected by the given velocity field
-
integrate
(dt, update_owners=True)[source]¶ Integrate the associated swarm in time, by dt, using the velocityfield that is associated with this class
Parameters: - dt (double) – The timestep to use in the intergration
- update_owners (bool) – If this is set to False, particle ownership (which element owns a particular particle) is not updated after advection. This is often necessary when both the mesh and particles are advecting simutaneously.
Example
>>> import underworld as uw >>> import numpy as np >>> import numpy.testing as npt >>> from underworld import function as fn >>> dim=2; >>> elementMesh = uw.mesh.FeMesh_Cartesian(elementType="Q1/dQ0", elementRes=(9,9), minCoord=(-1.,-1.), maxCoord=(1.,1.)) >>> velocityField = uw.mesh.MeshVariable( mesh=elementMesh, nodeDofCount=dim ) >>> swarm = uw.swarm.Swarm(mesh=elementMesh) >>> particle = np.zeros((1,2)) >>> particle[0] = [0.2,-0.2] >>> swarm.add_particles_with_coordinates(particle) array([0], dtype=int32) >>> velocityField.data[:]=[1.0,1.0] >>> swarmAdvector = uw.systems.SwarmAdvector(velocityField=velocityField, swarm=swarm, order=2) >>> dt=1.0 >>> swarmAdvector.integrate(dt) >>> npt.assert_allclose(swarm.particleCoordinates.data[0], [1.2,0.8], rtol=1e-8, verbose=True)