[p,l,s,c]=gmbasics;
The four return variables are four breps, which are
the following: a point embedded in R^0,
a line segment embedded in R^1, a square embedded
in R^2 and a cube embedded in R^3. The segment is [0,1], the
square is [0,1] x [0,1] and the cube is [0,1] x [0,1] x [0,1].
brep = gmmouse;
This routine allows the user to click the mouse in the drawing
window to draw the boundary a 2-D brep. The boundary is composed
of one or more closed paths. A brep is returned.
The user should not click a path in which boundary segments cross through each other. The program does not currently check for this error condition.
brep = gmpolygon(k);
Here, k is an integer 3 or larger. This routine constructs a regular
polygon (a 2D brep) with k sides centered at the origin. It is scaled
so that the inscribed circle of the polygon has radius=1.
There are also routines to construct affine transformation matrices. These transformation matrices are meant to be used by the gmapply function. A transformation matrix must be d-by-(d+1).
mat = gmdilate(arg_1,arg_2,..,arg_d);
Each arg is a nonzero real number. The result is a d-by-(d+1) transformation
that does dilation of the the coordinates. Negative entries do reflection.
The number of arguments d can be as large as 6.
mat = gmrotate(angle, d, coord1, coord2);
This routine constructs a d-dimensional rotation in plane
coord1-coord2. The amount of the rotation is angle radians in the
counterclockwise direction.
If there is only one argument to this routine, then d=2 is
assumed, and the plane rotation is in the (x_1,x_2) plane.
mat = gmtranslate(arg_1,arg_2,...,arg_d);
Each arg is a real number. This routine constructs a d-dimensional
translation operator, where the translation amounts in each
dimension are the arguments. The number of arguments can be
as large as 6.
chunk2 = gmapply(mat,chunk1);
Here, chunk1 is a brep or simplicial complex, and mat is an
affine transformation of the correct dimension. The affine
transformation is applied to the object, and the transformed
object is returned. This is the only unary operation that can
be applied to simplicial complexes as well as breps.
The gmapply function, when applied to breps, assigns the same property-value lists to faces of the output brep that were assigned to the input brep.
The mat
argument can be any matlab matrix of
size d-by-(d+1), but most often this matrix is constructed
using some combination of gmrotate
, gmtranslate
,
gmdilate
, and gmcompose
.
brep2=gmcone(brep1,point);
If brep1 has embedded dimension d, this creates a new brep
with embedded dimension d+1. The new brep is a cone whose base
is a copy of brep1 embedded in the x_{d+1}=0 plane. The second
argument is the apex of the cone. It must be a d+1-vector whose
last entry is not zero.
The gmcone function ignores the property-value lists of faces of the input brep, and its output brep has no property-value pairs initially.
brep2=gmembed(brep1);
If brep1 has embedded dimension d, this creates a new brep
with embedded dimension d+1. The new brep is simply a copy
of the old one, embedded in the x_{d+1}=0 plane of the
higher dimensional space.
The gmembed function, like gmapply, preserves property-value pairs.
brep2=gmprism(brep1);
If brep1 has embedded dimension d, this creates a new brep
with embedded dimension d+1. The new brep is the prism
formed by joining a copy of brep1 embedded in the x_{d+1}=0
plane to a parallel copy in the x_{d+1}=1 plane.
The gmprism function, like gmcone, clears property-value pairs.
If some of the faces of A or B have property-value pairs, then if this same face appears (either partly or in whole) in C, then its appearance in C is also marked with the same property-value pairs.
The algorithm for boolean operations is a fairly naive looping over all pairs of faces. The complexity of a boolean operation is O(mn(m+n)) where m is the combinatorial complexity of the first object and n is the combinatorial complexity of the second object. The standard reference for boolean operations on geometric objects is C. Hoffmann, ``Geometric and Solid Modeling: An Introduction,'' Morgan Kaufmann, 1989. The algorithms used in gm_bool for processing degeneracies (see below) are not the algorithms in Hoffmann's book, however. The optimizations described in the book, such as face boxing, have not been implemented.
If the intrinsic dimension is some other integer, this is called a degeneracy. The set boolean routines check for degeneracies, and if so, carefully work around them in order to produce a correct brep. We believe that, in the presence of exact arithmetic, the set-boolean operations correctly handle degeneracies of any dimension. However, because the operations are implemented in double-precision floating point, there are cases when the degeneracy routines might fail. For instance
c = gmunite(a,gmapply(gmrotate(1e-14),a));
(that is, the union of a two-dimensional brep with a copy of it
rotated by 10^{-14} radians) will very likely cause the routine
gm_bool to generate an error.
If there are degeneracies, then gm_bool also tries to ``simplify'' the resulting brep. For example, consider uniting a unit square with a translation of the unit square by 1 unit in the x dimension. The result is a 2-by-1 rectangle. Without simplification, it will be stored as a 6-sided polygon. The simplification routine notices that the two top edges of this 6-sided polygon are colinear, as are the bottom edges. So it merges these pairs of edges into two edges of length 2. Next, the simplification routine notices that there is a useless vertex in the middle of each of these two edges, so it deletes the two vertices. Thus, the user gets back a 2-by-1 rectangle with four edges and four vertices, which is presumable what he/she wanted.
In the second variant of the gmsubtract routine described above, the ``simplify'' routine is less aggressive about cleaning up the brep (else it would clean up the SL internal boundary that was just created!)
The gm_bool routine (mainly source file bool.C) is the second most complicated routine in all of QMG. (The most complicated is of course the mesh generator.)
mat = gmcompose(mat1,mat2);
This function composes two affine transformations, and returns
the composition. Thus, result = gmapply(mat1, gmapply(mat2, obj))
should be equivalent to
result = gmapply(gmcompose(mat1, mat2), obj)
for any obj, but the latter is more efficient.
gmget(obj)
prints a list of properties and values of an object. For a brep, the
properties are Type (which is the string ``brep''), EmbeddedDim,
IntrinsicDim, Level0Size, and so on up to LeveldSize, where
d is an integer that is the intrinsic dimension. There are
no values returned. Please see the documentation on breps.
For a simplicial complex, the properties are Type (which is the string ``simpcomplex''), EmbeddedDim, IntrinsicDim, NumVertices, NumSimplices, Vertices, Simplices, and SourceSpec. Please see the documentation on simplicial complexes.
The second way of calling gmget is
val = gmget(obj,prop);
Here, prop is a string that is the name of a valid property
for the object (e.g. 'IntrinsicDim'
) in question.
The return value is the value of this property.
For this call, the names of the properties are not case sensitive, i.e., intrinsicdim is the same as IntrinsicDim.
gmgetf(brep,fdim,faceind)
displays all the properties and values associated with face
of dimension fdim, index number faceind, of the brep.
Note that fdim must lie between 0 and the intrinsic dim of the brep
and faceind must lie between 0 and LeveldSize-1, where
d stands for fdim,
These properties always include the built-in properties of Children,
InternalBoundary, and Matrix. The remaining properties are
user-defined and can be inserted with the gm_addpropval function. A
property that is meaningful for the mesh generator is
size_control
. Properties that are meaningful for the
finite-element solver include bc
and
conductivity
. The color
property
is meaningful for the graphics routines. The user may define additional
properties at his/her discretion.
The second way to call the routine is
val = gmgetf(brep,fdim,faceind,prop);
Here, prop is a string that is a property name. If the brep face
in question does not have the property, then either an empty string
or an empty matrix [] is returned. In either case, this can
be tested by checking the length of val.
Property names for faces are case-sensitive.
newbrep=gm_addpropval(brep, fdim, faceind, prop, val);
This function adds a new property-value pair to a face of a
brep, and returns the so-modified brep. See the description
of gmgetf for the interpretation of fdim and faceind. Arguments
prop and val are both strings. The built-in properties (Children,
InternalBoundary, Matrix) cannot be changed.
If the face already has a property-value pair with the given name, then this function returns a brep with the old property-value pair replaced by the new one.
To delete a property-value pair from a face, just set the val to the empty string.
As mentioned above, property names of faces are case sensitive.
Note: this function does not affect the brep in the argument list. Thus, if you call gm_addpropval with no return argument, then it doesn't do anything. This is true of the entire QMG package, and of Matlab in general. Routines do not modify their calling arguments.
obj = gm_in(str);
where str is a string and obj is an object. This routine takes
a brep or simplicial complex in ascii form and returns the same
object in chunk form.
The various formats are documented
elsewhere.
str = gm_out(obj);
where str is a string and obj is an object, either a brep or
simplicial complex in chunk form.
The various formats are documented
elsewhere.
gm_write(obj,str);
where str is a string and obj is an object. This str argument
should be the name of a file in the current directory. This routine
writes a brep or simplicial complexes in Ascii format into the
file.
This routine is not the preferred way of saving breps; instead, the
preferred way is to save the chunk format with the ordinary matlab
save
command. But gm_write
is needed if you want to
edit property-value pairs directly, if you want to make an input for the standalone
mesh generator, or if you need to pass an object between different
hardware architectures.
obj = gm_rd(str);
where str is a string and obj is an object. This str argument
should be the name of a file in the current directory. This routine
reads in a brep/simplicial complex in Ascii format that was saved by
the gm_write function.
gmmeshgen
, which
is an m-file front-end for the mex file gm_meshgen. The mesh
generator is documented separately.
After every single run of the mesh generator,
you should immediate execute gmchecktri
, which is also documented separately.
A mesh can be refined; that is, each of its
simplices can be subdivided into 2^d subsimplices by gmrefine
.
This routine is documented separately.
The finite element solver is named gmfem. It is
documented separately. Before calling gmfem, you
must set boundary conditions on the brep using the
gm_addpropval
routine. You must also have generated a mesh with the
gmmeshgen function. Routine gmfem and its subsidiary routines
are all written in Matlab.
In order to figure out how QMG has numbered the faces of your brep (you must know the numbering to use gm_addproval) it is very helpful to use the graphics routines described below.
The package uses the convention that faces with no color property default to <0 0 0> unless some other default is specified. Color <0 0 0> is interpreted as `invisible.' (Usually <0 0 0> means black, but QMG interprets it as `invisible'). If you really need black, you can get very close to black with <0 0 0.001>.
Having a color that stands for invisible is very useful for 3D breps, because the only way to see internal boundaries or completely interior holes is to make some exterior faces invisible.
The graphics routines are as follows.
newbrep=gmrndcolor(brep, dim);
This routine produces a new brep that is the same as the old
brep, except that all faces of dimension dim are colored according
to a random ordering of 18 pre-designated colors. (Just like
gm_addpropval, this routine does not change its calling argument
and therefore is useless without a return variable.)
The dim argument is optional; if it is not specified, then it defaults to the intrinsic dimension of the brep minus 1 (i.e., facets).
If a face is already colored, then gmrndcolor leaves it with the color that it has.
gmshowcolor(brep, dim);
This routine displays a diagram showing the colors of the
faces of the brep of dimension dim. Using this routine in
conjunction with gmrndcolor and either gmviz or gmgeomview
is very helpful for figuring out the numbering of the faces
so that you can assign boundary conditions for the finite element package.
The dim argument is optional; if it is not specified, then it defaults to the intrinsic dimension of the brep minus 1 (i.e., facets).
gmviz(obj, colorspec, maxdim);
This creates a matlab window to display the object obj.
There is a column of buttons in the window which are described below.
Colorspec is used to specify the color to be used for faces that have no color assigned. If you omit this argument, it defaults to [0,0,0] which stands for `invisible.' Instead of a 3-vector, you can specify a one-character string like `r' for red.
Argument maxdim is the maximum dimension to show. For instance, if you display a 3D brep with maxdim=1 then you will get a line-diagram of its edges. If maxdim=2 you will get a rendering of its facets. This argument defaults to one less than the intrinsic dimension.
The object obj can be either a 2D brep, a 2D simplicial complex, or a 3D brep.
Note that if you specify a simplicial complex as the first argument, then you must specify a colorspec as the second argument because there is no way to assign colors to faces of a simplicial complex (and hence they all default to `invisible').
If you are displaying a 2D object, the meaning of the buttons in the GUI is as follows. `Zoom in' and `Zoom out' zoom in on the current center of the object. `Up', `Down', `Left', and `Right' move your viewpoint on the object.
If you are displaying a 3D brep, then `Zoom in' and `Zoom out' don't do anything. The buttons `Up', `Down', `Left' and `Right' change your viewpoint as if you are sitting on the surface of a circumscribing sphere around the object.
If the object is very complex, these buttons can take a long time to work,
In both 2D and 3D viewing modes, the `Clear' button clears objects plotted on the axis. The `Autoclear' check-button, if set, will cause the axis to get cleared for each new object displayed. Sometimes the axis is cleared even if this button is off (for instance, if you switch from 2D viewing to 3D viewing).
The rendering of 3D objects unfortunately doesn't work well in Matlab. There are two things that go wrong:
gmgeomview(brep, colorspec, dim, name);
This routine sends an RPC command to the geomview server to display
the brep, which must have embedded dimension of 3. The remaining
arguments are optional. Colorspec is the color to be used for
brep faces that have no color property-value pair. The default
is invisible. Dim is the dimension to show; dim=1 means show
a line diagram, and dim=2 means show the facets. The default for this
argument is the object's intrinsic dim minus one.
Geomview maintains its own namespace of objects. You can specify a name as the fourth argument (a string), in which case the name gets uploaded to geomview. The name should have only letters, digits and underscores. The default name is `matlab_obj.'
The advantage of specifying a name is that geomview can keep track of multiple objects at once if they have different names. If no name is specified, then each object is called `matlab_obj' and hence replaces the previous object called `matlab_obj.'
In order to use this function, you must have already initialized the geomview server. See the relevant section of the installation guide to learn how to initialize the geomview server.
If you have never used geomview before, you should read the manual to understand all the things you can do with displayed objects. In its default mode when it starts, if you drag the mouse across the object (which should displayed by geomview in a window marked `camera') while holding down the left mouse button, then the object will rotate in the direction that you dragged the mouse. The faster you drag, the faster it rotates. If you release the mouse button while you are still dragging, then it continues to rotate under its own `inertia.'
gmplot(scomplex, field);
where scomplex is the simplicial complex used for the finite element
method, and field is the solution vector. Thus, field should be
a column vector with real-number entries whose number of entries
is equal to the number of nodes in the simplicial complex. This is
the form of the return variable from gmfem. But you can also use
this function to show, for instance, the difference between two fem
solutions.
The routine uses the current colormap (jet is recommended). It uses interpolated colors on each triangle in the simplicial complex. However, it does not use the interpolated-color feature that is built in to matlab (I found that earlier versions of interpolated colors did not work on my color printer. I don't know if this has been fixed in the more recent Matlab releases). Instead, it divides each triangle into nine smaller triangles and uses flat colors (that are linearly interpolated by the gmplot program) on each of the nine subtriangles.
The matlab command colorbar
is useful in conjunction
with the gmplot routine.
This code is written in Matlab and involves a lot of scalar operations, plus it's a quadratic-time algorithm so it's somewhat slow. (The theoretically optimal algorithm for polygon triangulation requires linear time and is due to Chazelle, but that algorithm is probably too complex for practical use.)
This routine is necessary because neither Matlab nor Geomview can render polygons with holes. Since brep faces can be polygons with holes, each face must first be subdivided into polygons without holes before rendering. This is the purpose of gm_ftri.
The boundaries between the triangles on a given face are invisible in both gmviz and gmgeomview, but if you switch the default options in geomview, you will be able to see the individual triangles.
If you have your own 3D renderer that you want to connect to QMG, you might also find that gm_ftri is useful.
string = gm_tclexecute(filename, command);
This routine executes a TCL instruction. It creates an interpreter
and deletes it for just the one instruction, so no state is maintained
between consecutive calls to gm_tclexecute. The first argument
is name of a file that is `source'd by the TCL interpreter. (If
the string is empty, no file is read in.)
The second argument is the command itself. The return value,
a string, is whatever string is returned by the command.
The routine initializes the TCL interpreter and also the TCL-DP system, but not the Tk system (because Tk is very slow to initialize). For a quick overview of TCL and TCL-DP, see the installation guide.
This function is used by gmgeomview to send a command to the TCL interpreter that is the front-end to geomview.
This documentation is written by Stephen A. Vavasis and is copyright (c) 1995 by Cornell University. Permission to reproduce this documentation is granted provided this notice remains attached. There is no warranty of any kind on this software or its documentation. See the accompanying file 'Copyright' for a full statement of the copyright.
Back to mesh generation home page.
ref.html,v 1.1.1.1 1995/05/05 01:15:55 vavasis Exp
Stephen A. Vavasis, Computer Science Department, Cornell University, Ithaca, NY 14853, vavasis@cs.cornell.edu