underworld.mesh module

Implementation relating to meshing.

Classes

underworld.mesh.FeMesh The FeMesh class provides the geometry and topology of a finite element discretised domain.
underworld.mesh.FeMesh_Cartesian This class generates a finite element mesh which is topologically cartesian and geometrically regular.
underworld.mesh.FeMesh_IndexSet This class ties the FeMesh instance to an index set, and stores other metadata relevant to the set.
underworld.mesh.MeshVariable The MeshVariable class generates a variable supported by a finite element mesh.
class underworld.mesh.FeMesh(elementType, generator=None, **kwargs)[source]

Bases: underworld._stgermain.StgCompoundComponent, underworld.function._function.FunctionInput

The FeMesh class provides the geometry and topology of a finite element discretised domain. The FeMesh is implicitly parallel. Some aspects may be local or global, but this is generally handled automatically.

A number of element types are supported.

Parameters:
  • elementType (str) – Element type for FeMesh. See FeMesh.elementType docstring for further info.
  • generator (underworld.mesh.MeshGenerator) – Generator object which builds the FeMesh. See FeMesh.generator docstring for further info.
add_post_deform_function(function)[source]

Adds a function function to be executed after mesh deformation is applied.

Parameters:function (FunctionType) – Python (not underworld) function to be executed. Closures should be used where parameters are required.
add_pre_deform_function(function)[source]

Adds a function function to be executed before mesh deformation is applied.

Parameters:function (FunctionType) – Python (not underworld) function to be executed. Closures should be used where parameters are required.
add_variable(nodeDofCount, dataType='double', **kwargs)[source]

Creates and returns a mesh variable using the discretisation of the given mesh.

To set / read nodal values, use the numpy interface via the ‘data’ property.

Parameters:
  • dataType (string) – The data type for the variable. Note that only ‘double’ type variables are currently supported.
  • nodeDofCount (int) – Number of degrees of freedom per node the variable will have
Returns:

The newly created mesh variable.

Return type:

underworld.mesh.MeshVariable

Example

>>> linearMesh  = uw.mesh.FeMesh_Cartesian( elementType='Q1/dQ0', elementRes=(16,16), minCoord=(0.,0.), maxCoord=(1.,1.) )
>>> scalarFeVar = linearMesh.add_variable( nodeDofCount=1, dataType="double" )
>>> q0field     = linearMesh.subMesh.add_variable( 1 )  # adds variable to secondary elementType discretisation
data

Numpy proxy array proxy to underlying object vertex data. Note that the returned array is a proxy for all the local vertices, and it is provided as 1d list.

As these arrays are simply proxies to the underlying memory structures, no data copying is required.

Note that this property returns a read-only numpy array as default. If you wish to modify mesh vertex locations, you are required to use the deform_mesh context manager.

If you are modifying the mesh, remember to modify any submesh associated with it accordingly.

Returns:The data proxy array.
Return type:numpy.ndarray

Example

>>> import underworld as uw
>>> someMesh = uw.mesh.FeMesh_Cartesian( elementType='Q1', elementRes=(2,2), minCoord=(-1.,-1.), maxCoord=(1.,1.) )
>>> someMesh.data.shape
(9, 2)

You can retrieve individual vertex locations

>>> someMesh.data[1]
array([ 0., -1.])

You can modify these locations directly, but take care not to tangle the mesh! Mesh modifications must occur within the deform_mesh context manager.

>>> with someMesh.deform_mesh():
...    someMesh.data[1] = [0.1,-1.1]
>>> someMesh.data[1]
array([ 0.1, -1.1])
data_elementNodes
Returns:Array specifying the nodes (global node id) for a given element (local element id). NOTE: Length is local size.
Return type:numpy.ndarray
data_elgId
Returns:Array specifying global element ids. Length is domain size, (local+shadow).
Return type:numpy.ndarray
data_nodegId
Returns:Array specifying global node ids. Length is domain size, (local+shadow).
Return type:numpy.ndarray
deform_mesh(isRegular=False, remainsRegular=None)[source]

