proteus.MeshTools module

Tools for creating and manipulating 1,2, and 3D meshes.

Inheritance diagram of proteus.MeshTools

class proteus.MeshTools.Node(nodeNumber=0, x=0.0, y=0.0, z=0.0)[source]

Bases: object

A numbered point in 3D Euclidean space

Variables
  • N – node number

  • p – Euclidean coordinates

Comparison operators and a hash value are defined using the 3-tuple of coordinates. This allows using Node objects and tuples of node objects as dictionary keys, but in that use case one should be careful not to modify the node coordinates.

>>> n0 = Node(nodeNumber=0,x=0.0,y=0.0,z=0.0)
>>> n1 = Node(nodeNumber=1,x=1.0,y=1.0,z=1.0)
>>> n1 >= n0
True
xUnitVector = array([1., 0., 0.])[source]
yUnitVector = array([0., 1., 0.])[source]
zUnitVector = array([0., 0., 1.])[source]
computeGeometricInfo()[source]
class proteus.MeshTools.Element(elementNumber=0, nodes=[])[source]

Bases: object

An numbered polytope in R^n

Variables
  • N – element number

  • nodes – sorted tuple of nodes defining the polytope

The nodes data member can be used as a dictionary key for the polytope as long as the nodes aren’t later modified.

class proteus.MeshTools.Edge(edgeNumber=0, nodes=[])[source]

Bases: proteus.MeshTools.Element

xUnitVector = array([1., 1., 0.])[source]
yUnitVector = array([0., 1., 0.])[source]
zUnitVector = array([0., 0., 1.])[source]

1D Element–a line connecting two Nodes

The nodes are stored as a lexicographically sorted node list.

computeGeometricInfo()[source]
proteus.MeshTools.getNodesFromEdges(edges)[source]

Extract the subset of nodes from a list of edges.

class proteus.MeshTools.Polygon(polygonNumber=0, nodes=[])[source]

Bases: proteus.MeshTools.Element

An abstract 2D element–a closed set of Edges connecting a set of Nodes.

The nodes and edges are stored as lexicographically sorted lists.

proteus.MeshTools.getEdgesFromPolygons(polygons)[source]

Extract the subset of edges from a list of polygons

class proteus.MeshTools.Triangle(triangleNumber=0, nodes=[], edgeDict=None)[source]

Bases: proteus.MeshTools.Polygon

A 2D triangular element

edgeMap = {(0, 1): 2, (0, 2): 1, (1, 2): 0}[source]
zUnitVector = array([0., 0., 1.])[source]
computeGeometricInfo()[source]
class proteus.MeshTools.Quadrilateral(quadrilateralNumber=0, edges=[], simple=True)[source]

Bases: proteus.MeshTools.Polygon

A 2D quadrilateral element

sortNodes(nodeList)[source]
computeGeometricInfo()[source]
class proteus.MeshTools.Polyhedron(polyhedronNumber=0, nodes=[])[source]

Bases: proteus.MeshTools.Element

An abstract 3D Element–a closed set of Polygons connecting a set of Edges.

The nodes and edges are stored as lexicographically sorted lists.

class proteus.MeshTools.Tetrahedron(tetrahedronNumber, nodes, edgeDict=None, triangleDict=None)[source]

Bases: proteus.MeshTools.Polyhedron

A 3D tetrahedral element

triangleMap = {(0, 1, 2): 3, (0, 1, 3): 2, (0, 2, 3): 1, (1, 2, 3): 0}[source]
edgeMap = {(0, 1): 0, (0, 2): 1, (0, 3): 2, (1, 2): 3, (1, 3): 4, (2, 3): 5}[source]
computeGeometricInfo()[source]
class proteus.MeshTools.Hexahedron(HN, quadrilaterals)[source]

Bases: proteus.MeshTools.Polyhedron

A 3D hexahedral element

class proteus.MeshTools.MeshParallelPartitioningTypes[source]

Bases: object

fake an enum for parallel partitioning options

element = 0[source]
node = 1[source]
class proteus.MeshTools.Mesh[source]

Bases: object

A partition of a domain in R^n into elements.

This is the base class for meshes. Contains routines for plotting the edges of the mesh in Matlab

Variables

