underworld.utils module

Various utility classes & functions.

Functions

underworld.utils.is_kernel Function to determine if the script is being run in an ipython or jupyter notebook or in a regular python interpreter.
underworld.utils.matplotlib_inline This function simply enables Jupyter Notebook inlined matplotlib results.

Classes

underworld.utils.Integral The Integral class constructs the volume integral
underworld.utils.MeshVariable_Projection This class provides functionality for projecting data from any underworld function onto a provided mesh variable.
underworld.utils.SavedFileData A class used to define saved data.
underworld.utils.SolveLinearSystem To solve for x in the equation Ax=b
underworld.utils.is_kernel()[source]

Function to determine if the script is being run in an ipython or jupyter notebook or in a regular python interpreter.

Return true if in ipython or Jupyter notebook, False otherwise.

underworld.utils.matplotlib_inline()[source]

This function simply enables Jupyter Notebook inlined matplotlib results. This function should be called at the start of your notebooks as a replacement for the Jupyter Notebook %matplotlib inline magic. It provides the same functionality, however it allows notebooks to be converted to python without having to explicitly remove these calls.

class underworld.utils.Integral(fn, mesh=None, integrationType='volume', surfaceIndexSet=None, integrationSwarm=None, **kwargs)[source]

Bases: underworld._stgermain.StgCompoundComponent

The Integral class constructs the volume integral

\[F_{i} = \int_V \, f_i(\mathbf{x}) \, \mathrm{d} V\]

for some function \(f_i\) (specified by a Function object), over some domain \(V\) (specified by an FeMesh object), or the surface integral

\[F_{i} = \oint_{\Gamma} \, f_i(\mathbf{x}) \, \mathrm{d}\Gamma\]

for some surface \(\Gamma\) (specified via an IndexSet object on the mesh).

Parameters:
  • fn (uw.function.Function) – Function to be integrated.
  • mesh (uw.mesh.FeMesh) – The mesh over which integration is performed.
  • integrationType (str) – Type of integration to perform. Options are “volume” or “surface”.
  • surfaceIndexSet (uw.mesh.FeMesh_IndexSet) – Must be provided where integrationType is “surface”. This IndexSet determines which surface is to be integrated over. Note that surface integration over interior nodes is not currently supported.
  • integrationSwarm (uw.swarm.IntegrationSwarm (optional)) – User provided integration swarm.

Notes

Constructor must be called by collectively all processes.

Example

Calculate volume of mesh:

>>> import underworld as uw
>>> mesh = uw.mesh.FeMesh_Cartesian(minCoord=(0.,0.), maxCoord=(1.,1.))
>>> volumeIntegral = uw.utils.Integral(fn=1.,mesh=mesh)
>>> np.allclose( 1., volumeIntegral.evaluate(), rtol=1e-8)
True

Calculate surface area of mesh:

>>> surfaceIntegral = uw.utils.Integral(fn=1., mesh=mesh, integrationType='surface', surfaceIndexSet=mesh.specialSets["AllWalls_VertexSet"])
>>> np.allclose( 4., surfaceIntegral.evaluate(), rtol=1e-8)
True
evaluate()[source]

Perform integration.

Notes

Method must be called collectively by all processes.

Returns:result – Integration result. For vector integrals, a vector is returned.
Return type:list of floats
maskFn

The integration mask used where surface integration is performed.

class underworld.utils.MeshVariable_Projection(meshVariable=None, fn=None, voronoi_swarm=None, type=0, gauss_swarm=None, **kwargs)[source]

Bases: underworld._stgermain.StgCompoundComponent

This class provides functionality for projecting data from any underworld function onto a provided mesh variable.

For the variable \(\bf{U} = \bf{u}_a N_a\) and function \(F\), we wish to determine appropriate values for \(\bf{u}_a\) such that \(\bf{U} \simeq F\).