Any mesh deformation must occur within this python context manager. Note that certain algorithms may be switched to their irregular mesh equivalents (if not already set this way). This may have performance implications.

Any submesh will also be appropriately updated on return from the context manager, as will various mesh metrics.

Note that this method must be called collectively by all processes, irrespective of whether any given process does or does not need to deform any mesh nodes.

Parameters:isRegular (bool) – The general assumption is that the deformed mesh will no longer be regular (orthonormal), and more general but less efficient algorithms will be selected via this context manager. To over-ride this behaviour, set this parameter to True.

Notes

This method must be called collectively by all processes.

Example

>>> import underworld as uw
>>> someMesh = uw.mesh.FeMesh_Cartesian()
>>> with someMesh.deform_mesh():
...     someMesh.data[0] = [0.1,0.1]
>>> someMesh.data[0]
array([ 0.1,  0.1])
elementType
Returns:Element type for FeMesh. Supported types are “Q2”, “Q1”, “dQ1”, “dPc1” and “dQ0”.
Return type:str
elementsDomain
Returns:Returns the number of domain (local+shadow) elements on the mesh
Return type:int

Example

>>> someMesh = uw.mesh.FeMesh_Cartesian( elementType='Q1', elementRes=(2,2), minCoord=(-1.,-1.), maxCoord=(1.,1.) )
>>> someMesh.elementsDomain
4
elementsGlobal
Returns:Returns the number of global elements on the mesh
Return type:int

Example

>>> someMesh = uw.mesh.FeMesh_Cartesian( elementType='Q1', elementRes=(2,2), minCoord=(-1.,-1.), maxCoord=(1.,1.) )
>>> someMesh.elementsGlobal
4
elementsLocal
Returns:Returns the number of local elements on the mesh
Return type:int

Example

>>> someMesh = uw.mesh.FeMesh_Cartesian( elementType='Q1', elementRes=(2,2), minCoord=(-1.,-1.), maxCoord=(1.,1.) )
>>> someMesh.elementsLocal
4
generator

Getter/Setter for the mesh MeshGenerator object.

Returns:Object which builds the mesh. Note that the mesh itself may be a generator, in which case this property will return the mesh object iself.
Return type:underworld.mesh.MeshGenerator
load(filename)[source]

Load the mesh from disk.

Parameters:filename (str) – The filename for the saved file. Relative or absolute paths may be used, but all directories must exist.

Notes

This method must be called collectively by all processes.

If the file data array is the same length as the current mesh global size, it is assumed the file contains compatible data. Note that this may not be the case where for example where you have saved using a 2*8 resolution mesh, then loaded using an 8*2 resolution mesh.

Provided files must be in hdf5 format, and use the correct schema.

Example

Refer to example provided for ‘save’ method.

nodesDomain
Returns:Returns the number of domain (local+shadow) nodes on the mesh.
Return type:int
nodesGlobal
Returns:Returns the number of global nodes on the mesh
Return type:int

Example

>>> someMesh = uw.mesh.FeMesh_Cartesian( elementType='Q1', elementRes=(2,2), minCoord=(-1.,-1.), maxCoord=(1.,1.) )
>>> someMesh.nodesGlobal
9
nodesLocal
Returns:Returns the number of local nodes on the mesh.
Return type:int
reset()[source]

Reset the mesh.

Templated mesh (such as the DQ0 mesh) will be reset according to the current state of their geometryMesh.

Other mesh (such as the Q1 & Q2) will be reset to their post-construction state.

Notes

This method must be called collectively by all processes.

save(filename, units=None, **kwargs)[source]

Save the mesh to disk

Parameters:
  • filename (string) – The name of the output file.
  • units (pint unit object (optional)) – Define the units that must be used to save the data. The data will be dimensionalised and saved with the defined units. The units are saved as a HDF attribute.
  • keyword arguments are saved as string attributes. (Additional) –
Returns:

Data object relating to saved file. This only needs to be retained if you wish to create XDMF files and can be ignored otherwise.

Return type:

underworld.utils.SavedFileData

Notes

This method must be called collectively by all processes.

Example