elementBoundariesArray (array type) – This array lists the global edge number associated with every edge or face of an element.

partitionMesh(nLayersOfOverlap=1, parallelPartitioningType=1)[source]
partitionMeshFromFiles(filebase, base, nLayersOfOverlap=1, parallelPartitioningType=1)[source]
writeMeshXdmf(ar, name='', t=0.0, init=False, meshChanged=False, Xdmf_ElementTopology='Triangle', tCount=0, EB=False)[source]
buildFromC(cmesh)[source]
buildFromCNoArrays(cmesh)[source]
buildNodeStarArrays()[source]
buildArraysFromLists()[source]
computeGeometricInfo()[source]
buildMatlabMeshDataStructures(meshFileBase='meshMatlab', writeToFile=True)[source]

build array data structures for matlab finite element mesh representation and write to a file to view and play with in matlatb. The current matlab support is mostly for 2d, but this will return basic arrays for 1d and 3d too

in matlab can then print mesh with

pdemesh(p,e,t)

if one has pdetoolbox where

p is the vertex or point matrix e is the edge matrix, and t is the element matrix

e will be the elementBoundary matrix in 1d and 3d, but perhaps should remain the edge array?

points matrix is [nd x num vertices]
format :

row 1 = x coord, row 2 = y coord for nodes in mesh row 3 = z coord for nodes in mesh …

edge matrix is [2*nd+3 x num faces]
format:

row 1 = start vertex number … row nd = end vertex number row nd+1 = start value in edge parameterization, should be 0 row nd+2 = next value in edge parameterization, should be 1 or 2 row nd+nd= end value in edge parameterization, should be 2 or 1 row 2*nd+1 = global face id, base 1 row 2*nd+2 = subdomain on left? always 1 for now row 2*nd+3 = subdomain on right? always 0 for now

element matrix is [nd+2 x num elements]

row 1 = vertex 1 global number row 2 = vertex 2 global number … row nd+1 = vertex 3 global number row 4 = triangle subdomain number

where 1,2,3 is a local counter clockwise numbering of vertices in

triangle

writeEdgesMatlab(filename)[source]

store coordinates in files formatted for Matlab

viewTetrahedraMatlab(filename)[source]

plot the edges

viewMeshMatlab(filename)[source]

plot the edges

writeEdgesGnuplot(filename)[source]

store coordinates in files formatted for Matlab

writeEdgesGnuplot2(filename)[source]

store coordinates in files formatted for Matlab

viewMeshGnuplot(filename)[source]
viewMeshGnuplotPipe(filename)[source]
viewMeshGnuplotPipePar(filenames)[source]
convertFromPUMI(domain, MeshAdapt, faceList, regList, parallel=False, dim=3)[source]
class proteus.MeshTools.MultilevelMesh(levels=1)[source]

Bases: proteus.MeshTools.Mesh

A hierchical multilevel mesh

buildFromC(cmultilevelMesh)[source]
refine()[source]
locallyRefine(elementTagArray)[source]
buildArrayLists()[source]
calculateElementParents(recalculate=False)[source]
get array elementParents[l,e] = e_c, where element e_c is the parent of element e

elementParents[0,:] = -1

class proteus.MeshTools.PointMesh(points)[source]

Bases: proteus.MeshTools.Mesh

0D mesh

class proteus.MeshTools.EdgeGrid(nx=2, Lx=1.0)[source]

Bases: proteus.MeshTools.Mesh

A 1D regular grid on an interval

class proteus.MeshTools.QuadrilateralGrid(nx=2, ny=2, Lx=1.0, Ly=1.0)[source]

Bases: proteus.MeshTools.Mesh

A 2D regular grid of quadrilateral cells

class proteus.MeshTools.RectangularGrid(nx=1, ny=1, nz=1, Lx=1.0, Ly=1.0, Lz=1.0)[source]

Bases: proteus.MeshTools.Mesh

A regular partition into rectangles.

Nodes, edges, and faces can be indexed by (i,j,k) as follows. The edges and faces are divided according to orientation (i.e. x-edge…). An (i,j,k) index is associated with the type of edge or face having node (i,j,k) as the first node in a lexicographically sorted list of nodes corresponding to the edge or face.