Two projection methods are supported; weighted averages and weighted residuals. Generally speaking, weighted averages provide robust low order results, while weighted residuals give higher accuracy but spurious results for difficult functions \(F\).

The weighted average method is defined as:

\[\bf{u}_a = \frac{\int_{\Omega} \bf{F} N_a \partial\Omega }{\int_{\Omega} N_a \partial\Omega }\]

The weighted residual method constructs an SLE which is then solved to determine the unknowns:

\[\bf{u}_a\int_{\Omega} N_a N_b \partial\Omega = \int_{\Omega} \bf{F} N_b \partial\Omega\]
Parameters:
  • meshVariable (underworld.mesh.MeshVariable) – The variable you wish to project the function onto.
  • fn (underworld.function.Function) – The function you wish to project.
  • voronoi_swarm (underworld.swarm.Swarm) – Optional. If a voronoi_swarm is provided, voronoi type integration is utilised to integrate across elements. The provided swarm is used as the basis for the voronoi integration. If no swarm is provided, Gauss integration is used.
  • gauss_swarm (underworld.swarm.GaussIntegrationSwarm) – If provided this gauss_swarm will be used for (gaussian) numerical integration rather than a default gauss integration swarm that is automatically generated and dependent on the element order of the mesh. NB: if a voronoi_swarm is defined it OVERRIDES this gauss_swarm as the preferred integration swarm (quadrature method).
  • type (int, default=0) – Projection type. 0:weighted average, 1:weighted residual

Notes

Constructor must be called collectively by all processes.

Examples

>>> import underworld as uw
>>> import numpy as np
>>> mesh = uw.mesh.FeMesh_Cartesian()
>>> U = uw.mesh.MeshVariable( mesh, 1 )

Lets cast a constant value onto this mesh variable

>>> const = 1.23456
>>> projector = uw.utils.MeshVariable_Projection( U, const, type=0 )
>>> np.allclose(U.data, const)
False
>>> projector.solve()
>>> np.allclose(U.data, const)
True

Now cast mesh coordinates onto a vector variable

>>> U_coord = uw.mesh.MeshVariable( mesh, 2 )
>>> projector = uw.utils.MeshVariable_Projection( U_coord, uw.function.coord(), type=1 )
>>> projector.solve()
>>> np.allclose(U_coord.data, mesh.data)
True

Project one mesh variable onto another

>>> U_copy = uw.mesh.MeshVariable( mesh, 2 )
>>> projector = uw.utils.MeshVariable_Projection( U_copy, U_coord, type=1 )
>>> projector.solve()
>>> np.allclose(U_copy.data, U_coord.data)
True

Project the coords to the submesh (usually the constant mesh)

>>> U_submesh = uw.mesh.MeshVariable( mesh.subMesh, 2 )
>>> projector = uw.utils.MeshVariable_Projection( U_submesh, U_coord, type=1 )
>>> projector.solve()
>>> np.allclose(U_submesh.data, mesh.subMesh.data)
True

Create swarm, then project particle owning elements onto mesh

>>> U_submesh = uw.mesh.MeshVariable( mesh.subMesh, 1 )
>>> swarm = uw.swarm.Swarm(mesh)
>>> swarm.populate_using_layout(uw.swarm.layouts.PerCellSpaceFillerLayout(swarm,4))
>>> projector = uw.utils.MeshVariable_Projection( U_submesh, swarm.owningCell, type=1 )
>>> projector.solve()
>>> np.allclose(U_submesh.data[:,0], mesh.data_elgId)
True
class underworld.utils.SavedFileData(pyobj, filename)[source]

Bases: object

A class used to define saved data.

Parameters:
  • pyobj (object) – python object saved data relates to.
  • filename (str) – filename for saved data, full path
class underworld.utils.SolveLinearSystem(AMat, bVec, xVec, **kwargs)[source]

Bases: underworld._stgermain.StgCompoundComponent

To solve for x in the equation Ax=b

solve()[source]

Solve system