First create the mesh:

>>> mesh = uw.mesh.FeMesh_Cartesian( elementType='Q1/dQ0', elementRes=(16,16), minCoord=(0.,0.), maxCoord=(1.,1.) )

Save to a file (note that the ‘ignoreMe’ object isn’t really required):

>>> ignoreMe = mesh.save("saved_mesh.h5")

Now let’s try and reload. First create new mesh (note the different spatial size):

>>> clone_mesh = uw.mesh.FeMesh_Cartesian( elementType='Q1/dQ0', elementRes=(16,16), minCoord=(0.,0.), maxCoord=(1.5,1.5) )

Confirm clone mesh is different from original mesh:

>>> import numpy as np
>>> np.allclose(mesh.data,clone_mesh.data)
False

Now reload using saved file:

>>> clone_mesh.load("saved_mesh.h5")

Now check for equality:

>>> np.allclose(mesh.data,clone_mesh.data)
True
>>> # clean up:
>>> if uw.mpi.rank == 0:
...     import os;
...     os.remove( "saved_mesh.h5" )
specialSets
Returns:This dictionary stores a set of special data sets relating to mesh objects.
Return type:dict

Example

>>> import underworld as uw
>>> someMesh = uw.mesh.FeMesh_Cartesian( elementType='Q1', elementRes=(2,2), minCoord=(0.,0.), maxCoord=(1.,1.) )
>>> sorted(someMesh.specialSets.keys())    # NOTE THAT WE ONLY SORT THIS LIST SO THAT RESULTS ARE DETERMINISTIC FOR DOCTESTS
['AllWalls_VertexSet', 'Bottom_VertexSet', 'Empty', 'Left_VertexSet', 'MaxI_VertexSet', 'MaxJ_VertexSet', 'MinI_VertexSet', 'MinJ_VertexSet', 'Right_VertexSet', 'Top_VertexSet']
>>> someMesh.specialSets["MinJ_VertexSet"]
FeMesh_IndexSet([0, 1, 2])
class underworld.mesh.FeMesh_Cartesian(elementType='Q1/dQ0', elementRes=(4, 4), minCoord=(0.0, 0.0), maxCoord=(1.0, 1.0), periodic=None, partitioned=True, **kwargs)[source]

Bases: underworld.mesh._mesh.FeMesh, underworld.mesh._mesh.CartesianMeshGenerator

This class generates a finite element mesh which is topologically cartesian and geometrically regular. It is possible to directly build a dual mesh by passing a pair of element types to the constructor.

Refer to parent classes for parameters beyond those below.

Parameters:
  • elementType (str) – Mesh element type. Note that this class allows the user to (optionally) provide a pair of elementTypes for which a dual mesh will be created. The submesh is accessible through the ‘subMesh’ property. The primary mesh itself is the object returned by this constructor.
  • elementRes (list,tuple) – List or tuple of ints specifying mesh resolution. See CartesianMeshGenerator.elementRes docstring for further information.
  • minCoord (list, tuple) – List or tuple of floats specifying minimum mesh location. See CartesianMeshGenerator.minCoord docstring for further information.
  • maxCoord (list, tuple) – List or tuple of floats specifying maximum mesh location. See CartesianMeshGenerator.maxCoord docstring for further information.
  • periodic (list, tuple) – List or tuple of bools, specifying mesh periodicity in each direction.
  • partitioned (bool) – If false, the mesh is not partitioned across entire processor pool. Instead mesh is entirely owned by processor which generated it.

Examples

To create a linear mesh:

>>> import underworld as uw
>>> someMesh = uw.mesh.FeMesh_Cartesian( elementType='Q1', elementRes=(16,16), minCoord=(0.,0.), maxCoord=(1.,1.) )
>>> someMesh.dim
2
>>> someMesh.elementRes
(16, 16)

Alternatively, you can create a linear/constant dual mesh

>>> someDualMesh = uw.mesh.FeMesh_Cartesian( elementType='Q1/dQ0', elementRes=(16,16), minCoord=(0.,0.), maxCoord=(1.,1.) )
>>> someDualMesh.elementType
'Q1'
>>> subMesh = someDualMesh.subMesh
>>> subMesh.elementType
'DQ0'