getNodeNumber(i, j, k)[source]
getNode(i, j, k)[source]
getXEdgeNumber(ie, je, ke)[source]
getYEdgeNumber(ie, je, ke)[source]
getZEdgeNumber(ie, je, ke)[source]
getXEdge(ie, je, ke)[source]
getYEdge(ie, je, ke)[source]
getZEdge(ie, je, ke)[source]
getXYQuadrilateralNumber(ih, jh, kh)[source]
getXZQuadrilateralNumber(ih, jh, kh)[source]
getYZQuadrilateralNumber(ih, jh, kh)[source]
getXYQuadrilateral(ih, jh, kh)[source]
getXZQuadrilateral(ih, jh, kh)[source]
getYZQuadrilateral(ih, jh, kh)[source]
getHexahedronNumber(iH, jH, kH)[source]
getHexahedron(iH, jH, kH)[source]
getColor(i, j, k)[source]
refine(oldMesh, refineFactorX=2, refineFactorY=2, refineFactorZ=2)[source]
class proteus.MeshTools.MultilevelRectangularGrid(levels, nx, ny=1, nz=1, Lx=1.0, Ly=1.0, Lz=1.0, refinementLevels=1)[source]

Bases: proteus.MeshTools.MultilevelMesh

A hierarchical multilevel grid

refine()[source]
class proteus.MeshTools.TetrahedralMesh[source]

Bases: proteus.MeshTools.Mesh

A mesh of tetrahedra.

The nodes, edges, triangles, and tetrahedra are indexed by their node tuples. The corresponding lists are derived from the dictionaries, and sorted lexicographically. The global node numbers are redefined to give a lexicographic ordering.

The mesh can be generated from a rectangular grid and refined using either 4T or Freudenthal-Bey global refinement.

Variables
  • elementNodesArray (array_like) – A list of lists storing the node values associated with each element in the triangulation. The first index refers to the element number, while the second index refers to the global node value.

  • nodeArray (array_like) – A list of lists storing node coordinates. The first index referes to the global node number, while the second index refers to the x, y and z coordinates of the node respectively.

meshType()[source]
computeGeometricInfo()[source]
generateTetrahedralMeshFromRectangularGrid(nx, ny, nz, Lx, Ly, Lz)[source]
rectangularToTetrahedral6T(grid)[source]
rectangularToTetrahedral5T(grid)[source]
rectangularToTetrahedral(grid)[source]
fixLocalNumbering()[source]
finalize()[source]
buildLists()[source]
buildListsNodes()[source]
buildListsEdges()[source]
buildListsTriangles()[source]
buildListsTetrahedra()[source]
buildBoundaryMaps()[source]

Extract a mapping tn -> list((TN,tnLocal)) that provides all elements with the boundary face (triangle) tn and the local triangle number for that triangle. Likewise build mappings for edges and nodes Also extract a list of the triangles with only one associate element; these are the external boundary triangles. Then extract the edges and nodes from the boundary triangles.

newTetrahedron(nodes)[source]
registerEdges(t)[source]
registerTriangles(T)[source]
registerNode(node)[source]
readMeshADH(filename, adhBase=1)[source]
writeMeshXdmf(ar, name='', t=0.0, init=False, meshChanged=False, tCount=0, EB=False)[source]
writeMeshEnsight(filename, description=None)[source]
appendMeshEnsight(filename, description=None)[source]
writeMeshADH(filename, adhBase=1)[source]
writeBoundaryFacesADH(filename, adhBase=1)[source]
writeBoundaryNodesADH(filename, adhBase=1)[source]
refine4T(oldMesh)[source]
refineFreudenthalBey(oldMesh)[source]
refine(oldMesh)[source]
generateFromTetgenFiles(filebase, base, skipGeometricInit=False, parallel=False)[source]
generateFrom3DMFile(filebase, base=1)[source]
writeTetgenFiles(filebase, base)[source]
meshInfo()[source]
class proteus.MeshTools.HexahedralMesh[source]

Bases: proteus.MeshTools.Mesh

A mesh of hexahedra.

meshType()[source]
computeGeometricInfo()[source]
generateHexahedralMeshFromRectangularGrid(nx, ny, nz, Lx, Ly, Lz)[source]
finalize()[source]
buildLists()[source]
buildListsNodes()[source]
buildListsEdges()[source]
buildListsFaces()[source]
buildListsElems()[source]
buildBoundaryMaps()[source]

