3. Triangular Mesh Generation - TriMesh
The TriMesh module provides for the automatic generation of linear, parabolic or cubic triangular elements given a boundary mesh of linear, parabolic or cubic line segments. By default, the mesh generation procedure is guaranteed to not alter the boundary mesh in any way. The node numbering of the boundary mesh is preserved in the node numbering and connectivity of the resultant triangle mesh. The input boundary mesh is checked for closure.
The functions associated with a TriMesh object are the following.
Begin and end an instance of an object, return object error flag
msh_TriMeshBegin()
- create an instance of a TriMesh objectmsh_TriMeshEnd()
- destroy an instance of a TriMesh objectmsh_TriMeshError()
- return TriMesh object error flag
Define boundary mesh, mesh parameters and generate mesh
msh_TriMeshAbort()
- set abort flagmsh_TriMeshDef()
- define number of boundary nodes, lines.msh_TriMeshInq()
- inquire number of boundary nodes, lines.msh_TriMeshComputeArea()
- compute areamsh_TriMeshGenerate()
- generate finite element meshmsh_TriMeshGetConnect()
- get internal Connect objectmsh_TriMeshGetInteger()
- get integer mesh generation informationmsh_TriMeshRefine()
- refine a triangle meshmsh_TriMeshSetElem()
- set an orphan mesh for refinementmsh_TriMeshSetGeomSizing()
- specify sizing within geometric shapemsh_TriMeshSetLine()
- define a boundary linemsh_TriMeshSetLineAssoc()
- define a line associationmsh_TriMeshSetLineSizing()
- set line and region sizingmsh_TriMeshSetLineDepthSizing()
- set sizing within a depth in a regionmsh_TriMeshSetLineStat()
- set line statusmsh_TriMeshSetFunction()
- set a call back functionmsh_TriMeshSetParami()
- set mesh generation parametersmsh_TriMeshGetParami()
- get mesh generation parametersmsh_TriMeshSetParamd()
- set mesh generation parametersmsh_TriMeshGetParamd()
- get mesh generation parametersmsh_TriMeshSetNode()
- define a nodemsh_TriMeshSetNodeAssoc()
- define a node associationmsh_TriMeshSetNodeSizing()
- set node sizingmsh_TriMeshWrite()
- write input description to filemsh_TriMeshRead()
- read input description from file
Instance a TriMesh object using msh_TriMeshBegin()
.
Once a TriMesh object is instanced,
define the number of nodes and lines which define the boundary mesh
as well as the order of the lines (linear, parabolic or cubic)
using msh_TriMeshDef()
.
The boundary mesh is then input as
a set of nodes in space and a set of line
connectivities.
The nodes are defined using msh_TriMeshSetNode()
,
the line connectivities are defined using msh_TriMeshSetLine()
.
By default, the connectivity of the
lines is such that the normal to the line is pointed out of the area.
The input segments must form closed loops.
Loops will either mark the boundaries
of area or the boundaries of holes depending upon the connectivities.
A special case exists for loops consisting of two “back-to-back” segments.
These special loops are considered to be line constraints.
The node connectivity and normal conventions for input boundary lines is
illustrated in Figure 26-1. For linear input lines only two nodes are input,
for parabolic lines three nodes are input and for cubic lines four nodes
are input.
The function msh_TriMeshSetNodeAssoc()
is used to set associations at nodes.
The function msh_TriMeshSetLineAssoc()
is used to set associations at
input lines which, depending upon the type of association which is set,
will generate associations at generated triangles and triangle
edges. These associations are useful for tagging triangles generated
with an enclosed boundary, etc.
The function msh_TriMeshSetParami()
is used to set integer parameters which
affect the mesh generation process.
By default, all input nodes must be connected to input boundary lines.
If unconnected nodes are to be inserted in the interior
then the user must call msh_TriMeshSetParami()
to enable parameter
type TRIMESH_UNCONNECT.
It is possible to specify
parabolic or cubic triangles to be produced if the input lines are linear.
By default the output triangles are the same order as the input lines.
Also the sense of the normals of the input lines may be changed so that
normals are expected to point into the area.
An option exists to generate a mesh with no obtuse triangles.
The function msh_TriMeshSetParamd()
is used to set real double precision
parameters which control the quality and size of the generated triangles.
A global target edge length for interior elements may be prescribed.
Specific, location dependent, sizes may be set at input nodes
using msh_TriMeshSetNodeSizing()
and at input lines and regions bound
by input lines using msh_TriMeshSetLineSizing()
.
The specified growth rate is used
to regulate the element size transition rate between the location dependent
sizing and the overall target edge length.
Note that the element size at any
point in the mesh is taken to be the minimum of all user sizing input.
Finally the function msh_TriMeshGenerate()
generates the
nodes and finite elements and enters them into a Connect object.
The mesh generation process may fail due to errors in the boundary lines
or failure of the triangulation procedure. Use msh_TriMeshGetInteger()
to query for detailed information concerning errors in the boundary lines,
etc.
An existing mesh may be refined using msh_TriMeshRefine()
.
Mesh refinement is performed on the mesh resulting from the
call to msh_TriMeshGenerate()
or the last previous call to msh_TriMeshRefine()
.
A State object of element edge length
is input and the mesh is refined to satisfy the new element edge lengths.
Only mesh refinement is performed, all previously existing nodes are
retained in the refined mesh. By default the previous node locations
are unchanged.
If an element size is set to zero, the element will not be refined.
By default, all input boundary lines may be subdivided during mesh
generation or refinement. Use msh_TriMeshSetLineStat()
to control
subdivision of input lines.
Mesh refinement is supported both for meshes generated
with msh_TriMeshGenerate()
or for orphan meshes.
In the case of orphan meshes, in addition to specifying
the line discretized boundary and its associations, the orphan mesh is
specified by entering each original element with msh_TriMeshSetElem()
.
The function msh_TriMeshWrite()
is provided to write a complete description
of the defined input meshing parameters to a file. The primary
use of this file is to encapsulate cases
in which the TriMesh module fails
in some respect. This output file can then be
made available to Tech Soft 3D for failure diagnosis of TriMesh
meshing algorithms.
The TriMesh module provides for a callback function to be specified which
is called intermittently during the mesh generation process in
msh_TriMeshGenerate()
or msh_TriMeshRefine()
.
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_TriMeshAbort()
. This call will set an internal
flag which will cause msh_TriMeshGenerate()
or msh_TriMeshRefine()
to terminate the mesh generation process and return.
Set callback functions using msh_TriMeshSetFunction()
.
The following code fragment illustrates the basic framework of using the TriMesh module to generate linear triangular finite elements in an area bounded by a set of straight line segments.
/* declare objects */
msh_TriMesh *trimesh;
vis_Connect *connect;
/* input boundary lines */
Vint numnode, numline;
Vdouble xnod[][3] = { ... };
Vint ixlin[][2] = { ... };
Vdouble edgelength;
Vint ierr;
/* create trimesh object */
trimesh = msh_TriMeshBegin ();
msh_TriMeshDef (trimesh,numnode,numline,2);
/* define nodes and lines */
for (i = 0; i < numnode; i++) {
msh_TriMeshSetNode (trimesh,i+1,xnod[i]);
}
for (i = 0; i < numline; i++) {
msh_TriMeshSetLine (trimesh,i+1,ixlin[i]);
}
/* set meshing parameters */
msh_TriMeshSetParamd (trimesh,VIS_MESH_EDGELENGTH,edgelength);
/* create empty connect object to hold generated mesh */
connect = vis_ConnectBegin ();
vis_ConnectPre (connect,SYS_DOUBLE);
vis_ConnectDef (connect,0,0);
/* generate mesh */
msh_TriMeshGenerate (trimesh,connect);
/* process any errors */
ierr = msh_TriMeshError (trimesh);
if(ierr) {
...
}
/* access generated nodes and elements from Connect */
...
/* delete objects */
vis_ConnectEnd (connect);
msh_TriMeshEnd (trimesh);
3.1. Function Descriptions
The currently available TriMesh functions are described in detail in this section.
-
msh_TriMesh *msh_TriMeshBegin(void)
create an instance of a TriMesh object
Create an instance of a TriMesh object. Memory is allocated for the object private data and the pointer to the data is returned.
Destroy an instance of a TriMesh object using
void msh_TriMeshEnd (msh_TriMesh *trimesh)
Return the current value of a TriMesh object error flag using
Vint msh_TriMeshError (msh_TriMesh *trimesh)
- Returns:
The function returns a pointer to the newly created TriMesh object. If the object creation fails, NULL is returned.
-
void msh_TriMeshEnd(msh_TriMesh *p)
destroy an instance of a TriMesh object
-
Vint msh_TriMeshError(msh_TriMesh *p)
return the current value of a TriMesh object error flag
-
void msh_TriMeshDef(msh_TriMesh *p, Vint numnode, Vint numline, Vint maxi)
define number of boundary nodes, lines.
Define the number of nodes and lines which define the input boundary mesh. The number of nodes on a boundary line segment is specified by maxi . Define node coordinates using
msh_TriMeshSetNode()
and define line connectivities which reference the input nodes usingmsh_TriMeshSetLine()
.Inquire of defined numnode, numline and maxi as output arguments using
void msh_TriMeshInq (msh_TriMesh *trimesh, Vint *numnode, Vint *numline, Vint *maxi)
- Errors
VIS_ERROR_VALUE
is generated if an improper maxi is input or numnode is less than or equal to zero or numline is less than zero.
- Parameters:
p – Pointer to TriMesh object.
numnode – Number of nodes on input boundary mesh
numline – Number of lines on input boundary mesh
maxi – Number of nodes on line segment
=2 Linear boundary line segments =3 Parabolic boundary line segments =4 Cubic boundary line segments
-
void msh_TriMeshInq(msh_TriMesh *p, Vint *numnode, Vint *numline, Vint *maxi)
inquire of defined numnode, numline and maxi as output arguments
See
msh_TriMeshDef()
-
void msh_TriMeshSetFunction(msh_TriMesh *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. If an anisotropic sizing function is set, it takes precedence over any isotropic sizing function which has been set.
The monitor callback function prototype is
The first argument is the TriMesh object, p, and the second is a user defined object, object.void function (msh_TriMesh *trimesh, Vobject *object)
The isotropic sizing callback function prototype is
The first argument is the TriMesh object, trimesh, 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_TriMesh *trimesh, Vobject *object, Vdouble x[3], Vdouble *s)
The anisotropic sizing callback function prototype is
The first argument is the TriMesh object, trimesh, 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_TriMesh *trimesh, Vobject *object, Vdouble x[3], Vdouble s[3][3])
- Parameters:
p – Pointer to TriMesh object.
funtype – Type of callback function to set
x=SYS_FUNCTION_MONITOR Monitor callback =SYS_FUNCTION_SIZING Isotropic 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_TriMeshSetNode(msh_TriMesh *p, Vint id, Vdouble x[3])
define a node
Set a node location. The node id must be in the interval 1 <= id <= numnode where numnode is defined by
msh_TriMeshDef()
. The node location is projected to the x,y plane, the x[2] entry is ignored. Input boundary lines are defined usingmsh_TriMeshSetLine()
.- Errors
VIS_ERROR_VALUE
is generated if an improper id is specified.
- Parameters:
p – Pointer to TriMesh object.
id – Node id
x – Node coordinate location
-
void msh_TriMeshSetNodeAssoc(msh_TriMesh *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 TriMesh object.
type – Type of association
id – Node id
aid – Association value
-
void msh_TriMeshSetNodeSizing(msh_TriMesh *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 TriMesh object.
id – Node id
h – Sizing.
-
void msh_TriMeshSetLine(msh_TriMesh *p, Vint id, Vint ix[])
define a boundary line
Set a boundary line connectivity. The line id must be in the interval 1 <= id <= numline where numline is defined by
msh_TriMeshDef()
. The line connectivity, ix, must contain 2 nodes for linear boundary lines, 3 nodes for parabolic boundary lines and 4 nodes for cubic lines. By default, the connectivity of the lines is be such that the normal to the lines is pointed out of the area. The direction of the normal is determined by the right-hand rule applied to the line connectivity. The outward normal convention may be changed to an inward normal convention usingmsh_TriMeshSetParami()
. Note that all lines must be entered sequentially.- Errors
VIS_ERROR_VALUE
is generated if an improper id or ix is specified.VIS_ERROR_OPERATION
is generated if an id is entered out of order.
- Parameters:
p – Pointer to TriMesh object.
id – Line id
ix – Line connectivity
-
void msh_TriMeshSetLineStat(msh_TriMesh *p, Vint id, Vint stat)
set line status
Set line status. By default, boundary lines are split during generation and refinement.
- Errors
VIS_ERROR_VALUE
is generated if an improper id is specified.
- Parameters:
p – Pointer to TriMesh object.
id – Line id
stat – Status, bitwise OR
x=VIS_MESH_BIT_GENSPLIT Split during generation =VIS_MESH_BIT_REFSPLIT Split during refinement
-
void msh_TriMeshSetLineSizing(msh_TriMesh *p, Vint id, Vint enttype, Vdouble h)
set line and region sizing
Set element sizing to be enforced along an input line or a region bounded by an input line. if enttype is
SYS_EDGE
the specified size is applied to the input line. if enttype isSYS_ELEM
the specified size is enforced within a region. The region is defined as the material bound by input lines that exists to the right of the specified line. If a size is specified on more than one line bounding a region, the minimum size is used for the region.- Errors
VIS_ERROR_VALUE
is generated if an improper id is specified.VIS_ERROR_ENUM
is generated if an improper enttype is specified.
- Parameters:
p – Pointer to TriMesh object.
id – Line id
enttype – Entity type to which association is assigned
=SYS_EDGE Edge sizing =SYS_ELEM Region sizing
h – Sizing.
-
void msh_TriMeshSetLineDepthSizing(msh_TriMesh *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 line.
- 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 TriMesh object.
id – Line id
d – Depth
h – Sizing
-
void msh_TriMeshSetLineAssoc(msh_TriMesh *p, Vint type, Vint id, Vint enttype, Vint aid)
define a line association
Set associations to be assigned to generated triangles and triangle edges. The entity type of
SYS_EDGE
will propagate the association to triangle edges. The entity type ofSYS_ELEM
will propagate the association to all triangles adjacent to specified the line, bounded by adjacent input lines.- Errors
VIS_ERROR_ENUM
is generated if an improper type or enttype is specified.VIS_ERROR_VALUE
is generated if an improper id is specified.
- Parameters:
p – Pointer to TriMesh object.
type – Type of association
id – Line id
enttype – Entity type to which association is assigned
=SYS_EDGE Element edge association =SYS_ELEM Element association
aid – Association value
-
void msh_TriMeshSetElem(msh_TriMesh *p, Vint id, Vint shape, Vint maxi, Vint ix[])
set an orphan mesh for refinement
Input elements for an orphan triangle mesh to be refined. The API has been defined in a general way to support other element shapes. However, only triangle elements, linear, parabolic or cubic, are currently supported.
- Errors
VIS_ERROR_OPERATION
is generated if TriMesh was previously used to generate a mesh.msh_TriMeshSetElem()
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_SHAPETRI
.VIS_ERROR_VALUE
is generated if maxi is not equal to 0, 2, 3 or 4.
- Parameters:
p – Pointer to TriMesh object.
id – Element index
shape – Element shape. Only SYS_SHAPETRI is supported.
maxi – Order of the elements in the initial mesh.
ix – Element connectivity
-
void msh_TriMeshSetParami(msh_TriMesh *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 triangles.
The parameter
VIS_MESH_MAXI
sets the order of the generated triangular elements by specifying the number of nodes to generate along an edge. A value of 2 generates linear triangles, a value of 3 generates parabolic triangles and a value of 4 generates cubic triangles. By defaultVIS_MESH_MAXI
is set to the order of the input boundary lines defined bymsh_TriMeshDef()
.The parameter
VIS_MESH_INWARD
sets the normal convention of the input lines to point into the area. By defaultVIS_MESH_INWARD
is set toSYS_OFF
, that is, the line normals are expected to point out of the area.The parameter
VIS_MESH_BOUNDREFINE
enables the refinement of the boundary. By defaultVIS_MESH_BOUNDREFINE
is set toSYS_ON
.The parameter
VIS_MESH_INTERREFINE
enables the refinement of the interior. By defaultVIS_MESH_INTERREFINE
is set toSYS_ON
.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_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_NONOBTUSE
requests that the final mesh contain no obtuse triangles. By defaultVIS_MESH_NONOBTUSE
is set toSYS_OFF
.The parameter
VIS_MESH_POINTSET
requests a Delaunay triangulation 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
TRIMESH_SMOOTH
enables node location smoothing duringmsh_TriMeshGenerate()
. By defaultTRIMESH_SMOOTH
is set toSYS_OFF
.The parameter
VIS_MESH_REFINESMOOTH
enables node location smoothing duringmsh_TriMeshRefine()
. By defaultVIS_MESH_REFINESMOOTH
is set toSYS_OFF
.The parameter
TRIMESH_BOUNDCUT
enables the generation of rectilinear boundary lines in the interior. The generated boundary lines will connect input boundary lines which are nearly parallel to one another. The generated boundary lines will always connect an input node to a point on an opposing parallel boundary line. The point on the opposing boundary line may be an existing node or a generated point. The deviation from parallel in degrees is specified usingmsh_TriMeshSetParamd()
with parameterTRIMESH_CUTANGLE
. Use the parameterTRIMESH_CUTMINEDGELENGTH
to inhibit cuts which produce small edges. By defaultTRIMESH_BOUNDCUT
is set toSYS_OFF
.The parameters
TRIMESH_CUTFLAG
andTRIMESH_LINEFLAG
enable generation of specified element associations which flag the edges of each triangle which lie on a rectilinear boundary cut or boundary line respectively. The value of these parameters is an association type and must be zero or one ofVIS_MISCID0
throughVIS_MISCID8
. The value of each parameter may be the same association value, in which case the association will flag both hard edges and boundary lines. A value of zero disables the association generation. By defaultTRIMESH_CUTFLAG
andTRIMESH_LINEFLAG
are set to zero.By default, all input nodes must be connected to input lines. If unconnected nodes are to be inserted in the interior then parameter
TRIMESH_UNCONNECT
must be enabled. By defaultTRIMESH_UNCONNECT
is set to SYS_OFF.Inquire of defined integer parameters as output arguments using
void msh_TriMeshGetParami (msh_TriMesh *trimesh, 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 TriMesh object.
ptype – Type of display parameter to set
x=VIS_MESH_BOUNDREFINE Boundary refinement flag =VIS_MESH_DELAUNAY Delaunay mesh flag =VIS_MESH_FUNSIZINGMIN Size is min'ed with sizing callback =VIS_MESH_FUNSIZINGLEGACY Legacy sizing callback =VIS_MESH_INWARD Inward normal line flag =VIS_MESH_INTERREFINE Interior refinement flag =VIS_MESH_MAXI Order of generated elements =VIS_MESH_NONOBTUSE Nonobtuse mesh flag =VIS_MESH_POINTSET Delaunay triangulation of point set =TRIMESH_SMOOTH Smooth generated node locations =VIS_MESH_REFINESMOOTH Smooth refined node locations =TRIMESH_BOUNDCUT Toggle boundary cutting =TRIMESH_CUTFLAG Boundary cut flag associations =TRIMESH_LINEFLAG Element line flag associations =TRIMESH_UNCONNDECT Toggle unconnected interior nodes
iparam – Specifies the integer value that ptype will be set to.
x=SYS_OFF Disable =SYS_ON Enable
-
void msh_TriMeshGetParami(msh_TriMesh *p, Vint type, Vint *iparam)
get mesh generation parameters
-
void msh_TriMeshSetParamd(msh_TriMesh *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 triangles.
The parameter
VIS_MESH_EDGELENGTH
specifies a target edge length of generated triangles. The edge length constraint may not be satisfied for triangles generated near the boundary mesh due to boundary mesh sizing. By defaultVIS_MESH_EDGELENGTH
is set to 0. and is ignored.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 triangles. The maximum edge length is strictly enforced. By defaultVIS_MESH_MAXEDGELENGTH
is set to 0. and is ignored.The parameter
TRIMESH_CUTANGLE
specifies a deviation from parallel in degrees for use with the rectilinear boundary cut option enabled usingmsh_TriMeshSetParami()
with parameterTRIMESH_BOUNDCUT
. By defaultTRIMESH_CUTANGLE
is set to 2. degrees.The parameter
TRIMESH_CUTMINEDGELENGTH
specifies a minimum edge length which will be generated along an existing edge by a point introduced by a boundary cut. By defaultTRIMESH_CUTMINEDGELENGTH
is set to 0.Inquire of defined double precision parameters as output arguments using
void msh_TriMeshGetParamd (msh_TriMesh *trimesh, Vint ptype, Vdouble *dparam)
- Errors
VIS_ERROR_ENUM
is generated if an improper ptype is specified.VIS_ERROR_VALUE
is generated if an improper target edge length is specified.
- Parameters:
p – Pointer to TriMesh 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 =TRIMESH_CUTANGLE Deviation from parallel angle =TRIMESH_CUTMINEDGELENGTH Minimum edge length generated
dparam – Specifies the value that ptype will be set to.
-
void msh_TriMeshGetParamd(msh_TriMesh *p, Vint type, Vdouble *param)
get mesh generation parameters
-
void msh_TriMeshAbort(msh_TriMesh *p)
set abort flag
Set mesh generation abort flag. During the mesh generation process in
msh_TriMeshGenerate()
ormsh_TriMeshRefine()
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_TriMeshSetFunction()
.- Parameters:
p – Pointer to TriMesh object.
-
void msh_TriMeshSetGeomSizing(msh_TriMesh *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 TriMesh 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_TriMeshGetConnect(msh_TriMesh *p, vis_Connect **connect)
get internal Connect object
Get the internal Connect object which contains the user input nodes and lines. The Connect object must not be changed by the user.
-
void msh_TriMeshWrite(msh_TriMesh *p, Vint type, const Vchar *path)
write input description to file
Write a mesh description to a file which has been input to TriMesh. The contents of the file is intended for internal use only. This function should be called after all input and meshing parameters have been set to ensure that an accurate input 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 mesh description in the case that the TriMesh 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 .tri be used for ASCII files and .btri 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 has not been defined.
- Parameters:
p – Pointer to TriMesh object.
type – File type
x=SYS_ASCII ASCII format =SYS_BINARY Binary format
path – File path
-
void msh_TriMeshRead(msh_TriMesh *p, Vint type, const Vchar *path)
read input description from file
Read a mesh description from a file which has been previously written by
msh_TriMeshWrite()
. 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 TriMesh is read.
- Parameters:
p – Pointer to TriMesh object.
type – File type
x=SYS_ASCII ASCII format =SYS_BINARY Binary format
path – File path
-
void msh_TriMeshGetInteger(msh_TriMesh *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_TriMeshGenerate()
ormsh_TriMeshRefine()
.Use
TRIMESH_NUMFREEEND
to determine the number of lines which contain a free end. The queryTRIMESH_FREEEND
will return the list of all free end lines detected.Use
TRIMESH_NUMEXNODE
to determine the number of input nodes which lie outside the domain. The queryTRIMESH_EXNODE
will return the list of all nodes outside the domain.The query type
TRIMESH_NUMNONMAN
returns the number of non-manifold lines detected in the input. The queryTRIMESH_NONMAN
will return the list of non-manifold lines detected.The query type
TRIMESH_NUMINCONS
returns the number of inconsistent lines detected in the input. The queryTRIMESH_INCONS
will return the list of inconsistent lines detected.Use
TRIMESH_NUMINTERSECT
to determine the number of intersecting input lines. The queryTRIMESH_INTERSECT
will return the list of all intersecting input lines.The query type
TRIMESH_NUMUNCONN
returns the number of unconnected nodes detected in the input. The queryTRIMESH_UNCONN
will return the list of unconnected nodes detected.Use
TRIMESH_NUMREPEATEDNODE
to determine the number of input lines with repeated nodes. The queryTRIMESH_REPEATEDNODE
will return the list of all input lines with repeated nodes.Use
TRIMESH_NUMDUPLICATELINE
to determine the number of duplicate input lines. The queryTRIMESH_DUPLICATELINE
will return the list of all input lines duplicated.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 TriMesh object.
type – Type of integer information to query
x=TRIMESH_NUMNONMAN Number of non-manifold lines =TRIMESH_NONMAN Non-manifold lines =TRIMESH_NUMEXNODE Number of nodes outside domain =TRIMESH_EXNODE Nodes outside domain =TRIMESH_NUMINCONS Number of inconsistent lines =TRIMESH_INCONS Inconsistent lines =TRIMESH_NUMINTERSECT Number of intersecting lines =TRIMESH_INTERSECT Intersecting lines =TRIMESH_NUMFREEEND Number of free end lines =TRIMESH_FREEEND Free end lines =TRIMESH_NUMUNCONN Number of unconnected nodes =TRIMESH_UNCONN Unconnected nodes =TRIMESH_NUMREPEATEDNODE Number of lines with repeated nodes =TRIMESH_REPEATEDNODE Lines with repeated nodes =TRIMESH_NUMDUPLICATELINE Number of duplicate lines =TRIMESH_DUPLICATELINE Duplicate lines =TRIMESH_NUMNEGJAC Number of negative Jacobians =TRIMESH_NUMNONQUAL Number of bad jacobian ratio =VIS_MESH_PROGRESS Current meshing progress
iparam – [out] Returned integer information
-
void msh_TriMeshComputeArea(msh_TriMesh *p, Vdouble *area)
compute area
Compute area enclosed by input boundary.
- Parameters:
p – Pointer to TriMesh object.
area – [out] Computed area
-
void msh_TriMeshGenerate(msh_TriMesh *p, vis_Connect *connect)
generate finite element mesh
Generate interior nodes and triangle elements. The generated nodes and elements are appended to any existing nodes and elements in the input Connect object. The input boundary lines are checked for free ends, connectivity consistency, enclosing negative or zero area, etc. Use
msh_TriMeshGetInteger()
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 boundary lines.VIS_ERROR_COMPUTE
is generated if mesh generation fails during triangulation.
-
void msh_TriMeshRefine(msh_TriMesh *p, vis_State *state, vis_Connect *connect)
refine a triangle mesh
Refine a previously generated or refined mesh. The input State object must be an element scalar state of target element edge lengths. The number of elements in State must be exactly equal to the number of elements resulting from the last call to
msh_TriMeshGenerate()
ormsh_TriMeshRefine()
. If an element edge length is zero, it is ignored. The refinement process does not coarsen the mesh, all previously existing element corner nodes are preserved and by default their coordinate locations are unaltered. Usemsh_TriMeshSetParami()
with parameterVIS_MESH_REFINESMOOTH
to control smoothing node locations.If an orphan mesh is to be refined - i.e., a mesh not generated with TriMeshGenerate, use
msh_TriMeshSetElem()
to set the initial triangle mesh. The boundary of the orphan mesh must be specified withmsh_TriMeshSetLine()
and, if needed,msh_TriMeshSetLineAssoc()
, etc.The generated nodes and elements are appended to any existing nodes and elements in the input Connect object. Use
msh_TriMeshGetInteger()
is query for detailed information about the mesh generation process including errors.- Errors
VIS_ERROR_OPERATION
is generated if TriMesh has not been used to generate the original mesh or an orphan mesh is to be refined andmsh_TriMeshSetElem()
has not been called.VIS_ERROR_OPERATION
is generated if the input State object is not an element scalar state with the proper number of element entities.VIS_ERROR_OPERATION
is generated if the mesh generation fails.