To set / read vertex coords, use the numpy interface via the ‘data’ property.

integrate(fn)[source]

Perform a domain integral of the given underworld function over this mesh

Parameters:mesh (uw.mesh.FeMesh_Cartesian) – Domain to perform integral over.

Examples

>>> mesh = uw.mesh.FeMesh_Cartesian(minCoord=(0.0,0.0), maxCoord=(1.0,2.0))
>>> fn_1 = uw.function.misc.constant(2.0)
>>> np.allclose( mesh.integrate( fn_1 )[0], 4 )
True
>>> fn_2 = uw.function.misc.constant(2.0) * (0.5, 1.0)
>>> np.allclose( mesh.integrate( fn_2 ), [2,4] )
True
subMesh
Returns:Returns the submesh where the object is a dual mesh, or None otherwise.
Return type:underworld.mesh.FeMesh
class underworld.mesh.FeMesh_IndexSet(object, topologicalIndex=None, *args, **kwargs)[source]

Bases: underworld.container._indexset.ObjectifiedIndexSet, underworld.function._function.FunctionInput

This class ties the FeMesh instance to an index set, and stores other metadata relevant to the set.

Parameters:
  • object (underworld.mesh.FeMesh) – The FeMesh instance from which the IndexSet was extracted.
  • topologicalIndex (int) – Mesh topological index for which the IndexSet relates. See docstring for further info.

Example

>>> amesh = uw.mesh.FeMesh_Cartesian( elementType='Q1/dQ0', elementRes=(4,4), minCoord=(0.,0.), maxCoord=(1.,1.) )
>>> iset = uw.libUnderworld.StgDomain.RegularMeshUtils_CreateGlobalMaxISet( amesh._mesh )
>>> uw.mesh.FeMesh_IndexSet( amesh, topologicalIndex=0, size=amesh.nodesGlobal, fromObject=iset )
FeMesh_IndexSet([ 4,  9, 14, 19, 24])
topologicalIndex
Returns:
The topological index for the indices. The mapping is:
0 - vertex 1 - edge 2 - face 3 - volume
Return type:int
class underworld.mesh.MeshVariable(mesh, nodeDofCount, dataType='double', **kwargs)[source]

Bases: underworld._stgermain.StgCompoundComponent, underworld.function._function.Function, underworld._stgermain.Save, underworld._stgermain.Load

The MeshVariable class generates a variable supported by a finite element mesh.

To set / read nodal values, use the numpy interface via the ‘data’ property.

Parameters:
  • mesh (underworld.mesh.FeMesh) – The supporting mesh for the variable.
  • dataType (string) – The data type for the variable. Note that only ‘double’ type variables are currently supported.
  • nodeDofCount (int) – Number of degrees of freedom per node the variable will have.

See property docstrings for further information.

Example

For example, to create a scalar meshVariable:

>>> linearMesh = uw.mesh.FeMesh_Cartesian( elementType='Q1/dQ0', elementRes=(16,16), minCoord=(0.,0.), maxCoord=(1.,1.) )
>>> scalarFeVar = uw.mesh.MeshVariable( mesh=linearMesh, nodeDofCount=1, dataType="double" )

or a vector meshvariable can be created:

>>> vectorFeVar = uw.mesh.MeshVariable( mesh=linearMesh, nodeDofCount=3, dataType="double" )
copy(deepcopy=False)[source]

This method returns a copy of the meshvariable.

Parameters:deepcopy (bool) – If True, the variable’s data is also copied into new variable.
Returns:The mesh variable copy.
Return type:underworld.mesh.MeshVariable

Example

>>> mesh = uw.mesh.FeMesh_Cartesian()
>>> var = uw.mesh.MeshVariable(mesh,2)
>>> import math
>>> var.data[:] = (math.pi,math.exp(1.))
>>> varCopy = var.copy()
>>> varCopy.mesh == var.mesh
True
>>> varCopy.nodeDofCount == var.nodeDofCount
True
>>> import numpy as np
>>> np.allclose(var.data,varCopy.data)
False
>>> varCopy2 = var.copy(deepcopy=True)
>>> np.allclose(var.data,varCopy2.data)
True
data