Extract a mapping tn -> list((TN,tnLocal)) that provides all elements with the boundary face tn and the local triangle number for that face Likewise build mappings for edges and nodes Also extract a list of the triangles with only one associate element; these are the external boundary triangles. Then extract the edges and nodes from the boundary triangles.

registerEdges(t)[source]
registerFaces(T)[source]
registerNode(node)[source]
meshInfo()[source]
writeMeshXdmf(ar, name='', t=0.0, init=False, meshChanged=False, tCount=0, EB=False)[source]
generateFromHexFile(filebase, base=0)[source]
class proteus.MeshTools.Mesh2DM(filename, adhBase=1)[source]

Bases: proteus.MeshTools.Mesh

A triangular mesh based on an ADH 3dm file

buildEdgeArrays()[source]
writeBoundaryMeshADH(filename, adhBase=1)[source]
writeMeshEnsight(filename, description=None)[source]
writeMeshXdmf(ar, name='', t=0.0, init=False, meshChanged=False, Xdmf_ElementTopology='Triangle', tCount=0)[source]
class proteus.MeshTools.Mesh3DM(filename, adhBase=1)[source]

Bases: proteus.MeshTools.Mesh

A Mesh for reading in tetrahedral meshes in the .3dm format

buildTriangleArrays()[source]
buildEdgeArray()[source]
writeBoundaryMeshADH(filename, adhBase=1)[source]
writeMeshEnsight(filename, description=None)[source]
writeBoundaryMeshEnsight(filename, description=None)[source]
writeMeshXdmf(ar, name='', t=0.0, init=False, meshChanged=False, Xdmf_ElementTopology='Tetrahedron', tCount=0)[source]
class proteus.MeshTools.MultilevelTetrahedralMesh(nx, ny, nz, x=0.0, y=0.0, z=0.0, Lx=1.0, Ly=1.0, Lz=1.0, refinementLevels=1, skipInit=False, nLayersOfOverlap=1, parallelPartitioningType=1)[source]

Bases: proteus.MeshTools.MultilevelMesh

A hierarchical multilevel mesh with tetrahedral cells

generateFromExistingCoarseMesh(mesh0, refinementLevels, nLayersOfOverlap=1, parallelPartitioningType=1)[source]
generatePartitionedMeshFromPUMI(mesh0, refinementLevels, nLayersOfOverlap=1)[source]
generatePartitionedMeshFromTetgenFiles(filebase, base, mesh0, refinementLevels, nLayersOfOverlap=1, parallelPartitioningType=1)[source]
refine()[source]
computeGeometricInfo()[source]
class proteus.MeshTools.MultilevelHexahedralMesh(nx, ny, nz, px=0, py=0, pz=0, x=0.0, y=0.0, z=0.0, Lx=1.0, Ly=1.0, Lz=1.0, refinementLevels=1, skipInit=False, nLayersOfOverlap=1, parallelPartitioningType=1)[source]

Bases: proteus.MeshTools.MultilevelMesh

A hierarchical multilevel mesh with hexahedral cells

generateFromExistingCoarseMesh(mesh0, refinementLevels, nLayersOfOverlap=1, parallelPartitioningType=1)[source]
refine()[source]
computeGeometricInfo()[source]
proteus.MeshTools.buildReferenceSimplex(nd=2)[source]

Create and return a Proteus mesh object for the reference element.

Parameters

nd (int) – Dimension of reference element

Returns

mesh – Simplex mesh

Return type

proteus.MeshTools.TriangularMesh

class proteus.MeshTools.TriangularMesh[source]

Bases: proteus.MeshTools.Mesh

A mesh of triangles

The nodes, edges, and triangles are indexed by their node tuples. The corresponding lists are derived from the dictionaries, and sorted lexicographically. The global node numbers are redefined to give a lexicographic ordering.

The mesh can be generated from a rectangular grid and refined using either 3t or Freudenthal-Bey global refinement.

