5. 3D Volume Mesh Generation - TetMesh, VolMesh
CEETRON Mesh supports automatic, generally unstructured 3D volume mesh generation using the TetMesh and VolMesh objects. The TetMesh object is strictly designed for tetrahedral meshing in which only triangular input faces (boundary and internal) are allowed. The VolMesh object is designed to allow for both triangular and quadrilateral faces to be input and allow pyramids to be generated from the boundary. If VolMesh is unable to completely fill a volume with pyramids the remainder of the volume is filled with tetrahedra.
5.1. 3D Tetrahedral Mesh Generation - TetMesh
The TetMesh module provides for the automatic generation of linear, parabolic or cubic tetrahedral elements given a surface mesh of linear, parabolic or cubic triangles. By default, the mesh generation procedure is guaranteed to not alter the surface mesh in any way. The node numbering of the surface mesh is preserved in the node numbering and connectivity of the resultant tetrahedral mesh. The input surface mesh is checked for closure, that is, all input surface triangles must be connected to at least one triangle across each edge.
The functions associated with a TetMesh object are the following.
Begin and end an instance of an object, return object error flag
msh_TetMeshBegin()
- create an instance of a TetMesh objectmsh_TetMeshEnd()
- destroy an instance of a TetMesh objectmsh_TetMeshError()
- return TetMesh object error flag
Define surface mesh, mesh parameters and generate mesh
msh_TetMeshAbort()
- set abort flagmsh_TetMeshDef()
- define number of nodes, tris.msh_TetMeshInq()
- inquire number of nodes, tris.msh_TetMeshConnect()
- set the contents of Connect objectmsh_TetMeshComputeVolume()
- compute volumemsh_TetMeshGenerate()
- generate finite element meshmsh_TetMeshGetConnect()
- get internal Connect objectmsh_TetMeshGetInteger()
- get integer mesh generation informationmsh_TetMeshSetEdge()
- set interior edgemsh_TetMeshSetEdgeAssoc()
- set interior edge associationsmsh_TetMeshSetEdgeStat()
- set interior edge statusmsh_TetMeshSetFunction()
- set a call back functionmsh_TetMeshSetParami()
- set mesh generation parametersmsh_TetMeshGetParami()
- get mesh generation parametersmsh_TetMeshSetParamd()
- set mesh generation parametersmsh_TetMeshGetParamd()
- get mesh generation parametersmsh_TetMeshRefine()
- refine a tetrahedral meshmsh_TetMeshSetElem()
- set an orphan mesh for refinementmsh_TetMeshSetGeomSizing()
- specify sizing within geometric shapemsh_TetMeshSetNode()
- define a nodemsh_TetMeshSetNodeAssoc()
- define a node associationmsh_TetMeshSetNodeSizing()
- set node sizingmsh_TetMeshSetNodeStat()
- set node statusmsh_TetMeshSetTri()
- define a trianglemsh_TetMeshSetTriAssoc()
- define a triangle associationmsh_TetMeshSetTriSizing()
- set sizing within a regionmsh_TetMeshSetTriDepthSizing()
- set sizing within a depth in a regionmsh_TetMeshSetTriStat()
- set triangle statusmsh_TetMeshWrite()
- write mesh description to filemsh_TetMeshRead()
- read mesh description from file
Instance a TetMesh object using msh_TetMeshBegin()
.
Once a TetMesh object is instanced,
define the number of nodes and triangles which define the surface mesh
as well as the order of the triangles (linear, parabolic or cubic)
using msh_TetMeshDef()
.
The surface mesh is then input as
a set of nodes in 3D space and a set of triangle connectivities.
The nodes are defined using msh_TetMeshSetNode()
,
the boundary triangle connectivities are defined using msh_TetMeshSetTri()
.
Each edge in the boundary triangulation must be connected to an even
number of boundary triangles - usually two.
Interior edges may be specified using msh_TetMeshSetEdge()
.
The interior edges may connect to the
exterior boundary triangulation.
Interior faces may be defined by specifying two “double backed” triangles
with opposite connectivity sense.
By default, the connectivity of the
triangles is such that the normal to the triangle is pointed out of the
volume. The direction of the normal is determined by the right-hand rule
applied to the triangle connectivity.
As a convenience, the function msh_TetMeshConnect()
will take
the contents of a Connect object and set it as the geometric input
to TetMesh. The typical case is using the Connect object which
contains the mesh generated by a previous call to msh_SurfMeshGenerate()
or
msh_SurfMeshRefine()
.
The function assumes certain “best practices” conventions have been
used with respect to the node and element associations which have been
set to identify the underlying geometry entities.
The function msh_TetMeshSetParami()
is used to set integer parameters which
affect the mesh generation process.
It is possible to specify
parabolic or cubic tetrahedra to be produced if the input triangles are linear.
In this case nodes are added to the edges lying on the
input surface triangulation to create the higher order tetrahedra throughout.
By default the output tetrahedra are the same order as the input triangles.
If parabolic or cubic triangles are input, then only parabolic or cubic
tetrahedra respectively are generated.
Also the sense of the normals of the input triangles may be changed so that
normals are expected to point into the volume.
The function msh_TetMeshSetParamd()
is used to set real double precision
parameters which control the quality and size of the generated tetrahedra.
A global target edge length for interior elements may be prescribed.
Specific, location dependent, sizes may be set at input nodes
using msh_TetMeshSetNodeSizing()
and for regions using
msh_TetMeshSetTriSizing()
.
By default, all input nodes must be connected to input boundary triangles
and interior edges and faces. If unconnected nodes are to be inserted
in the interior then the user must call msh_TetMeshSetParami()
to enable parameter type TETMESH_UNCONNECT.
The user is able to define integer associations at nodes and triangles
which will be assigned to the node and elements generated on the
respective entity.
These associations are useful for identifying nodes and elements
in the output mesh for the application of loads, boundary conditions, material
properties, etc.
The function msh_TetMeshSetNodeAssoc()
is used to set associations at
nodes generated at input points. The function msh_TetMeshSetNodeAssoc()
may be called any number of times for a single node with the same
association type. The association values are accumulated as multiple
association on the node.
The function msh_TetMeshSetTriAssoc()
is used to
set associations which will be attached to nodes generated along
triangle edges or on triangle faces. These associations will also
be generated as element entity associations along the respective
tetrahedron edges and faces connected to the input triangles.
The element associations are useful for tagging tetrahedra generated
with an enclosed boundary, etc.
Finally the function msh_TetMeshGenerate()
generates the
nodes and elements and enters them into a Connect object.
The mesh generation process may fail due to errors in the surface triangles
or failure of the tetrahedralization procedure. Use msh_TetMeshGetInteger()
to query for detailed information concerning errors in the surface triangles,
etc.
During the mesh generation process, the input surface triangulation and
triangle corner node locations (by default) remain unchanged.
However, the midside node locations
may be perturbed under certain circumstances to satisfy quality constraints
in the output tetrahedra.
The function msh_TetMeshSetParami()
is
used to set the parameter
VIS_MESH_MIDSIDEQUAL which controls the overall movement of midside
nodes. By default all midside nodes are susceptible to movement.
The function msh_TetMeshWrite()
is provided to write a complete description
of the defined input surface mesh and meshing parameters to a file. The primary
use of this file is to encapsulate cases
in which the TetMesh module should fail in some respect.
This output file can then be
made available to Tech Soft 3D for failure diagnosis of TetMesh
meshing algorithms.
The TetMesh module provides for a monitor
callback function to be specified which
is called intermittently during the mesh generation process in
msh_TetMeshGenerate()
. The primary purpose of this function is to
allow the user to interrupt or interrogate
the on-going mesh generation process for any reason. If mesh generation is to
be terminated, call msh_TetMeshAbort()
. This call will set an internal
flag which will cause msh_TetMeshGenerate()
to terminate the mesh generation
process and return.
In addition, a user defined sizing callback function may be set.
Set callback functions using msh_TetMeshSetFunction()
.
The following code fragment illustrates the basic framework of using the TetMesh module to generate linear tetrahedral finite elements in a volume bounded by the surface defined by a set of linear triangular facets.
/* declare objects */
msh_TetMesh *tetmesh;
vis_Connect *connect;
/* input surface triangles */
Vint numnode, numtris;
Vdouble xpnt[][3] = { ... };
Vint ixtri[][3] = { ... };
Vint ierr;
/* create tetmesh object */
tetmesh = msh_TetMeshBegin ();
msh_TetMeshDef (tetmesh,numnode,numtris,2);
/* define nodes and triangles */
for (i = 0; i < numnode; i++) {
msh_TetMeshSetNode (tetmesh,i+1,xpnt[i]);
}
for (i = 0; i < numtris; i++) {
msh_TetMeshSetTri (tetmesh,i+1,ixtri[i]);
}
/* create empty connect object to hold generated mesh */
connect = vis_ConnectBegin ();
vis_ConnectDef (connect,0,0);
/* generate mesh */
msh_TetMeshGenerate (tetmesh,connect);
/* process any errors */
ierr = msh_TetMeshError (tetmesh);
if(ierr) {
...
}
/* access generated nodes and elements from Connect */
...
/* delete objects */
vis_ConnectEnd (connect);
msh_TetMeshEnd (tetmesh);
5.2. Non-manifold Geometry, Voids, Inclusions and Hard Nodes and Edges
The TetMesh module supports the solid meshing of non-manifold geometries containing voids and inclusions. Any number of internal hard nodes and edges may be enforced in the final mesh.
The geometry input to the TetMesh module consists of nodes, boundary
triangles and enforced edges and faces. All input nodes are considered
to be hard points and by default are not moved. All boundary
triangles are input using msh_TetMeshSetTri()
and may include
external boundary triangles (including the boundaries of voids)
and internal boundary triangles representing inclusions
and other internal surfaces. All triangles representing internal
surfaces must be “double backed”, ie for each triangle there must
be an associated coincident triangle with opposite connectivity sense.
Boundary triangles are subject to change during refinement operations.
Internal hard edges which are enforced in the final
mesh are input using msh_TetMeshSetEdge()
.
5.3. Mesh Refinement and Adaptivity
A tetrahedral mesh can be refined using msh_TetMeshRefine()
.
Mesh refinement is performed on the mesh resulting from the
call to msh_TetMeshGenerate()
or the last previous call to
msh_TetMeshRefine()
.
The refinement process usually occurs as a result of a need to change
element sizing due to the results of a solution performed on the mesh.
A new mesh size is specified using a vis_State object of sizes per element
in the mesh to be refined. If an element size is set to zero, the element
will not be refined.
During refinement existing nodes are not moved. Refinement is performed
by splitting edges and performing allowable swaps.
All corner node identifiers
are retained, midside node identifiers are, in general, reassigned.
Additional nodes are generated to satisfy the refinement requirements.
Any input triangle or interior edge is not swapped,
any triangle edge with a specified
association is not swapped. By default all boundary triangles and edges
may be subdivided by the refinement process. Use the functions
msh_TetMeshSetTriStat()
and msh_TetMeshSetEdgeStat()
.
to mark entities which are not to be subdivided.
Mesh refinement is supported both for meshes generated
with msh_TetMeshGenerate()
or for orphan meshes.
In the case of orphan meshes, in addition to specifying
the triangle-discretized boundary and its associations, the orphan mesh is
specified by entering each original element with msh_TetMeshSetElem()
.
5.4. Function Descriptions
The currently available TetMesh functions are described in detail in this section.
-
msh_TetMesh *msh_TetMeshBegin(void)
create an instance of a TetMesh object
Create an instance of an TetMesh object. Memory is allocated for the object private data and the pointer to the data is returned.
Destroy an instance of a TetMesh object using
void msh_TetMeshEnd (msh_TetMesh *tetmesh)
Return the current value of a TetMesh object error flag using
Vint msh_TetMeshError (msh_TetMesh *tetmesh)
- Returns:
The function returns a pointer to the newly created TetMesh object. If the object creation fails, NULL is returned.
-
void msh_TetMeshEnd(msh_TetMesh *p)
destroy an instance of a TetMesh object
-
Vint msh_TetMeshError(msh_TetMesh *p)
return the current value of a TetMesh object error flag
-
void msh_TetMeshDef(msh_TetMesh *p, Vint numnode, Vint numtris, Vint maxi)
define number of nodes, tris.
Define the number of nodes and triangles which define the input surface mesh. The number of nodes on a surface triangle edge is specified by maxi. Define node coordinates using
msh_TetMeshSetNode()
and define triangle connectivities which reference the input nodes usingmsh_TetMeshSetTr()
.Inquire of defined numnode, numtris and maxi as output arguments using
void msh_TetMeshInq (msh_TetMesh *tetmesh, Vint *numnode, Vint *numtris, Vint *maxi)
- Errors
VIS_ERROR_VALUE
is generated if an improper maxi is input or numnode is less than or equal to zero or numtris is less than zero.
- Parameters:
p – Pointer to TetMesh object.
numnode – Number of nodes on input surface mesh
numtris – Number of triangles on input surface mesh
maxi – Number of nodes on triangle edge
=2 Linear surface triangles =3 Parabolic surface triangles =4 Cubic surface triangles
-
void msh_TetMeshInq(msh_TetMesh *p, Vint *numnode, Vint *numtris, Vint *maxi)
inquire of defined numnode, numtris and maxi as output arguments
See
msh_TetMeshDef()
-
void msh_TetMeshSetFunction(msh_TetMesh *p, Vint funtype, Vfunc *function, Vobject *object)
set a call back function
Set callback functions. By default the callback functions are NULL. A callback is not invoked if it is NULL.
The monitor callback function prototype is
The first argument is the TetMesh object, tetmesh, and the second is a user defined object, object.void function (msh_TetMesh *tetmesh, Vobject *object)
The sizing callback function prototype is
The first argument is the TetMesh object, tetmesh, the second is a user defined object, object, the third is the coordinate location x and the fourth is the returned size s.void function (msh_TetMesh *tetmesh, Vobject *object, Vdouble x[3], Vdouble *s)
The anisotropic sizing callback function prototype is
The first argument is the TetMesh object, tetmesh, the second is a user defined object, object, the third is the coordinate location x and the fourth is the returned anisotropic scaled orthogonal direction vectors. The first 3 components are the first size scaled direction, the next 3 components are the second size scaled direction, the next 3 components are the third size scaled direction.void function (msh_TetMesh *tetmesh, Vobject *object, Vdouble x[3], Vdouble s[3][3])
- Parameters:
p – Pointer to TetMesh object.
funtype – Type of callback function to set
x=SYS_FUNCTION_MONITOR Monitor callback =SYS_FUNCTION_SIZING Sizing callback =SYS_FUNCTION_ASIZING Anisotropic sizing callback
function – Pointer to callback function
object – Pointer to the object to be returned as function argument
-
void msh_TetMeshSetParami(msh_TetMesh *p, Vint ptype, Vint iparam)
set mesh generation parameters
Specify mesh generation parameters. These parameters specify the order and connectivity convention to be used for the generated tetrahedra.
The parameter
VIS_MESH_MAXI
sets the order of the generated tetrahedral elements by specifying the number of nodes to generate along an edge. A value of 2 generates linear tetrahedra, a value of 3 generates parabolic tetrahedra and a value of 4 generates cubic tetrahedra. By defaultVIS_MESH_MAXI
is set to the order of the input surface triangles defined bymsh_TetMeshDef()
. Note that the order of the generated tetrahedra can only be changed from the default if linear input surface triangles are input.The parameter
VIS_MESH_INWARD
sets the normal convention of the input triangles to point into the volume. By defaultVIS_MESH_INWARD
is set toSYS_OFF
, that is, the triangle normals are expected to point out of the volume.The parameter
VIS_MESH_DELAUNAY
requests that the final mesh satisfy the Delaunay criterion. By defaultVIS_MESH_DELAUNAY
is set toSYS_OFF
.The parameter
VIS_MESH_INTERREFINE
enables the refinement of the interior. If this parameter is disabled, the generation of zero interior points is not guaranteed. The placement of zero interior points is impossible in some cases. By defaultVIS_MESH_INTERREFINE
is set toSYS_ON
.The parameter
VIS_MESH_MIDSIDEQUAL
toggles the adjustment of midside nodes of parabolic and cubic elements to satisfy specified quality constraints. The midside nodes are adjusted to ensure that all node Jacobians are positive and that the Jacobian ratio in any element is less than the specified upper Jacobian ratio bound. This bound is specified usingmsh_TetMeshSetParamd()
. This procedure normally moves the midside node location of interior edges to satisfy quality constraints. It may move midside nodes on the surface triangulation in rare cases. By defaultVIS_MESH_MIDSIDEQUAL
is set toSYS_ON
.The parameter
VIS_MESH_FUNSIZINGMIN
is used to control the way final element sizing is computed when a sizing callback is installed. By default, the sizing returned by the sizing callback is used. IfVIS_MESH_FUNSIZINGMIN
is enabled, then the sizing is the minimum of the sizing callback and normal sizing. By defaultVIS_MESH_FUNSIZINGMIN
is set toVIS_OFF
.The parameter
VIS_MESH_FUNSIZINGLEGACY
is used to revert the sizing callback behavior to its legacy implementation. In this case the sizing is only applied to the nodes, and not to elements nor edges. By defaultVIS_MESH_FUNSIZINGLEGACY
is set toVIS_OFF
.The parameter
VIS_MESH_POINTSET
requests a Delaunay tetrahedralization of the input points. The mesh boundary is the convex hull of the point set. By defaultVIS_MESH_POINTSET
is set toSYS_OFF
.The parameter
VIS_MESH_REPAIRINCONS
toggles the repair of inconsistent connectivity of input triangles. If enabled, all input triangle connectivities are checked for consistent sense and enclosing positive total volume. The connectivities are then selectively reversed to enforce consistency and total positive volume. By defaultVIS_MESH_REPAIRINCONS
is set toSYS_OFF
.By default, all input nodes must be connected to input boundary triangles and interior edges and faces. If unconnected nodes are to be inserted in the interior then parameter
TETMESH_UNCONNECT
must be enabled. By defaultTETMESH_UNCONNECT
is set toSYS_OFF
.Inquire of defined integer parameters as output arguments using
void msh_TetMeshGetParami (msh_TetMesh *tetmesh, Vint ptype, Vint *iparam)
- Errors
VIS_ERROR_ENUM
is generated if an improper ptype is specified.VIS_ERROR_VALUE
is generated if an improper order of generated elements is specified.
- Parameters:
p – Pointer to TetMesh object.
ptype – Type of display parameter to set
x=VIS_MESH_MAXI Order of generated elements =VIS_MESH_FUNSIZINGMIN Size is min'ed with sizing callback =VIS_MESH_FUNSIZINGLEGACY Legacy sizing callback =VIS_MESH_INWARD Inward normal triangle flag. =VIS_MESH_DELAUNAY Delaunay mesh flag =VIS_MESH_INTERREFINE Interior refinement flag =VIS_MESH_MIDSIDEQUAL Quality midside node adjustment =VIS_MESH_POINTSET Delaunay tetrahedralization of point set =VIS_MESH_REPAIRINCONS Repair inconsistent input =TETMESH_UNCONNECT Toggle unconnected interior nodes
iparam – Specifies the integer value that ptype will be set to.
x=SYS_OFF Disable =SYS_ON Enable
-
void msh_TetMeshGetParami(msh_TetMesh *p, Vint type, Vint *iparam)
get mesh generation parameters
-
void msh_TetMeshSetParamd(msh_TetMesh *p, Vint ptype, Vdouble dparam)
set mesh generation parameters
Specify mesh generation parameters. These parameters are used to control the quality and size of the generated tetrahedra.
The parameter
VIS_MESH_EDGELENGTH
specifies a target edge length of generated tetrahedra. The edge length constraint may not be satisfied for tetrahedra generated near the surface mesh due to surface mesh sizing. By defaultVIS_MESH_EDGELENGTH
is set to 1.The parameter
VIS_MESH_GROWTHRATE
specifies a maximum growth rate. The growth rate governs the rate at which the size of adjacent elements can grow. The ratio of edge lengths of unconstrained adjacent elements will not exceed the specified growth rate. The growth rate must be greateer than 1.001. By defaultVIS_MESH_GROWTHRATE
is set to 2.0The parameter
VIS_MESH_MAXEDGELENGTH
specifies a maximum edge length of generated elements. The maximum edge length constraint may not be satisfied for tetrahedra generated near the surface mesh due to surface mesh sizing otherwise the maximum edge length is strictly enforced. By defaultVIS_MESH_MAXEDGELENGTH
is set to 0. and is ignored.The parameter
VIS_MESH_MAXEDGEALT
specifies an upper bound on the edge-altitude ratio of an element. The edge-altitude ratio is defined as the maximum edge length over the smallest altitude. A series of operations are performed in an attempt to satisfy the edge-altitude ratio quality. By defaultVIS_MESH_MAXEDGEALT
is set to 100.The parameter
VIS_MESH_REFINEFACTOR
controls the number of elements generated duringmsh_TetMeshRefine()
. It specifies a factor to be applied to the number of starting elements to yield the maximum number of elements to be allowed in the refinement step. By defaultVIS_MESH_REFINEFACTOR
is set to 10.Inquire of defined double precision parameters as output arguments using
void msh_TetMeshGetParamd (msh_TetMesh *tetmesh, Vint ptype, Vdouble *dparam)
- Errors
VIS_ERROR_ENUM
is generated if an improper ptype is specified.VIS_ERROR_VALUE
is generated if an improper value is specified.
- Parameters:
p – Pointer to TetMesh object.
ptype – Type of display parameter to set
x=VIS_MESH_EDGELENGTH Target edge length =VIS_MESH_GROWTHRATE Maximum growth rate =VIS_MESH_MAXEDGELENGTH Maximum edge length =VIS_MESH_MAXEDGEALT Maximum edge altitude ratio =VIS_MESH_MAXANGLE Maximum dihedral angle =VIS_MESH_REFINEFACTOR Control refine element count
dparam – Specifies the value that ptype will be set to.
-
void msh_TetMeshGetParamd(msh_TetMesh *p, Vint type, Vdouble *param)
get mesh generation parameters
-
void msh_TetMeshAbort(msh_TetMesh *p)
set abort flag
Set the contents of the Connect object as geometry input. All nodes and elements in the Connect object are set as input. All
VIS_GEOVERT
,VIS_GEOEDGE
andVIS_GEOFACE
associations on nodes are set. AllVIS_GEOEDGE
andVIS_GEOFACE
associations on triangles are set as element entity associations. AllVIS_GEOBODY
associations are set as element associations. AllVIS_GEOEDGE
associations on edges are set.- Errors
VIS_ERROR_OPERATION
is generated if the Connect object does not contain proper element types to serve as boundary input.
- Parameters:
p – Pointer to TetMesh object.
-
void msh_TetMeshSetNode(msh_TetMesh *p, Vint id, Vdouble x[3])
define a surface node
Set a node location. The node id must be in the interval 1 <= id <= numnode where numnode is defined by
msh_TetMeshDef()
. Input surface triangles are defined usingmsh_TetMeshSetTri()
.- Errors
VIS_ERROR_VALUE
is generated if an improper id is specified.
- Parameters:
p – Pointer to TetMesh object.
id – Node id
x – Node coordinate location
-
void msh_TetMeshSetNodeStat(msh_TetMesh *p, Vint id, Vint stat)
set surface node status
Set a surface corner or midside node movement status. By default all corner nodes may not be moved and all midside nodes may be moved. Set the status to 1 for any surface corner or midsides nodes which must not be moved for any reason. The default movement status of surface corner and midside nodes may be specified using
msh_TetMeshSetParami()
.- Errors
VIS_ERROR_VALUE
is generated if an improper id is specified.
- Parameters:
p – Pointer to TetMesh object.
id – Node id
stat – Surface node status
=0 Surface node may be moved =1 Surface node may not be moved
-
void msh_TetMeshSetNodeAssoc(msh_TetMesh *p, Vint type, Vint id, Vint aid)
define a node association
Set a node association. The association will be propagated to the node in the generated mesh. This function may be called any number of times for a particular node and association type, all association values are accumulated.
- Errors
VIS_ERROR_VALUE
is generated if an improper id is specified.
- Parameters:
p – Pointer to TetMesh object.
type – Type of association
id – Node id
aid – Association value
-
void msh_TetMeshSetNodeSizing(msh_TetMesh *p, Vint id, Vdouble h)
set node sizing
Set element size to be enforced at a node.
- Errors
VIS_ERROR_VALUE
is generated if an improper id is specified.
- Parameters:
p – Pointer to TetMesh object.
id – Node id
h – Sizing.
-
void msh_TetMeshSetTri(msh_TetMesh *p, Vint id, Vint ix[])
define a triangle
Set a surface triangle connectivity. The triangle id must be in the interval 1 <= id <= numtris where numtris is defined by
msh_TetMeshDef()
. The triangle connectivity, ix, must contain 3 nodes for linear surface triangles, 6 nodes for parabolic surface triangles and 9 nodes for cubic surface triangles. The connectivity convention for surface triangles is Serendipity (3 corner nodes followed by 3 or 6 midside nodes in counter clockwise direction). By default, the connectivity of the triangles is be such that the normal to the triangle is pointed out of the volume. The direction of the normal is determined by the right-hand rule applied to the triangle connectivity. The outward normal convention may be changed to an inward normal convention usingmsh_TetMeshSetParami()
. Note that all triangles must be entered sequentially.- Errors
VIS_ERROR_VALUE
is generated if an improper id is specified.VIS_ERROR_OPERATION
is generated if an id is entered out of order.
- Parameters:
p – Pointer to TetMesh object.
id – Triangle id
ix – Triangle connectivity
-
void msh_TetMeshSetTriAssoc(msh_TetMesh *p, Vint type, Vint id, Vint enttype, Vint no, Vint aid)
define a triangle association
Set associations to be assigned to generated tetrahedron elements, edges and faces. The entity type of
SYS_EDGE
will propagate the association to tetrahedron edges and all nodes on the triangle edge. The entity type ofSYS_FACE
will propagate the association to tetrahedron faces and all nodes on the triangle face. The entity type ofSYS_ELEM
will propagate the association to all tetrahedra adjacent to the specified triangle, bounded by adjacent input triangles. The value of no is ignored.- Errors
VIS_ERROR_ENUM
is generated if an improper type or enttype is specified.VIS_ERROR_VALUE
is generated if an improper id or no is specified.
- Parameters:
p – Pointer to TetMesh object.
type – Type of association
id – Triangle id
enttype – Entity type to which association is assigned
=SYS_EDGE Element edge association =SYS_FACE Element face association =SYS_ELEM Element association
no – Triangle edge or face number
aid – Association value
-
void msh_TetMeshSetTriStat(msh_TetMesh *p, Vint id, Vint no, Vint stat)
set triangle status
Set triangle status. By default triangle edges are not swapped or split during generation and may be split but not swapped during refinement.
- Errors
VIS_ERROR_VALUE
is generated if an improper id or no is specified.
- Parameters:
p – Pointer to TetMesh object.
id – Triangle id
no – Triangle edge number
stat – Status, bitwise OR
x=VIS_MESH_BIT_GENSPLIT Split during generation =VIS_MESH_BIT_GENSWAP Swap during generation =VIS_MESH_BIT_REFSPLIT Split during refinement =VIS_MESH_BIT_REFSWAP Swap during refinement
-
void msh_TetMeshSetTriSizing(msh_TetMesh *p, Vint id, Vdouble h)
set sizing within a region
Set element size to be enforced within a region. The region is defined as the material bound by input triangles that exists to the left of the specified triangle, ie. the normal of the triangle points out of the material region. If a size is specified on more than one triangle bounding a region, the minimum size is used for the region.
- Errors
VIS_ERROR_VALUE
is generated if an improper id is specified.
- Parameters:
p – Pointer to TetMesh object.
id – Triangle id
h – Sizing.
-
void msh_TetMeshSetTriDepthSizing(msh_TetMesh *p, Vint id, Vdouble d, Vdouble h)
set sizing within a depth in a region
Set element size to be enforced within a specified depth from the input boundary triangle.
- Errors
VIS_ERROR_VALUE
is generated if an improper id is specified.VIS_ERROR_VALUE
is generated if either d or h is less than zero.
- Parameters:
p – Pointer to TetMesh object.
id – Triangle id
d – Depth
h – Sizing
-
void msh_TetMeshSetEdge(msh_TetMesh *p, Vint id, Vint ix[])
set interior edge
Set an interior edge connectivity. The edge connectivity, ix, must contain 2 nodes for linear input surface triangles and 3 nodes for parabolic input surface triangles and 4 nodes for cubic input surface triangles. The nodes of an interior edge may be connected to the boundary.
- Errors
VIS_ERROR_VALUE
is generated if an improper id or ix is specified.
- Parameters:
p – Pointer to TetMesh object.
id – Edge id
ix – Edge connectivity
-
void msh_TetMeshSetEdgeStat(msh_TetMesh *p, Vint id, Vint stat)
set interior edge status
Set an interior edge status. By default interior edges are not swapped or split during generation and may be split but not swapped during refinement.
- Errors
VIS_ERROR_VALUE
is generated if an improper id is specified.
- Parameters:
p – Pointer to TetMesh object.
id – Edge id
stat – Status, bitwise OR
x=VIS_MESH_BIT_GENSPLIT Split during generation =VIS_MESH_BIT_GENSWAP Swap during generation =VIS_MESH_BIT_REFSPLIT Split during refinement =VIS_MESH_BIT_REFSWAP Swap during refinement
-
void msh_TetMeshSetEdgeAssoc(msh_TetMesh *p, Vint type, Vint id, Vint aid)
set interior edge associations
Set associations at interior input edges which will propagate to generated tetrahedron edges.
- Errors
VIS_ERROR_ENUM
is generated if an improper type is specified.VIS_ERROR_VALUE
is generated if an improper id is specified.
- Parameters:
p – Pointer to TetMesh object.
type – Type of association
id – Edge id
aid – Association value
-
void msh_TetMeshSetGeomSizing(msh_TetMesh *p, Vint type, Vdouble xo[3], Vdouble xa[3], Vdouble xb[3], Vdouble d[3], Vdouble h)
specify sizing within geometric shape
Set element size to be enforced within a geometric shape.
- Errors
VIS_ERROR_ENUM
is generated if an improper type is specified.
- Parameters:
p – Pointer to TetMesh object.
type – Geometry type
x=VIS_MESH_GEOMSPHERE Ellipsoid =VIS_MESH_GEOMCYLINDER Ellipsoidal cylinder =VIS_MESH_GEOMBOX Box
xo – Coordinates of geometry origin
xa – Coordinates of point along x’ axis
xb – Coordinates of point in x’y’ plane
d – Half length of geometric shape along each axis
h – Sizing
-
void msh_TetMeshGetConnect(msh_TetMesh *p, vis_Connect **connect)
get internal Connect object
Get the internal Connect object which contains the user input nodes and triangles. The Connect object must not be changed by the user.
-
void msh_TetMeshWrite(msh_TetMesh *p, Vint type, const Vchar *path)
write mesh description from file
Write a surface mesh description to a file which has been input to TetMesh. The contents of the file is intended for internal use only. This function should be called after all surface mesh input and meshing parameters have been set to ensure that an accurate surface mesh description including correct meshing parameters is written. The binary file format is preferred for accuracy and generally results in smaller file sizes. The primary use of this function is to create a file of the complete surface mesh description in the case that the TetMesh module fails to properly generate a mesh. This file can then be transmitted to Visual Kinematics for debugging purposes. It is suggested that a file extension of .tet be used for ASCII files and .btet be used for BINARY files.
- Errors
SYS_ERROR_ENUM
is generated if an improper type is specified.SYS_ERROR_FILE
is generated if the file can not be opened.SYS_ERROR_OPERATION
is generated if the input surface mesh has not been defined.
- Parameters:
p – Pointer to TetMesh object.
type – File type
x=SYS_ASCII ASCII format =SYS_BINARY Binary format
path – [out] File path
-
void msh_TetMeshRead(msh_TetMesh *p, Vint type, const Vchar *path)
read mesh description from file
Read a surface mesh description from a file which has been previously written by
msh_TetMeshWrite()
. The format of the file is not published and is intended for internal use only.- Errors
SYS_ERROR_ENUM
is generated if an improper type is specified.SYS_ERROR_FILE
is generated if the file can not be opened.SYS_ERROR_FORMAT
is generated if the file contents are not correctly formatted or if a file written by a later version of TetMesh is read.
- Parameters:
p – Pointer to TetMesh object.
type – File type
x=SYS_ASCII ASCII format =SYS_BINARY Binary format
path – [out] File path
-
void msh_TetMeshComputeVolume(msh_TetMesh *p, Vdouble *volume)
compute volume
Compute volume enclosed by input boundary.
- Parameters:
p – Pointer to TetMesh object.
volume – [out] Computed volume
-
void msh_TetMeshGenerate(msh_TetMesh *p, vis_Connect *connect)
generate finite element mesh
Generate interior nodes and tetrahedral elements. The generated nodes and elements are appended to any existing nodes and elements in the input Connect object. The input surface triangles are checked for free edges, connectivity consistency, degeneracy, enclosing negative or zero volume, etc. Use
msh_TetMeshGetInteger()
is query for detailed information about the mesh generation process including errors.- Errors
VIS_ERROR_OPERATION
is generated if mesh generation fails due to errors detected in the surface triangles.VIS_ERROR_COMPUTE
is generated if mesh generation fails during tetrahedralization.
-
void msh_TetMeshGetInteger(msh_TetMesh *p, Vint type, Vint iparam[])
get integer mesh generation information
Query for integer mesh generation information. This function is useful for determining where errors in the mesh generation procedure occurred during a call to
msh_TetMeshGenerate()
. UseTETMESH_NUMFREEEDGE
to determine the number of triangles which contained a free edge. The queryTETMESH_FREEEDGE
will return the list of all free edge triangles detected.The query type
TETMESH_NUMBOUNDINSERT
returns the number of nodes inserted on the boundary during mesh generation. This number is generally zero and will only be non-zero in the rare case that a node needed to be inserted on the boundary as a last resort during the boundary recovery phase. The queryTETMESH_BOUNDINSERT
will return the list of nodes inserted on the boundary.The query type
TETMESH_NUMNONMAN
returns the number of improper non-manifold boundary triangles detected in the input. This number will only include those triangles detected on the first non-manifold edge encountered, it will not generally include all non-manifold triangles. The queryTETMESH_NONMAN
will return the list of non-manifold triangles detected.The query type
TETMESH_NUMINCONS
returns the number of inconsistent triangles detected in the input. This number will only include those triangles detected on the first inconsistent edge encountered, it will not generally include all inconsistent triangles. The queryTETMESH_INCONS
will return the list of inconsistent triangles detected.The query type
TETMESH_NUMNEGJAC
returns the number of elements with negative corner Jacobians. This check is only performed if it is enabled usingmsh_TetMeshSetParami()
.The query type
TETMESH_NUMUNCONN
returns the number of unconnected nodes detected in the input. The queryTETMESH_UNCONN
will return the list of unconnected nodes.The query type
TETMESH_NUMTRIZEROANG
returns the number of input triangles with zero angle detected in the input. The queryTETMESH_TRIZEROANG
will return the list of zero angle triangles.The query type
TETMESH_NUMINTERSECT
returns the number of input faces which intersect. The queryTETMESH_INTERSECT
will return the list of intersecting triangles.The query type
TETMESH_NUMINTERSECTPAIR
returns the number of pairs of intersecting triangles. The queryTETMESH_INTERSECTPAIR
will return the list of pairs of intersecting triangles. Note that the array of triangles returned needs large enough to hold 2 times the number of triangle pairs triangles.The query type
TETMESH_NUMUNREC
returns the number of input faces which were unable to be recovered. The queryTETMESH_UNREC
will return the list of unrecovered triangles.The query type
TETMESH_NUMNODEINSERT
returns the number of input nodes which were unable to be inserted. The queryTETMESH_NODEINSERT
will return the list of nodes which were not inserted.The query type
VIS_MESH_PROGRESS
returns information about the current state of the meshing process. This function is most useful when called from the monitor function. The following four integer values are returned.iparam[0] Meshing phase Phase 1 - Pre-processing boundary data Phase 2 - Boundary recovery Phase 3 - Generating interior points Phase 4 - Mesh cleanup Phase 5 - Mesh quality improvement Phase 101 - Refine Pre-processing Phase 102 - Refine point insertion Phase 103 - Refine Mesh quality improvement iparam[1] Current number of nodes iparam[2] Current number of elements iparam[3] Estimate percent completion
- Errors
VIS_ERROR_ENUM
is generated if an improper type is specified.
- Parameters:
p – Pointer to TetMesh object.
type – Type of integer information to query
x=TETMESH_NUMBOUNDINSERT Number of nodes inserted on boundary =TETMESH_BOUNDINSERT Nodes inserted on boundary =TETMESH_NUMFREEEDGE Number of surface triangle free edges =TETMESH_FREEEDGE Free edge surface triangles =TETMESH_NUMINCONS Number of inconsistent triangles =TETMESH_INCONS Inconsistent surface triangles =TETMESH_NUMINTERSECT Number of intersecting triangles =TETMESH_INTERSECT Intersecting triangles =TETMESH_NUMINTERSECTPAIR No. of intersecting triangle pairs =TETMESH_INTERSECTPAIR Intersecting triangle pairs =TETMESH_NUMNONMAN Number of non manifold triangles =TETMESH_NONMAN Non manifold surface triangles =TETMESH_NUMNEGJAC Number of negative Jacobians =TETMESH_NUMUNCONN Number of unconnected nodes =TETMESH_UNCONN Unconnected nodes =TETMESH_NUMTRIZEROANG Number of zero angle triangles =TETMESH_TRIZEROANG Zero angle triangles =TETMESH_NUMUNREC Number of unrecovered triangles =TETMESH_UNREC Unrecovered triangles =TETMESH_NUMNODEINSERT Number of nodes not inserted =TETMESH_NODEINSERT Not inserted nodes =VIS_MESH_PROGRESS Current meshing progress
iparam – [out] Returned integer information
-
void msh_TetMeshSetElem(msh_TetMesh *p, Vint id, Vint shape, Vint maxi, Vint ix[])
set an orphan mesh for refinement
Input elements for an orphan tetrahedron mesh to be refined. The API has been defined in a general way to support other element shapes. However, only tetrahedron elements, linear, parabolic or cubic, are currently supported.
- Errors
VIS_ERROR_OPERATION
is generated if TetMesh was previously used to generate a mesh.msh_TetMeshSetElem()
should only be used for orphan meshes.VIS_ERROR_VALUE
is generated if id’s are entered in a non-sequential manner.VIS_ERROR_VALUE
is generated if shape is notVIS_SHAPETET
.VIS_ERROR_VALUE
is generated if maxi is not equal to 0, 2, 3 or 4.
- Parameters:
p – Pointer to TetMesh object.
id – Element index
shape – Element shape. Only SYS_SHAPETET is supported.
maxi – Order of the elements in the initial mesh.
ix – Element connectivity
-
void msh_TetMeshRefine(msh_TetMesh *p, vis_State *state, vis_Connect *connect)
refine a tetrahedral mesh
Refine a tetrahedral mesh. Use the State object to specify the new mesh size for each element. As the refinement takes place TetMesh will use the initial mesh and state to interpolate new mesh sizes at any point in space. Any element with a size of zero in the state object will not be refined.
If an orphan mesh is to be refined - i.e., a mesh not generated with TetMesh, use
msh_TetMeshSetElem()
to set the initial tetrahedron mesh. The boundary of the orphan mesh must be specified withmsh_TetMeshSetTri()
and, if needed,msh_TetMeshSetTriAssoc()
.Since TetMesh has no knowledge of preserved edges it uses all edges with an edge association as an indication of a preserved edge.
- Errors
VIS_ERROR_OPERATION
is generated if TetMesh has not been used to generate the original mesh or an orphan mesh is to be refined andmsh_TetMeshSetElem()
has not been called.VIS_ERROR_OPERATION
is generated if state is not defined as an element state or if the state’s size differs from the existing mesh.VIS_ERROR_OPERATION
is generated if state is not defined with a data typeVIS_SCALAR
.VIS_ERROR_OPERATION
can also be generated if auxiliary, internal tools cannot be created.
5.5. 3D General Volume Mesh Generation - VolMesh
The VolMesh module provides for the automatic generation of complex hybrid surface and volumetric elements given a surface mesh of linear, parabolic or cubic triangles and/or quadrilaterals. By default, the mesh generation procedure is guaranteed to not alter the surface mesh in any way. The node numbering of the surface mesh is preserved in the node numbering and connectivity of the resultant mesh. The input surface mesh is checked for closure, that is, all input surface faces must be connected to at least one face across each edge.
The functions associated with a VolMesh object are the following.
Begin and end an instance of an object, return object error flag
msh_VolMeshBegin()
- create an instance of a VolMesh objectmsh_VolMeshEnd()
- destroy an instance of a VolMesh objectmsh_VolMeshError()
- return VolMesh object error flag
Define surface mesh, mesh parameters and generate mesh
msh_VolMeshAbort()
- set abort flagmsh_VolMeshDef()
- define number of nodes, faces.msh_VolMeshInq()
- inquire number of nodes, faces.msh_VolMeshConnect()
- set the contents of Connect objectmsh_VolMeshGenerate()
- generate finite element meshmsh_VolMeshGetConnect()
- get internal Connect objectmsh_VolMeshGetInteger()
- get integer mesh generation informationmsh_VolMeshSetEdge()
- set interior edgemsh_VolMeshSetEdgeStat()
- set interior edge statusmsh_VolMeshSetFunction()
- set a call back functionmsh_VolMeshSetParami()
- set mesh generation parametersmsh_VolMeshGetParami()
- get mesh generation parametersmsh_VolMeshSetParamd()
- set mesh generation parametersmsh_VolMeshGetParamd()
- get mesh generation parametersmsh_VolMeshRefine()
- refine a meshmsh_VolMeshSetGeomSizing()
- specify sizing within geometric shapemsh_VolMeshSetNode()
- define a nodemsh_VolMeshSetNodeAssoc()
- define a node associationmsh_VolMeshSetNodeSizing()
- set node sizingmsh_VolMeshSetNodeStat()
- set node statusmsh_VolMeshSetFace()
- define a facemsh_VolMeshSetFaceAssoc()
- define a face associationmsh_VolMeshSetFaceFill()
- define a bounding face fill typemsh_VolMeshSetFaceSizing()
- set sizing within a regionmsh_VolMeshSetFaceDepthSizing()
- set sizing within a depth in a regionmsh_VolMeshSetFaceStat()
- set face statusmsh_VolMeshWrite()
- write mesh description to filemsh_VolMeshRead()
- read mesh description from file
Instance a VolMesh object using msh_VolMeshBegin()
.
Once a VolMesh object is instanced,
define the number of nodes and faces which define the surface mesh
using msh_VolMeshDef()
.
The surface mesh is then input as
a set of nodes in 3D space and a set of face (triangle and/or quadrilateral)
connectivities.
The nodes are defined using msh_VolMeshSetNode()
,
the boundary face connectivities are defined using msh_VolMeshSetFace()
.
Each edge in the boundary must be connected to an even
number of boundary faces - usually two.
Interior edges may be specified using msh_VolMeshSetEdge()
.
The interior edges may connect to the
exterior boundary faces.
Interior faces which bound material on both sides
may be defined by specifying two “double backed” faces
with opposite connectivity sense.
The connectivity of the
faces is such that the normal to the faces is pointed out of the
volume. The direction of the normal is determined by the right-hand rule
applied to the faces connectivity.
Each bounded region may be assigned a fill type. Currently two fill
types are supported, no fill and fill with tetrahedra. The fill type
is specified using the function msh_VolMeshSetFaceFill()
. If the boundary
contains quadrilateral elements then a pyramid element is generated
to transition conformally from the quadrilateral face to the interior
tetrahedra. It is possible to input mixed order triangle and quadrilateral
elements. For example, the input may contain linear quadrilaterals and
parabolic or cubic and mixed order linear/parabolic-cubic triangles.
The interior
tetrahedron are parabolic or cubic and all pyramids and tetrahedra connecting
to linear edges in the boundary mesh will be of mixed
linear/parabolic-cubic order. Note that the order of the input surface
triangles and quadrilaterals is unchanged by msh_VolMeshGenerate()
, ie.
no nodes will be inserted on the edges of the input surface.
The function msh_VolMeshSetParami()
is used to set integer parameters which
affect the mesh generation process.
The function msh_VolMeshSetParamd()
is used to set real double precision
parameters which control the quality and size of the generated elements.
A global target edge length for interior elements may be prescribed.
Specific, location dependent, sizes may be set at input nodes
using msh_VolMeshSetNodeSizing()
and for regions using
msh_VolMeshSetFaceSizing()
.
By default, all input nodes must be connected to input boundary faces
and interior edges and faces. If unconnected nodes are to be inserted
in the interior then the user must call msh_VolMeshSetParami()
to enable parameter type VOLMESH_UNCONNECT.
The user is able to define integer associations at nodes and faces
which will be assigned to the node and elements generated on the
respective entity.
These associations are useful for identifying nodes and elements
in the output mesh for the application of loads, boundary conditions, material
properties, etc.
The function msh_VolMeshSetNodeAssoc()
is used to set associations at
nodes generated at input points. The function msh_VolMeshSetNodeAssoc()
may be called any number of times for a single node with the same
association type. The association values are accumulated as multiple
association on the node.
The function msh_VolMeshSetFaceAssoc()
is used to
set associations which will be attached to nodes generated on faces and along
face edges. These associations will also
be generated as element entity associations along the respective
element edges and faces connected to the input faces.
The element associations are useful for tagging elements generated
with an enclosed boundary, etc.
Finally the function msh_VolMeshGenerate()
generates the
nodes and elements and enters them into a Connect object.
The mesh generation process may fail due to errors in the surface faces
or failure of the meshing procedure. Use msh_VolMeshGetInteger()
to query for detailed information concerning errors in the surface faces,
etc.
During the mesh generation process, the input surface faces and
face corner node locations (by default) remain unchanged.
However, the midside node locations
may be perturbed under certain circumstances to satisfy quality constraints
in the output elements.
The function msh_VolMeshSetParami()
is
used to set the parameter
VIS_MESH_MIDSIDEQUAL which controls the overall movement of midside
nodes. By default all midside nodes are susceptible to movement.
The function msh_VolMeshWrite()
is provided to write a complete description
of the defined input surface mesh and meshing parameters to a file. The primary
use of this file is to encapsulate cases
in which the VolMesh module should fail in some respect.
This output file can then be
made available to Tech Soft 3D for failure diagnosis of VolMesh
meshing algorithms.
The VolMesh module provides for a monitor
callback function to be specified which
is called intermittently during the mesh generation process in
msh_VolMeshGenerate()
. The primary purpose of this function is to
allow the user to interrupt or interrogate
the on-going mesh generation process for any reason. If mesh generation is to
be terminated, call msh_VolMeshAbort()
. This call will set an internal
flag which will cause msh_VolMeshGenerate()
to terminate the mesh generation
process and return.
In addition, a user defined sizing callback function may be set.
Set callback functions using msh_VolMeshSetFunction()
.
5.6. Function Descriptions
The currently available VolMesh functions are described in detail in this section.
-
msh_VolMesh *msh_VolMeshBegin(void)
create an instance of a VolMesh object
Create an instance of an VolMesh object. Memory is allocated for the object private data and the pointer to the data is returned.
Destroy an instance of a VolMesh object using
void msh_VolMeshEnd (msh_VolMesh *volmesh)
Return the current value of a VolMesh object error flag using
Vint msh_VolMeshError (msh_VolMesh *volmesh)
- Returns:
The function returns a pointer to the newly created VolMesh object. If the object creation fails, NULL is returned.
-
void msh_VolMeshEnd(msh_VolMesh *p)
destroy an instance of a VolMesh object
-
Vint msh_VolMeshError(msh_VolMesh *p)
destroy an instance of a VolMesh object
-
void msh_VolMeshDef(msh_VolMesh *p, Vint numnode, Vint numface)
define number of nodes, faces.
Define the number of nodes and faces which define the input surface mesh. Define node coordinates using
msh_VolMeshSetNode()
and define face connectivities which reference the input nodes usingmsh_VolMeshSetFace()
.Inquire of defined numnode and numface. as output arguments using
void msh_VolMeshInq (msh_VolMesh *volmesh, Vint *numnode, Vint *numface)
- Errors
VIS_ERROR_VALUE
is generated if an improper maxi is input or numnode or numface is less than or equal to zero.
- Parameters:
p – Pointer to VolMesh object.
numnode – Number of nodes on input surface mesh
numface – Number of faces on input surface mesh
-
void msh_VolMeshInq(msh_VolMesh *p, Vint *numnode, Vint *numface)
destroy an instance of a VolMesh object
See
msh_VolMeshDef()
-
void msh_VolMeshSetFunction(msh_VolMesh *p, Vint funtype, Vfunc *function, Vobject *object)
set a call back function
Set callback functions. By default the callback functions are NULL. A callback is not invoked if it is NULL.
The monitor callback function prototype is
The first argument is the VolMesh object, volmesh, and the second is a user defined object, object.void function (msh_VolMesh *volmesh, Vobject *object)
The sizing callback function prototype is
The first argument is the VolMesh object, volmesh, the second is a user defined object, object, the third is the coordinate location x and the fourth is the returned size s.void function (msh_VolMesh *volmesh, Vobject *object, Vdouble x[3], Vdouble *s)
The anisotropic sizing callback function prototype is
The first argument is the VolMesh object, volmesh, the second is a user defined object, object, the third is the coordinate location x and the fourth is the returned anisotropic scaled orthogonal direction vectors. The first 3 components are the first size scaled direction, the next 3 components are the second size scaled direction, the next 3 components are the third size scaled direction.void function (msh_VolMesh *volmesh, Vobject *object, Vdouble x[3], Vdouble s[3][3])
- Parameters:
p – Pointer to VolMesh object.
funtype – Type of callback function to set
x=SYS_FUNCTION_MONITOR Monitor callback =SYS_FUNCTION_SIZING Sizing callback =SYS_FUNCTION_ASIZING Anisotropic sizing callback
function – Pointer to callback function
object – Pointer to the object to be returned as function argument
-
void msh_VolMeshSetParami(msh_VolMesh *p, Vint ptype, Vint iparam)
set mesh generation parameters
Specify mesh generation parameters. These parameters specify the order and connectivity convention to be used for the generated elements.
The parameter
VIS_MESH_MAXI
sets the order of the generated tetrahedral elements by specifying the number of nodes to generate along an edge. A value of 2 generates linear tetrahedra, a value of 3 generates parabolic tetrahedra and a value of 4 generates cubic tetrahedra. By defaultVIS_MESH_MAXI
is set to 2.The parameter
VIS_MESH_MIDSIDEQUAL
toggles the adjustment of midside nodes of parabolic and cubic elements to satisfy specified quality constraints. The midside nodes are adjusted to ensure that all node Jacobians are positive and that the Jacobian ratio in any element is less than the specified upper Jacobian ratio bound. This bound is specified usingmsh_VolMeshSetParamd()
. This procedure normally moves the midside node location of interior edges to satisfy quality constraints. It may move midside nodes on the surface face in rare cases. By defaultVIS_MESH_MIDSIDEQUAL
is set toSYS_ON
.The parameter
VIS_MESH_FUNSIZINGMIN
is used to control the way final element sizing is computed when a sizing callback is installed. By default, the sizing returned by the sizing callback is used. IfVIS_MESH_FUNSIZINGMIN
is enabled, then the sizing is the minimum of the sizing callback and normal sizing. By defaultVIS_MESH_FUNSIZINGMIN
is set toVIS_OFF
.The parameter
VIS_MESH_FUNSIZINGLEGACY
is used to revert the sizing callback behavior to its legacy implementation. In this case the sizing is only applied to the nodes, and not to elements nor edges. By defaultVIS_MESH_FUNSIZINGLEGACY
is set toVIS_OFF
.By default, all input nodes must be connected to input boundary faces and interior edges and faces. If unconnected nodes are to be inserted in the interior then parameter
VOLMESH_UNCONNECT
must be enabled. By defaultVOLMESH_UNCONNECT
is set toSYS_OFF
.Inquire of defined integer parameters as output arguments using
void msh_VolMeshGetParami (msh_VolMesh *volmesh, Vint ptype, Vint *iparam)
- Errors
VIS_ERROR_ENUM
is generated if an improper ptype is specified.VIS_ERROR_VALUE
is generated if an improper order of generated elements is specified.
- Parameters:
p – Pointer to VolMesh object.
ptype – Type of display parameter to set
x=VIS_MESH_MAXI Order of generated elements =VIS_MESH_MIDSIDEQUAL Quality midside node adjustment =VIS_MESH_REPAIRINCONS Repair inconsistent input =VIS_MESH_FUNSIZINGMIN Size is min'ed with sizing callback =VIS_MESH_FUNSIZINGLEGACY Legacy sizing callback =VOLMESH_UNCONNECT Toggle unconnected interior nodes
iparam – Specifies the integer value that ptype will be set to.
x=SYS_OFF Disable =SYS_ON Enable
-
void msh_VolMeshGetParami(msh_VolMesh *p, Vint type, Vint *iparam)
get mesh generation parameters
-
void msh_VolMeshSetParamd(msh_VolMesh *p, Vint ptype, Vdouble dparam)
set mesh generation parameters
Specify mesh generation parameters. These parameters are used to control the quality and size of the generated elements.
The parameter
VIS_MESH_EDGELENGTH
specifies a target edge length of generated elements. The edge length constraint may not be satisfied for elements generated near the surface mesh due to surface mesh sizing. By defaultVIS_MESH_EDGELENGTH
is set to 1.The parameter
VIS_MESH_GROWTHRATE
specifies a maximum growth rate. The growth rate governs the rate at which the size of adjacent elements can grow. The ratio of edge lengths of adjacent elements will not exceed the specified growth rate. The growth rate must be greateer than 1.001. By defaultVIS_MESH_GROWTHRATE
is set to 2.0The parameter
VIS_MESH_MAXEDGELENGTH
specifies a maximum edge length of generated elements. The maximum edge length constraint may not be satisfied for elements generated near the surface mesh due to surface mesh sizing. By defaultVIS_MESH_MAXEDGELENGTH
is set to 0. and is ignored.The parameter
VIS_MESH_MAXEDGEALT
specifies an upper bound on the edge-altitude ratio of an element. The edge-altitude ratio is defined as the maximum edge length over the smallest altitude. A series of operations are performed in an attempt to satisfy the edge-altitude ratio quality. By defaultVIS_MESH_MAXEDGEALT
is set to 100.Inquire of defined double precision parameters as output arguments using
void msh_VolMeshGetParamd (msh_VolMesh *volmesh, Vint ptype, Vdouble *dparam)
- Errors
VIS_ERROR_ENUM
is generated if an improper ptype is specified.VIS_ERROR_VALUE
is generated if an improper value is specified.
- Parameters:
p – Pointer to VolMesh object.
ptype – Type of display parameter to set
dparam – Specifies the value that ptype will be set to.
x=VIS_MESH_EDGELENGTH Target edge length =VIS_MESH_GROWTHRATE Maximum growth rate =VIS_MESH_MAXEDGELENGTH Maximum edge length =VIS_MESH_MAXEDGEALT Maximum edge-altitude ratio =VIS_MESH_MAXANGLE Maximum dihedral angle
-
void msh_VolMeshGetParamd(msh_VolMesh *p, Vint type, Vdouble *param)
get mesh generation parameters
-
void msh_VolMeshAbort(msh_VolMesh *p)
set abort flag
Set mesh generation abort flag. During the mesh generation process in
msh_VolMeshGenerate()
this flag is intermittently checked. If it has been set, the mesh generation process terminates and returns. The abort flag is usually set in the user defined callback function specified bymsh_VolMeshSetFunction()
.- Parameters:
p – Pointer to VolMesh object.
-
void msh_VolMeshSetNode(msh_VolMesh *p, Vint id, Vdouble x[3])
define a surface node
Set a node location. The node id must be in the interval 1 <= id <= numnode where numnode is defined by
msh_VolMeshDef()
. Input surface faces are defined usingmsh_VolMeshSetFace()
.- Errors
VIS_ERROR_VALUE
is generated if an improper id is specified.
- Parameters:
p – Pointer to VolMesh object.
id – Node id
x – Node coordinate location
-
void msh_VolMeshSetNodeStat(msh_VolMesh *p, Vint id, Vint stat)
set surface node status
Set a surface corner or midside node movement status. By default all corner nodes may not be moved and all midside nodes may be moved. Set the status to 1 for any surface corner or midsides nodes which must not be moved for any reason. The default movement status of surface corner and midside nodes may be specified using
msh_VolMeshSetParami()
.- Errors
VIS_ERROR_VALUE
is generated if an improper id is specified.
- Parameters:
p – Pointer to VolMesh object.
id – Node id
stat – Surface node status
=0 Surface node may be moved =1 Surface node may not be moved
-
void msh_VolMeshSetNodeSizing(msh_VolMesh *p, Vint id, Vdouble h)
set node sizing
Set element size to be enforced at a node.
- Errors
VIS_ERROR_VALUE
is generated if an improper id is specified.
- Parameters:
p – Pointer to VolMesh object.
id – Node id
h – Sizing.
-
void msh_VolMeshSetNodeAssoc(msh_VolMesh *p, Vint type, Vint id, Vint aid)
define a node association
Set a node association. The association will be propagated to the node in the generated mesh. This function may be called any number of times for a particular node and association type, all association values are accumulated.
- Errors
VIS_ERROR_VALUE
is generated if an improper id is specified.
- Parameters:
p – Pointer to VolMesh object.
type – Type of association
id – Node id
aid – Association value
-
void msh_VolMeshSetFace(msh_VolMesh *p, Vint id, Vint shape, Vint maxi, Vint maxj, Vint ix[])
define a face
Set a surface face topology and connectivity. The face id must be in the interval 1 <= id <= numface where numface is defined by
msh_VolMeshDef()
. The connectivity of the faces is be such that the normal to the faces is pointed out of the volume. The direction of the normal is determined by the right-hand rule applied to the faces connectivity. Note that all faces must be entered sequentially.- Errors
VIS_ERROR_VALUE
is generated if an improper id or order is specified.VIS_ERROR_ENUM
is generated if an improper shape is specified.VIS_ERROR_OPERATION
is generated if an id is entered out of order.
- Parameters:
p – Pointer to VolMesh object.
id – Face id
shape – Face shape
maxi – Face order
maxj – Face order
ix – Face connectivity
-
void msh_VolMeshSetFaceAssoc(msh_VolMesh *p, Vint type, Vint id, Vint enttype, Vint no, Vint aid)
define a face association
Set associations to be assigned to generated elements, edges and faces. The entity type of
SYS_EDGE
will propagate the association to element edges and all nodes on the face edge. The entity type ofSYS_FACE
will propagate the association to element faces and all nodes on the face, the value of no is ignored. The entity type ofSYS_ELEM
will propagate the association to all elements adjacent to the specified face, bounded by adjacent input faces, the value of no is ignored.- Errors
VIS_ERROR_ENUM
is generated if an improper type or enttype is specified.VIS_ERROR_VALUE
is generated if an improper id or no is specified.
- Parameters:
p – Pointer to VolMesh object.
type – Type of association
id – Face id
enttype – Entity type to which association is assigned
=SYS_EDGE Element edge association =SYS_FACE Element face association =SYS_ELEM Element association
no – Face edge number
aid – Association value
-
void msh_VolMeshSetFaceStat(msh_VolMesh *p, Vint id, Vint no, Vint stat)
set face status
Set face status. By default face edges are not swapped or split during generation and may be split but not swapped during refinement.
- Errors
VIS_ERROR_VALUE
is generated if an improper id or no is specified.
- Parameters:
p – Pointer to VolMesh object.
id – Face id
no – Face edge number
stat – Status, bitwise OR
x=VIS_MESH_BIT_GENSPLIT Split during generation =VIS_MESH_BIT_GENSWAP Swap during generation =VIS_MESH_BIT_REFSPLIT Split during refinement =VIS_MESH_BIT_REFSWAP Swap during refinement
-
void msh_VolMeshSetFaceSizing(msh_VolMesh *p, Vint id, Vdouble h)
set sizing within a region
Set element size to be enforced within a region. The region is defined as the material bound by input faces that exists to the left of the specified face, ie. the normal of the face points out of the material region. If a size is specified on more than one face bounding a region, the minimum size is used for the region.
- Errors
VIS_ERROR_VALUE
is generated if an improper id is specified.
- Parameters:
p – Pointer to VolMesh object.
id – Face id
h – Sizing.
-
void msh_VolMeshSetFaceDepthSizing(msh_VolMesh *p, Vint id, Vdouble d, Vdouble h)
set sizing within a depth in a region
Set element size to be enforced within a specified depth from the input boundary face.
- Errors
VIS_ERROR_VALUE
is generated if an improper id is specified.VIS_ERROR_VALUE
is generated if either d or h is less than zero.
- Parameters:
p – Pointer to VolMesh object.
id – Face id
d – Depth
h – Sizing
-
void msh_VolMeshSetFaceFill(msh_VolMesh *p, Vint id, Vint fill)
define a bounding face fill type
Set type of element fill to be generated adjacent to a given boundary face. A fill type of
VOLMESH_FILL_TET
will fill the region adjacent to the face with tetrahedra. Any quadrilateral input faces will have pyramids created adjacent to the face to transition to the tetrahedral mesh. A fill type ofVOLMESH_FILL_BOUNDARY
will output the boundary face directly as a generated element. No volumetric elements will be generated in the region adjacent to the boundary face.- Errors
VIS_ERROR_VALUE
is generated if an improper id is specified.VIS_ERROR_ENUM
is generated if an improper fill is specified.
- Parameters:
p – Pointer to VolMesh object.
id – Face id
fill – Fill type
x=VOLMESH_FILL_TET Fill directly with tetrahedra =VOLMESH_FILL_BOUNDARY Do not volumetrically fill
-
void msh_VolMeshSetEdge(msh_VolMesh *p, Vint id, Vint maxi, Vint ix[])
set interior edge
Set an interior edge order, maxi and connectivity, ix. The nodes of an interior edge may be connected to the boundary.
- Errors
VIS_ERROR_VALUE
is generated if an improper id, maxi or ix is specified.
- Parameters:
p – Pointer to VolMesh object.
id – Edge id
maxi – Number of edge nodes
ix – Edge connectivity
-
void msh_VolMeshSetEdgeStat(msh_VolMesh *p, Vint id, Vint stat)
set interior edge status
Set an interior edge status. By default interior edges are not swapped or split during generation and may be split but not swapped during refinement.
- Errors
VIS_ERROR_VALUE
is generated if an improper id is specified.
- Parameters:
p – Pointer to VolMesh object.
id – Edge id
stat – Status, bitwise OR
x=VIS_MESH_BIT_GENSPLIT Split during generation =VIS_MESH_BIT_GENSWAP Swap during generation =VIS_MESH_BIT_REFSPLIT Split during refinement =VIS_MESH_BIT_REFSWAP Swap during refinement
-
void msh_VolMeshSetGeomSizing(msh_VolMesh *p, Vint type, Vdouble xo[3], Vdouble xa[3], Vdouble xb[3], Vdouble d[3], Vdouble h)
specify sizing within geometric shape
Set element size to be enforced within a geometric shape.
- Errors
VIS_ERROR_ENUM
is generated if an improper type is specified.
- Parameters:
p – Pointer to VolMesh object.
type – Geometry type
x=VIS_MESH_GEOMSPHERE Ellipsoid =VIS_MESH_GEOMCYLINDER Ellipsoidal cylinder =VIS_MESH_GEOMBOX Box
xo – Coordinates of geometry origin
xa – Coordinates of point along x’ axis
xb – Coordinates of point in x’y’ plane
d – Half length of geometric shape along each axis
h – Sizing
-
void msh_VolMeshGetConnect(msh_VolMesh *p, vis_Connect **connect)
get internal Connect object
Get the internal Connect object which contains the user input nodes and faces. The Connect object must not be changed by the user.
-
void msh_VolMeshWrite(msh_VolMesh *p, Vint type, const Vchar *path)
write mesh description to file
Write a surface mesh description to a file which has been input to VolMesh. The contents of the file is intended for internal use only. This function should be called after all surface mesh input and meshing parameters have been set to ensure that an accurate surface mesh description including correct meshing parameters is written. The binary file format is preferred for accuracy and generally results in smaller file sizes. The primary use of this function is to create a file of the complete surface mesh description in the case that the VolMesh module fails to properly generate a mesh. This file can then be transmitted to Visual Kinematics for debugging purposes. It is suggested that a file extension of .vol be used for ASCII files and .bvol be used for BINARY files.
- Errors
SYS_ERROR_ENUM
is generated if an improper type is specified.SYS_ERROR_FILE
is generated if the file can not be opened.SYS_ERROR_OPERATION
is generated if the input surface mesh has not been defined.
- Parameters:
p – Pointer to VolMesh object.
type – File type
x=SYS_ASCII ASCII format =SYS_BINARY Binary format
path – File path
-
void msh_VolMeshRead(msh_VolMesh *p, Vint type, const Vchar *path)
read mesh description from file
Read a surface mesh description from a file which has been previously written by
msh_VolMeshWrite()
. The format of the file is not published and is intended for internal use only.- Errors
SYS_ERROR_ENUM
is generated if an improper type is specified.SYS_ERROR_FILE
is generated if the file can not be opened.SYS_ERROR_FORMAT
is generated if the file contents are not correctly formatted or if a file written by a later version of VolMesh is read.
- Parameters:
p – Pointer to VolMesh object.
type – File type
x=SYS_ASCII ASCII format =SYS_BINARY Binary format
path – File path
-
void msh_VolMeshConnect(msh_VolMesh *p, vis_Connect *connect)
set the contents of Connect object
Set the contents of the Connect object as geometry input. All nodes and elements in the Connect object are set as input. All
VIS_GEOVERT
,VIS_GEOEDGE
andVIS_GEOFACE
associations on nodes are set. AllVIS_GEOEDGE
andVIS_GEOFACE
associations on faces are set as element entity associations. AllVIS_GEOBODY
associations are set as element associations.- Errors
VIS_ERROR_OPERATION
is generated if the Connect object does not contain proper element types to serve as boundary input.
-
void msh_VolMeshGenerate(msh_VolMesh *p, vis_Connect *connect)
generate finite element mesh
Generate interior nodes and elements. The generated nodes and elements are appended to any existing nodes and elements in the input Connect object. The input surface faces are checked for free edges, connectivity consistency, degeneracy, enclosing negative or zero volume, etc. Use
msh_VolMeshGetInteger()
is query for detailed information about the mesh generation process including errors.- Errors
VIS_ERROR_OPERATION
is generated if mesh generation fails due to errors detected in the surface faces.VIS_ERROR_COMPUTE
is generated if mesh generation fails.
-
void msh_VolMeshGetInteger(msh_VolMesh *p, Vint type, Vint iparam[])
get integer mesh generation information
Query for integer mesh generation information. This function is useful for determining where errors in the mesh generation procedure occurred during a call to
msh_VolMeshGenerate()
. UseVOLMESH_NUMFREEEDGE
to determine the number of faces which contained a free edge. The queryVOLMESH_FREEEDGE
will return the list of all free edge faces detected.The query type
VOLMESH_NUMNONMAN
returns the number of improper non-manifold boundary face detected in the input. This number will only include those face detected on the first non-manifold edge encountered, it will not generally include all non-manifold face. The queryVOLMESH_NONMAN
will return the list of non-manifold face detected.The query type
VOLMESH_NUMINCONS
returns the number of inconsistent faces detected in the input. This number will only include those faces detected on the first inconsistent edge encountered, it will not generally include all inconsistent faces. The queryVOLMESH_INCONS
will return the list of inconsistent faces detected.The query type
VOLMESH_NUMNEGJAC
returns the number of elements with negative corner Jacobians. This check is only performed if it is enabled usingmsh_VolMeshSetParami()
.The query type
VOLMESH_NUMUNCONN
returns the number of unconnected nodes detected in the input. The queryVOLMESH_UNCONN
will return the list of unconnected nodes.The query type
VOLMESH_NUMTRIZEROANG
returns the number of input faces with zero angle detected in the input. The queryVOLMESH_TRIZEROANG
will return the list of zero angle faces.The query type
VOLMESH_NUMINTERSECT
returns the number of input faces which intersect. The queryVOLMESH_INTERSECT
will return the list of intersecting faces.The query type
VOLMESH_NUMUNREC
returns the number of input faces which were unable to be recovered. The queryVOLMESH_UNREC
will return the list of unrecovered faces.The query type
VOLMESH_NUMNODEINSERT
returns the number of input nodes which were unable to be inserted. The queryVOLMESH_NODEINSERT
will return the list of nodes which were not inserted.The query type
VIS_MESH_PROGRESS
returns information about the current state of the meshing process. This function is most useful when called from the monitor function. The following four integer values are returned.iparam[0] Meshing phase Phase 1 - Pre-processing boundary data Phase 2 - Boundary recovery Phase 3 - Generating interior points Phase 4 - Mesh cleanup Phase 5 - Mesh quality improvement iparam[1] Current number of nodes iparam[2] Current number of elements iparam[3] Estimate percent completion
- Errors
VIS_ERROR_ENUM
is generated if an improper type is specified.
- Parameters:
p – Pointer to VolMesh object.
type – Type of integer information to query
x=VOLMESH_NUMFREEEDGE Number of surface face free edges =VOLMESH_FREEEDGE Free edge surface face =VOLMESH_NUMINCONS Number of inconsistent face =VOLMESH_INCONS Inconsistent surface face =VOLMESH_NUMINTERSECT Number of intersecting face =VOLMESH_INTERSECT Intersecting face =VOLMESH_NUMNONMAN Number of non manifold face =VOLMESH_NONMAN Non manifold surface face =VOLMESH_NUMNEGJAC Number of negative Jacobians =VOLMESH_NUMUNCONN Number of unconnected nodes =VOLMESH_UNCONN Unconnected nodes =VOLMESH_NUMTRIZEROANG Number of zero angle face =VOLMESH_TRIZEROANG Zero angle face =VOLMESH_NUMUNREC Number of unrecovered face =VOLMESH_UNREC Unrecovered face =VOLMESH_NUMNODEINSERT Number of nodes not inserted =VOLMESH_NODEINSERT Not inserted nodes =VIS_MESH_PROGRESS Current meshing progress
iparam – [out] Returned integer information