Numpy proxy array to underlying variable data. Note that the returned array is a proxy for all the local nodal data, and is provided as 1d list. It is possible to change the shape of this numpy array to reflect the cartesian topology (where appropriate), though again only the local portion of the decomposed domain will be available, and the shape will not necessarily be identical on all processors.

As these arrays are simply proxies to the underlying memory structures, no data copying is required.

Returns:The proxy array.
Return type:numpy.ndarray

Example

>>> linearMesh = uw.mesh.FeMesh_Cartesian( elementType='Q1/dQ0', elementRes=(16,16), minCoord=(0.,0.), maxCoord=(1.,1.) )
>>> scalarFeVar = uw.mesh.MeshVariable( mesh=linearMesh, nodeDofCount=1, dataType="double" )
>>> scalarFeVar.data.shape
(289, 1)

You can retrieve individual nodal values

>>> scalarFeVar.data[100]
array([ 0.])

Likewise you can modify nodal values

>>> scalarFeVar.data[100] = 15.333
>>> scalarFeVar.data[100]
array([ 15.333])
dataType
Returns:Data type for variable. Supported types are ‘double’.
Return type:str
fn_gradient

Returns a Function for the gradient field of this meshvariable.

Note that for a scalar variable T, the gradient function returns an array of the form:

\[[ \frac{\partial T}{\partial x}, \frac{\partial T}{\partial y}, \frac{\partial T}{\partial z} ]\]

and for a vector variable v:

\[[ \frac{\partial v_x}{\partial x}, \frac{\partial v_x}{\partial y}, \frac{\partial v_x}{\partial z}, \frac{\partial v_y}{\partial x}, \frac{\partial v_y}{\partial y}, \frac{\partial v_y}{\partial z}, \frac{\partial v_z}{\partial x}, \frac{\partial v_z}{\partial y}, \frac{\partial v_z}{\partial z} ]\]
Returns:The gradient function.
Return type:underworld.function.Function
load(filename, interpolate=False)[source]

Load the MeshVariable from disk.

Parameters:
  • filename (str) – The filename for the saved file. Relative or absolute paths may be used, but all directories must exist.
  • interpolate (bool) – Set to True to interpolate a file containing different resolution data. Note that a temporary MeshVariable with the file data will be build on each processor. Also note that the temporary MeshVariable can only be built if its corresponding mesh file is available. Also note that the supporting mesh mush be regular.

Notes

This method must be called collectively by all processes.

If the file data array is the same length as the current variable global size, it is assumed the file contains compatible data. Note that this may not be the case where for example where you have saved using a 2*8 resolution mesh, then loaded using an 8*2 resolution mesh.

Provided files must be in hdf5 format, and use the correct schema.

Example

Refer to example provided for ‘save’ method.

mesh
Returns:Supporting FeMesh for this MeshVariable.
Return type:underworld.mesh.FeMesh
nodeDofCount
Returns:Degrees of freedom on each mesh node that this variable provides.
Return type:int
save(filename, meshHandle=None, units=None, **kwargs)[source]

Save the MeshVariable to disk.

Parameters:
  • filename (string) – The name of the output file. Relative or absolute paths may be used, but all directories must exist.
  • meshHandle (uw.utils.SavedFileData , optional) – The saved mesh file handle. If provided, a link is created within the mesh variable file to this saved mesh file. Important for checkpoint when the mesh deforms.
  • units (pint unit object (optional)) – Define the units that must be used to save the data. The data will be dimensionalised and saved with the defined units. The units are saved as a HDF attribute. Note if units are in celsius (see scaling.pint_degc_labels) the data is scaled and save to degrees kelvin.
  • keyword arguments are saved as string attributes. (Additional) –

Notes

This method must be called collectively by all processes.

Returns:Data object relating to saved file. This only needs to be retained if you wish to create XDMF files and can be ignored otherwise.
Return type:underworld.utils.SavedFileData