meshType()[source]
computeGeometricInfo()[source]
generateTriangularMeshFromRectangularGrid(nx, ny, Lx, Ly, triangleFlag=1)[source]
rectangularToTriangularOriented(grid)[source]
rectangularToTriangularOrientedOtherWay(grid)[source]
rectangularToTriangularRedBlack(grid)[source]
rectangularToTriangular(grid)[source]
generateFromTriangleMesh(ctrirep, base)[source]
generateFromTriangleFiles(filebase, base)[source]
writeTriangleFiles(filebase, base)[source]
generateFrom2DMFile(filebase, base=1)[source]
constructTriangularMeshOnRectangle(Lx, Ly, nx, ny, writeMesh=0, meshFileBase='mesh2d')[source]

wrapper function for making a triangular mesh on the rectangle [0,Lx] x [0,Ly].

buildFromSets(triangleSet, edgeSet, nodeSet)[source]
fixLocalNumbering()[source]
finalize()[source]
buildLists()[source]
buildListsNodes()[source]
buildListsEdges()[source]
buildListsTriangles()[source]
newTriangle(nodes)[source]
registerEdges(t)[source]
registerNode(node)[source]
buildLevelSetMesh(value, nodalValues)[source]
refine3t(oldMesh)[source]
refineFreudenthalBey(oldMesh)[source]
refine(oldMesh)[source]
meshInfo()[source]
readMeshADH(filename, adhBase=1, suffix='3dm')[source]
writeMeshXdmf(ar, name='', t=0.0, init=False, meshChanged=False, tCount=0, EB=False)[source]
writeMeshEnsight(filename, description=None)[source]
appendMeshEnsight(filename, description=None)[source]
writeMeshADH(filename, adhBase=1)[source]
writeAsymptote(fileprefix, L, x, units='m')[source]

Write a representation of the triangular mesh in the Asymptote vector graphics language

class proteus.MeshTools.QuadrilateralMesh[source]

Bases: proteus.MeshTools.Mesh

A mesh of quads

The nodes, edges, and triangles are indexed by their node tuples. The corresponding lists are derived from the dictionaries, and sorted lexicographically. The global node numbers are redefined to give a lexicographic ordering.

The mesh can be generated from a rectangular grid and refined using either 3t or Freudenthal-Bey global refinement.

buildFromSets(faceSet, edgeSet, nodeSet)[source]
rectangularToQuadrilateral(grid, x=0.0, y=0.0, z=0.0)[source]

WIP - I think this is the first function that needs to be written so that MultilevelQuadrilateralMesh can work. This function does not call C functions.

generateQuadrilateralMeshFromRectangularGrid(nx, ny, Lx, Ly)[source]
generateFromQuadFileIFISS(meshfile)[source]

WIP - read a matlab.mat file containing IFISS vertices and elements

meshType()[source]
meshInfo()[source]
newQuadrilateral(edges)[source]
registerEdges(q)[source]

check if an edge is in the mesh dictionary if it is, point to existing entry otherwise, create a new entry

registerNode(node)[source]

check if a node is in the mesh dictionary if it is, point to existing entry otherwise, create a new entry

refine(oldMesh)[source]
finalize()[source]

WIP

buildLists()[source]

WIP

buildListsNodes()[source]
buildListsEdges()[source]
buildListsQuadrilaterals()[source]
writeMeshXdmf(ar, name='', t=0.0, init=False, meshChanged=False, tCount=0, EB=False)[source]
buildNodeDiameterArray()[source]
class proteus.MeshTools.MultilevelTriangularMesh(nx, ny, nz, x=0.0, y=0.0, z=0.0, Lx=1.0, Ly=1.0, Lz=1.0, refinementLevels=1, skipInit=False, nLayersOfOverlap=1, parallelPartitioningType=1, triangleFlag=0)[source]

Bases: proteus.MeshTools.MultilevelMesh

A hierarchical multilevel mesh of triangular cells

cmeshTools = <module 'proteus.cmeshTools' from '/home/travis/build/erdc/proteus/proteus/cmeshTools.cpython-310-x86_64-linux-gnu.so'>[source]
generateFromExistingCoarseMesh(mesh0, refinementLevels, nLayersOfOverlap=1, parallelPartitioningType=1)[source]
generatePartitionedMeshFromPUMI(mesh0, refinementLevels, nLayersOfOverlap=1)[source]
generatePartitionedMeshFromTriangleFiles(filebase, base, mesh0, refinementLevels, nLayersOfOverlap=1, parallelPartitioningType=1)[source]
refine()[source]
computeGeometricInfo()[source]
locallyRefine(elementTagArray, flagForRefineType=0)[source]

