sfepy.discrete.fem.mesh module

class sfepy.discrete.fem.mesh.Mesh(name='mesh', cmesh=None)[source]

The Mesh class is a light proxy to CMesh.

Input and output is handled by the MeshIO class and subclasses.

property coors
copy(name=None)[source]

Make a deep copy of the mesh.

Parameters:
namestr

Name of the copied mesh.

create_conn_graph(verbose=True)[source]

Create a graph of mesh connectivity.

Returns:
graphcsr_matrix

The mesh connectivity graph as a SciPy CSR matrix.

static from_data(name, coors, ngroups, conns, mat_ids, descs, nodal_bcs=None)[source]

Create a mesh from mesh IO data.

static from_file(filename=None, io='auto', prefix_dir=None, omit_facets=False, file_format=None)[source]

Read a mesh from a file.

Parameters:
filenamestring or function or MeshIO instance or Mesh instance

The name of file to read the mesh from. For convenience, a mesh creation function or a MeshIO instance or directly a Mesh instance can be passed in place of the file name.

io*MeshIO instance

Passing *MeshIO instance has precedence over filename.

prefix_dirstr

If not None, the filename is relative to that directory.

omit_facetsbool

If True, do not read cells of lower dimension than the space dimension (faces and/or edges). Only some MeshIO subclasses support this!

static from_region(region, mesh_in, localize=False, is_surface=False, tdim=None)[source]

Create a mesh corresponding to cells, or lower dimensional entities according tdim parameter, of a given region. If is_surface is True, then tdim = dim - 1.

get_bounding_box()[source]
get_cmesh(desc)[source]
get_conn(desc, ret_cells=False, tdim=None)[source]

Get the rectangular cell-vertex connectivity corresponding to desc. If ret_cells is True, the corresponding cells are returned as well.

transform_coors(mtx_t, ref_coors=None)[source]

Transform coordinates of the mesh by the given transformation matrix.

Parameters:
mtx_tarray

The transformation matrix T (2D array). It is applied depending on its shape:

  • (dim, dim): x = T * x

  • (dim, dim + 1): x = T[:, :-1] * x + T[:, -1]

ref_coorsarray, optional

Alternative coordinates to use for the transformation instead of the mesh coordinates, with the same shape as self.coors.

write(filename=None, io=None, out=None, float_format=None, file_format=None, **kwargs)[source]

Write mesh + optional results in out to a file.

Parameters:
filenamestr, optional

The file name. If None, the mesh name is used instead.

ioMeshIO instance or ‘auto’, optional

Passing ‘auto’ respects the extension of filename.

outdict, optional

The output data attached to the mesh vertices and/or cells.

float_formatstr, optional

The format string used to print floats in case of a text file format.

**kwargsdict, optional

Additional arguments that can be passed to the MeshIO instance.

sfepy.discrete.fem.mesh.find_map(x1, x2, allow_double=False, join=True)[source]

Find a mapping between common coordinates in x1 and x2, such that x1[cmap[:,0]] == x2[cmap[:,1]]

sfepy.discrete.fem.mesh.fix_double_nodes(coor, ngroups, conns)[source]

Detect and attempt fixing double nodes in a mesh.

The double nodes are nodes having the same coordinates w.r.t. precision given by eps.

sfepy.discrete.fem.mesh.get_min_vertex_distance(coor, guess)[source]

Can miss the minimum, but is enough for our purposes.

sfepy.discrete.fem.mesh.get_min_vertex_distance_naive(coor)[source]
sfepy.discrete.fem.mesh.make_mesh(coor, ngroups, conns, mesh_in)[source]

Create a mesh reusing mat_ids and descs of mesh_in.

sfepy.discrete.fem.mesh.merge_mesh(x1, ngroups1, conn1, mat_ids1, x2, ngroups2, conn2, mat_ids2, cmap)[source]

Merge two meshes in common coordinates found in x1, x2.

Notes

Assumes the same number and kind of element groups in both meshes!

sfepy.discrete.fem.mesh.set_accuracy(eps)[source]