Example

First create the mesh add a variable:

>>> mesh = uw.mesh.FeMesh_Cartesian( elementType='Q1/dQ0', elementRes=(16,16), minCoord=(0.,0.), maxCoord=(1.,1.) )
>>> var = uw.mesh.MeshVariable( mesh=mesh, nodeDofCount=1, dataType="double" )

Write something to variable

>>> import numpy as np
>>> var.data[:,0] = np.arange(var.data.shape[0])

Save to a file (note that the ‘ignoreMe’ object isn’t really required):

>>> meshHandle = mesh.save("saved_mesh.h5")
>>> ignoreMe = var.save("saved_mesh_variable.h5", meshHandle)
>>> ignoreMe = var.save("saved_mesh_variable.h5", meshHandle) #test dup

Now let’s try and reload.

>>> clone_var = uw.mesh.MeshVariable( mesh=mesh, nodeDofCount=1, dataType="double" )
>>> clone_var.load("saved_mesh_variable.h5")

Now check for equality:

>>> np.allclose(var.data,clone_var.data)
True

Now check the field can be loaded on a different mesh topology (interpolation)

>>> mesh19 = uw.mesh.FeMesh_Cartesian( elementType='Q1/dQ0', elementRes=(19,19), minCoord=(0.,0.), maxCoord=(1.,1.) )
>>> clone_var2 = mesh19.add_variable(1)
>>> clone_var2.load("saved_mesh_variable.h5", interpolate=True)
>>> np.allclose(mesh.integrate(var), mesh19.integrate(clone_var2))
True
>>> # clean up:
>>> if uw.mpi.rank == 0:
...     import os;
...     os.remove( "saved_mesh_variable.h5" )
...     os.remove( "saved_mesh.h5" )
syncronise()[source]

This method is often necessary when Underworld is operating in parallel.

It will syncronise the mesh variable so that it is consistent with it’s parallel neighbours. Specifically, the shadow space of each process obtains the required data from neighbouring processes.

xdmf(filename, fieldSavedData, varname, meshSavedData, meshname, modeltime=0.0)[source]

Creates an xdmf file, filename, associating the fieldSavedData file on the meshSavedData file

Notes

xdmf contain 2 files: an .xml and a .h5 file. See http://www.xdmf.org/index.php/Main_Page This method only needs to be called by the master process, all other processes return quiely.

Parameters:
  • filename (str) – The output path to write the xdmf file. Relative or absolute paths may be used, but all directories must exist.
  • varname (str) – The xdmf name to give the field
  • meshSavedData (underworld.utils.SaveFileData) – Handler returned for saving a mesh. underworld.mesh.save(xxx)
  • meshname (str) – The xdmf name to give the mesh
  • fieldSavedData (underworld.SavedFileData) – Handler returned from saving a field. underworld.mesh.save(xxx)
  • modeltime (float) – The time recorded in the xdmf output file

Example

First create the mesh add a variable:

>>> mesh = uw.mesh.FeMesh_Cartesian( elementType='Q1/dQ0', elementRes=(16,16), minCoord=(0.,0.), maxCoord=(1.,1.) )
>>> var = uw.mesh.MeshVariable( mesh=mesh, nodeDofCount=1, dataType="double" )

Write something to variable

>>> import numpy as np
>>> var.data[:,0] = np.arange(var.data.shape[0])

Save mesh and var to a file:

>>> meshDat = mesh.save("saved_mesh.h5")
>>> varDat = var.save("saved_mesh_variable.h5", meshDat)

Now let’s create the xdmf file

>>> var.xdmf("TESTxdmf", varDat, "var1", meshDat, "meshie" )

Does file exist?

>>> import os
>>> if uw.mpi.rank == 0: os.path.isfile("TESTxdmf.xdmf")
True

Clean up:

>>> if uw.mpi.rank == 0:
...     import os;
...     os.remove( "saved_mesh_variable.h5" )
...     os.remove( "saved_mesh.h5" )
...     os.remove( "TESTxdmf.xdmf" )