simple local refinement assuming elementTagArray[eN]=1 –> bisect

flagForRefineType = 0 – newest node, 1 – 4T, 2 – U4T

class proteus.MeshTools.MultilevelQuadrilateralMesh(nx, ny, nz, x=0.0, y=0.0, z=0.0, Lx=1.0, Ly=1.0, Lz=1.0, refinementLevels=1, skipInit=False, nLayersOfOverlap=1, parallelPartitioningType=1, triangleFlag=0, useC=True)[source]

Bases: proteus.MeshTools.MultilevelMesh

A heirarchical multilevel mesh of quadrilaterals WIP

refine()[source]
class proteus.MeshTools.InterpolatedBathymetryMesh(domain, triangleOptions, atol=0.0001, rtol=0.0001, maxElementDiameter=None, maxLevels=20, maxNodes=100000, bathyType='points', bathyAssignmentScheme='interpolation', errorNormType='L2', refineType=0)[source]

Bases: proteus.MeshTools.MultilevelTriangularMesh

A triangular mesh that interpolates bathymetry from a point cloud

setMeshBathymetry(mesh)[source]
setMeshBathymetry_interpolate(mesh)[source]
setMeshBathymetry_localAveraging(mesh)[source]

calculate the arithmetic mean bathymetry of points inside each triangle and then assign the area-weighted average of the element means to each node

locatePoints(mesh)[source]

locate the element containing each point

this should only be used on very coarse meshes

locatePoints_refined(mesh)[source]

locate the element containing each point

this should only be used on very coarse meshes

locatePoints_initial(mesh)[source]

locate the element containing each point

first find the nearest node, then loop over that node’s elements

interpolateBathymetry()[source]

interpolate bathymetry for the refinement from the parent mesh

tagElements(mesh)[source]

loop over points and calculate whether the interpolation error is within the tolerance

this should only be used on very coarse meshes

class proteus.MeshTools.EdgeMesh[source]

Bases: proteus.MeshTools.Mesh

A mesh of edges

The nodes, and edges are indexed by their node tuples. The corresponding lists are derived from the dictionaries, and sorted lexicographically. The global node numbers are redefined to give a lexicographic ordering.

computeGeometricInfo()[source]
generateEdgeMeshFromRectangularGrid(nx, Lx)[source]
rectangularToEdge(grid)[source]
finalize()[source]
buildLists()[source]
buildListsNodes()[source]
buildListsEdges()[source]
newEdge(nodes)[source]
registerNode(node)[source]
refine2e(oldMesh)[source]
refine(oldMesh)[source]
meshInfo()[source]
writeMeshADH(filename)[source]
writeMeshXdmf(ar, name='', t=0.0, init=False, meshChanged=False, tCount=0)[source]
writeMeshEnsight(filename, description=None)[source]
class proteus.MeshTools.MultilevelEdgeMesh(nx, ny, nz, x=0.0, y=0.0, z=0.0, Lx=1.0, Ly=1.0, Lz=1.0, refinementLevels=1, nLayersOfOverlap=1, parallelPartitioningType=1)[source]

Bases: proteus.MeshTools.MultilevelMesh

A hierarchical multilevel mesh of intervals (edges)

cmeshTools = <module 'proteus.cmeshTools' from '/home/travis/build/erdc/proteus/proteus/cmeshTools.cpython-310-x86_64-linux-gnu.so'>[source]
refine()[source]
computeGeometricInfo()[source]
locallyRefine(elementTagArray)[source]

simple local refinement assuming elementTagArray[eN]=1 –> bisect

class proteus.MeshTools.MultilevelSimplicialMesh(nd, nx, ny=1, nz=1, Lx=1.0, Ly=1.0, Lz=1.0, refinementLevels=1)[source]

Bases: proteus.MeshTools.MultilevelMesh

A wrapper for all the simplicial hierarchical meshes in 1,2, and 3D

refine()[source]
proteus.MeshTools.findXMLgridElement(xmf, MeshTag='Spatial_Domain', id_in_collection=- 1, verbose=0)[source]

Try to find the element of the xml tree xmf that holds a uniform grid with the name given in MeshTag by searching through Temporal Grid Collections and Grid Collections.

If MeshTag isn’t found, uses the first entry in the Domain

proteus.MeshTools.extractPropertiesFromXdmfGridNode(Grid)[source]

unpack the Topology, Geometry, NodeMaterials, and ElementMaterials nodes from xdmf node for a uniform grid

proteus.MeshTools.readUniformElementTopologyFromXdmf(elementTopologyName, Topology, hdf5, topologyid2name, topology2nodes)[source]

Read xmdf element topology information when there are uniform elements in the mesh Type of element given by elementTopologyName Heavy data stored in hdf5 topologyid2name – lookup for number of nodes in a given element type

returns

nElements_global – the number of elements in the mesh nNodes_element – number of nodes per element elementNodesArray – element –> node connectivity stored as flattened array accessed using elementNodes_offset elementNodes_offset – offsets into the elementNodesArray storage for element connectivity, element eN nodes are in elementNodesArray[elementNodes_offset[eN]:elementNodes_offset[eN+1]]

proteus.MeshTools.readMixedElementTopologyFromXdmf(elementTopologyName, Topology, hdf5, topologyid2name, topology2nodes)[source]

Read xmdf element topology information when there are mixed elements in the mesh Heavy data stored in hdf5 topologyid2name – lookup for number of nodes in a given element type

returns

nElements_global – the number of elements in the mesh elementNodesArray – element –> node connectivity stored as flattened array accessed using elementNodes_offset elementNodes_offset – offsets into the elementNodesArray storage for element connectivity, element eN nodes are inelementNodesArray[elementNodes_offset[eN]:elementNodes_offset[eN+1]]

proteus.MeshTools.readMeshXdmf(xmf_archive_base, heavy_file_base, MeshTag='Spatial_Domain', hasHDF5=True, verbose=0)[source]

Read in a mesh from XDMF, assuming heavy data is in hdf5

Returns

a BasicMeshInfo object with the minimal information read

proteus.MeshTools.writeHexMesh(mesh_info, hexfile_base, index_base=0)[source]

Write a hex mesh in Ido’s format with base numbering index_base HEX nNodes_global nElements_global x0 y0 z0 x1 y1 z1 … xN yN zN [n0 n1 n2 n3 n4 n5 n6 n7 mat0] [n0 n1 n2 n3 n4 n5 n6 n7 mat1]

class proteus.MeshTools.MultilevelNURBSMesh(nx, ny, nz, x=0.0, y=0.0, z=0.0, px=1, py=1, pz=1, Lx=1.0, Ly=1.0, Lz=1.0, refinementLevels=1, skipInit=False, nLayersOfOverlap=1, parallelPartitioningType=1)[source]

Bases: proteus.MeshTools.MultilevelMesh

generateFromExistingCoarseMesh(mesh0, refinementLevels, nLayersOfOverlap=1, parallelPartitioningType=1)[source]
refine()[source]
computeGeometricInfo()[source]
class proteus.MeshTools.NURBSMesh[source]

Bases: proteus.MeshTools.HexahedralMesh

A mesh consisting of NURBS.

generateHexahedralMeshFromRectangularGrid(nx, ny, nz, Lx, Ly, Lz)[source]
generateNURBSMeshFromRectangularGrid(nx, ny, nz, px, py, pz, Lx, Ly, Lz)[source]
proteus.MeshTools.distance(a, b)[source]
proteus.MeshTools.triangleVerticesToNormals(elementVertices)[source]

Given a set of vertices to a triangle, return normals and a point corresponding to each normal

proteus.MeshTools.tetrahedronVerticesToNormals(elementVertices)[source]

Given a set of vertices to a tetrahedron, return normals and a point corresponding to each normal

proteus.MeshTools.intersectPoints(line, points)[source]

Given a line segment (defined as two points), identify all points that the line segment intersects.

This hasn’t been vectorized.

proteus.MeshTools.intersectEdges(line, edges)[source]

Given a line segment (defined as two points), identify the locations of its intersections with all given edges (defined as line segments). If the line and an edge overlap, the furthest point along the line (closest to the second point) that is still on each edge is returned.

This hasn’t been vectorized.

proteus.MeshTools.intersectPolyhedron(line, polyhedron)[source]

Given a line (defined as two points), identify the locations that it enters and exits the polyhedron (defined as a collection of half-planes in three-space in normal, vertex form)

If the facets of the polyhedron are in edge form, the normal can be computed by taking the cross product of any two non-parallel edges of the facet (in three-space). Any vertex of the facet will work.

Implementation of algorithm described here: http://geomalgorithms.com/a13-_intersect-4.html

This hasn’t been vectorized.

proteus.MeshTools.getMeshIntersections(mesh, toPolyhedron, endpoints)[source]

Return all intersections between a line segment and a Proteus mesh

:param mesh - a Proteus mesh :param toPolyhedron - a method for converting Proteus element vertices to polyhedra in normal/point form :param endpoints - a pair of points in 3-space defining the line segment

:return a list of pairs of intersections through the mesh

proteus.MeshTools.runTriangle(polyfile, baseFlags='Yp', name='')[source]

Generate tetgen files from a polyfile.

Parameters
  • polyfile (str) – Filename with appropriate data for tengen.

  • baseFlags (str) – Standard Tetgen options for generation

  • name (str) –

proteus.MeshTools.runTetgen(polyfile, baseFlags='Yp', name='')[source]

Generate tetgen files from a polyfile.

Parameters
  • polyfile (str) – Filename with appropriate data for tengen.

  • baseFlags (str) – Standard Tetgen options for generation

  • name (str) –

proteus.MeshTools.genMeshWithTriangle(polyfile, nbase=1)[source]

Generate a mesh from a set of triangle files.

Parameters
  • polyfile (str) – Filename base for triangle files

  • nbase (int) –

Returns

mesh – Simplex mesh

Return type

proteus.MeshTools.TriangularMesh

proteus.MeshTools.genMeshWithTetgen(polyfile, nbase=1)[source]

Generate a mesh from a set of tetgen files.

Parameters
  • polyfile (str) – Filename base for tetgen files

  • nbase (int) –

Returns

mesh – Simplex mesh

Return type

proteus.MeshTools.TetrahedralMesh

class proteus.MeshTools.MeshOptions(nd=None)[source]

Bases: object

Mesh options for the domain

Parameters

nd (2 for 2D, 3 for 3D) –

setElementSize(he)[source]

Sets element size for uniform mesh.

Parameters

he (float) – mesh characteristic element size

setParallelPartitioningType(partitioning_type='node', layers_overlap=0)[source]

Changes parallel partitioning type

Parameters
  • partitioning_type (Optional[str, int]) – parallel partitioning type (default: ‘node’ (1))

  • layers (int) – layers of overlap for paralllel (default: 0)

setTriangleOptions(triangleOptions=None)[source]

Sets the trangle options

Parameters

triangle_options (Optional[str]) – string for triangle options. If not passed, it will be set with triangle_string attribute and ‘he’ value, with default for 2D: he**2/2; default for 3D: he**3/6

setMeshGenerator(generator)[source]

Indicates mesh generator to use

Parameters
  • generator (str) – options: ‘gmsh’, ‘triangle’, ‘tetgen’

  • current ((!) Only has an effect when setting to 'gmsh' in) –

  • 2D (implementation (triangle is default for) –

  • 3D) (tetgen for) –

setOutputFiles(name='mesh', poly=True, ply=False, asymptote=False, geo=False)[source]

Output files to be created

Parameters
  • name (Optional[str]) – name of the mesh files (prefix) (default: ‘mesh’)

  • poly (Optional[bool]) – create a poly file

  • ply (Optional[bool]) – create a ply file

  • asymptote (Optional[bool]) – create an asymptote file

  • geo – create a geofile

proteus.MeshTools.msh2simplex(fileprefix, nd)[source]

Converts a .msh file (Gmsh) to .ele .edge .node files (triangle). (!) Works only with triangle elements in 2D and tetrahedral elements in 3D.

Parameters

fileprefix (str) – prefix of the .msh file (e.g. ‘mesh’ if file called ‘mesh.msh’)

proteus.MeshTools.generateMesh(physics, numerics, generatePartitionedMeshFromFiles=False)[source]