Files | |
file | aprx.h |
Class Aprx: the APpRoXimation library object. | |
Classes | |
struct | sAprx |
Class Aprx: Parameters and datatypes Class Aprx: Definition. More... | |
Typedefs | |
typedef struct sAprx | Aprx |
Declaration of the Aprx class as the Aprx structure. | |
Functions | |
Aprx * | Aprx_ctor (Vmem *vmem, Gem *tgm, PDE *tpde) |
Construct a linear Approximation object. | |
Aprx * | Aprx_ctor2 (Vmem *vmem, Gem *tgm, PDE *tpde, int order) |
Approximation object constructor. | |
void | Aprx_dtor (Aprx **thee) |
Approximation object destructor. | |
int | Aprx_numB (Aprx *thee) |
Return the number of blocks. | |
Bmat * | Aprx_P (Aprx *thee) |
Return a pointer to the current prolongation matrix. | |
int | Aprx_read (Aprx *thee, int key, Vio *sock) |
Read in the user-specified initial mesh given in the "MCSF" or "MCEF" format, and transform into our internal datastructures. Do a little more than a "Gem_read", in that we also initialize the extrinsic and intrinsic spatial dimensions corresponding to the input mesh, and we also then build the reference elements. | |
void | Aprx_reset (Aprx *thee) |
Reset all of the datastructures. | |
void | Aprx_initNodes (Aprx *thee) |
Make node information. | |
void | Aprx_buildNodes (Aprx *thee) |
Determine the node types via a single fast traversal of the elements by looking at the element faces. | |
void | Aprx_evalTrace (Aprx *thee, Bvec *Wud, Bvec *Wui, Bvec *Wut) |
Evaluate the trace information. | |
int | Aprx_refine (Aprx *thee, int rkey, int bkey, int pkey) |
Refine the mesh. | |
int | Aprx_unRefine (Aprx *thee, int rkey, int pkey) |
Unrefine the mesh. | |
int | Aprx_deform (Aprx *thee, Bvec *def) |
Deform the mesh. | |
void | Aprx_buildProlong (Aprx *thee, int numRc, int numR, int pkey) |
Create storage for and build prolongation operator; the operator prolongates a vector FROM a coarser level TO this level. | |
void | Aprx_buildProlongBmat (Aprx *thee, Mat *p) |
Create block version of prolongation operator; also updates the Bnode object. | |
int | Aprx_nodeCount (Aprx *thee, Re *re) |
Count up the total number of basis functions (nodes). This is basically an a priori count of the number of rows in the stiffness matrix. | |
void | Aprx_interact (Aprx *thee, Re *re, Re *reT, int *numR, int *numO, int *numOYR, int **IJA, int **IJAYR) |
Count the number of basis (linear/quadratic/etc) functions which interact with a given basis function. This is basically an a priori count of the number of nonzeros per row in the stiffness matrix, and the number of nozeros per row above the diagonal (and thus the number of nonzeros per column below the diagonal if we have symmetric nonzero structure). By doing this count, we can one-time malloc storage for matrices before assembly, and avoid extra copies as well as avoid the use of linked list structures for representing matrices with a priori unknown nonzero structure. | |
void | Aprx_interactBlock (Aprx *thee, Re *re[MAXV], Re *reT[MAXV], MATsym sym[MAXV][MAXV], MATmirror mirror[MAXV][MAXV], MATformat frmt[MAXV][MAXV], int numR[MAXV], int numO[MAXV][MAXV], int *IJA[MAXV][MAXV]) |
Determine all basis function interations. | |
int | Aprx_nodeIndex (Aprx *thee, Re *re, int i, int dim, int q) |
Compute global index of node of dimension dim in the matrix. | |
int | Aprx_nodeIndexTT (Aprx *thee, Re *re, TT *t, int idf) |
Given an index of node DF within the simplex and its TT structure, compute global index for that node. | |
void | Aprx_buildBRC (Aprx *thee, Bmat *A, Bmat *M) |
Build boundary Bnode information into two Bmats. | |
void | Aprx_initSpace (Aprx *thee, int i, int r, double phi[], double phix[][3], double V[], double dV[][3]) |
Initialize vector test function "i" and its gradient, generated by scalar test function "r" for this element. | |
void | Aprx_initEmat (Aprx *thee, int bumpKey, TT *t, Emat *em) |
Initialize an element matrix. | |
void | Aprx_quadEmat (Aprx *thee, int evalKey, int energyKey, int residKey, int tangKey, int massKey, int bumpKey, int face, TT *t, Emat *em, Bvec *Wu, Bvec *Wud) |
Build an element matrix and an element load vector for a single given element, using quadrature. | |
void | Aprx_fanEmat (Aprx *thee, int evalKey, int energyKey, int residKey, int tangKey, int massKey, int bumpKey, Emat *em, Bmat *A, Bmat *M, Bvec *F) |
Fan out an element matrix to the global matrix, and fan out an element vector to the global vector. | |
void | Aprx_assemEmat (Aprx *thee, SS *sm, int evalKey, int energyKey, int residKey, int tangKey, int massKey, int bumpKey, Emat *em, Bvec *Wu, Bvec *Wud) |
Assemble one element matrix. | |
double | Aprx_assem (Aprx *thee, int evalKey, int energyKey, int residKey, int tangKey, int massKey, int bumpKey, int ip[], double rp[], Bmat *A, Bmat *M, Bvec *F, Bvec *Wu, Bvec *Wud) |
Assemble the tangent matrix and nonlinear residual vector (which reduces to the normal stiffness matrix and load vector in the case of a linear problem) for a general nonlinear system of m components in spatial dimension d. | |
int | Aprx_markRefine (Aprx *thee, int key, int color, int bkey, double elevel, Bvec *u, Bvec *ud, Bvec *f, Bvec *r) |
Mark simplices to be refined. | |
int | Aprx_estRefine (Aprx *thee, int key, int color, int bkey, double elevel, Bvec *u, Bvec *ud, Bvec *f, Bvec *r) |
A posteriori error estimation via several different approaches. | |
int | Aprx_markRefineFixed (Aprx *thee, int num2ref, int color) |
Selects a fixed number/fraction of simplicies to be refined. | |
void | Aprx_evalFunc (Aprx *thee, Bvec *u, int block, int numPts, double *pts, double *vals, int *marks) |
Evaluate a finite element solution at a set of arbitrary points. | |
void | Aprx_fe2fd (Aprx *thee, Bvec *u, int block, int nx, int ny, int nz, double x0, double y0, double z0, double dx, double dy, double dz, int derivs, int outTypeKey) |
Interpolates data onto finite difference grid. | |
void | Aprx_bndIntegral (Aprx *thee, Bvec *u, Bvec *ud, Bvec *ut) |
Evaluate a given boundary integral. | |
void | Aprx_admMass (Aprx *thee, Bvec *u, Bvec *ud, Bvec *ut, int block) |
Evaluate ADM mass from conformal factor. | |
double | Aprx_evalError (Aprx *thee, int pcolor, int key, Bvec *u, Bvec *ud, Bvec *ut) |
Evaluate error in solution in one of several norms. | |
double | Aprx_evalH1 (Aprx *thee, SS *sm, Bvec *u, Bvec *ud, Bvec *ut) |
Evaluate difference of two functions in H1 norm in one element. | |
void | Aprx_writeGEOM (Aprx *thee, Vio *sock, int defKey, int colKey, int chartKey, double gluVal, int fkey, Bvec *u, char *format) |
Write a finite element mesh or mesh function in some format. | |
void | Aprx_writeSOL (Aprx *thee, Vio *sock, Bvec *u, char *format) |
Write a finite element mesh or mesh function in some format. | |
int | Aprx_partSet (Aprx *thee, int pcolor) |
Clear the partitioning (set all colors to pcolor). | |
int | Aprx_partSmooth (Aprx *thee) |
Smooth the partitioning. | |
int | Aprx_partInert (Aprx *thee, int pcolor, int numC, double *evec, simHelper *simH) |
Partition the domain using inertial bisection. Partition sets of points in R^d (d=2 or d=3) by viewing them as point masses of a rigid body, and by then employing the classical mechanics ideas of inertia and Euler axes. | |
int | Aprx_partSpect (Aprx *thee, int pcolor, int numC, double *evec, simHelper *simH, int *ford, int *rord, int general) |
Partition the domain using spectral bisection. | |
int | Aprx_part (Aprx *thee, int pkey, int pwht, int ppow) |
Do recursive bisection. | |
int | Aprx_partOne (Aprx *thee, int pkey, int pwht, int pcolor, int poff) |
Do a single bisection step. |
Evaluate ADM mass from conformal factor.
thee | Pointer to an Aprx allocated memory location | |
u | block NODAL solution vector | |
ud | block NODAL dirichlet vector | |
ut | block NODAL analytical vector | |
block | index of a block |
double Aprx_assem | ( | Aprx * | thee, | |
int | evalKey, | |||
int | energyKey, | |||
int | residKey, | |||
int | tangKey, | |||
int | massKey, | |||
int | bumpKey, | |||
int | ip[], | |||
double | rp[], | |||
Bmat * | A, | |||
Bmat * | M, | |||
Bvec * | F, | |||
Bvec * | Wu, | |||
Bvec * | Wud | |||
) |
Assemble the tangent matrix and nonlinear residual vector (which reduces to the normal stiffness matrix and load vector in the case of a linear problem) for a general nonlinear system of m components in spatial dimension d.
thee | Pointer to an Aprx allocated memory location | |
evalKey | == 0 ==> Primal problem: evaluation at z=[0+ud]. == 1 ==> Primal problem: evaluation at z=[u+ud]. == 2 ==> Dual problem: evaluation at z=[0+ud]. == 3 ==> Dual problem: evaluation at z=[u+ud]. == ? ==> Same as 0. | |
energyKey | == 0 ==> DON'T assemble an energy. == 1 ==> Assemble the energy J(z). == 2 ==> Assemble the energy 0. == ? ==> Same as 0. | |
residKey | == 0 ==> DON'T assemble a residual. == 1 ==> Assemble the residual F(z)(v). == 2 ==> Assemble the residual 0. == ? ==> Same as 0. | |
tangKey | == 0 ==> DON'T assemble a tangent matrix. == 1 ==> Assemble the tangent matrix DF(z)(w,v). == 2 ==> Assemble the tangent matrix 0. == ? ==> Same as 0. | |
massKey | == 0 ==> DON'T assemble a mass matrix. == 1 ==> Assemble the mass matrix p(w,v). == 2 ==> Assemble the mass matrix 0. == ? ==> Same as 0. | |
bumpKey | == 0 ==> Do not assemble with bumps. == 1 ==> Assemble bilinear form with bumps. == 2 ==> Assemble residual form with bumps. == ? ==> Same as 1. | |
ip | index for assembled energies | |
rp | parameter for initially assembled PDE | |
A | Pointer to stiffness matrix | |
M | Pointer to Mass matrix | |
F | block NODAL residual vector | |
Wu | block NODAL solution work vector | |
Wud | block NODAL dirichlet work vector |
void Aprx_assemEmat | ( | Aprx * | thee, | |
SS * | sm, | |||
int | evalKey, | |||
int | energyKey, | |||
int | residKey, | |||
int | tangKey, | |||
int | massKey, | |||
int | bumpKey, | |||
Emat * | em, | |||
Bvec * | Wu, | |||
Bvec * | Wud | |||
) |
Assemble one element matrix.
thee | Pointer to an Aprx allocated memory location | |
sm | Pointer to a given simplex | |
evalKey | == 0 ==> Primal problem: evaluation at z=[0+ud]. == 1 ==> Primal problem: evaluation at z=[u+ud]. == 2 ==> Dual problem: evaluation at z=[0+ud]. == 3 ==> Dual problem: evaluation at z=[u+ud]. == ? ==> Same as 0. | |
energyKey | == 0 ==> DON'T assemble an energy. == 1 ==> Assemble the energy J(z). == 2 ==> Assemble the energy 0. == ? ==> Same as 0. | |
residKey | == 0 ==> DON'T assemble a residual. == 1 ==> Assemble the residual F(z)(v). == 2 ==> Assemble the residual 0. == ? ==> Same as 0. | |
tangKey | == 0 ==> DON'T assemble a tangent matrix. == 1 ==> Assemble the tangent matrix DF(z)(w,v). == 2 ==> Assemble the tangent matrix 0. == ? ==> Same as 0. | |
massKey | == 0 ==> DON'T assemble a mass matrix. == 1 ==> Assemble the mass matrix p(w,v). == 2 ==> Assemble the mass matrix 0. == ? ==> Same as 0. | |
bumpKey | == 0 ==> Do not assemble with bumps. == 1 ==> Assemble bilinear form with bumps. == 2 ==> Assemble residual form with bumps. == ? ==> Same as 1. | |
em | Pointer to element matrix | |
Wu | block NODAL solution work vector | |
Wud | block NODAL dirichlet work vector |
Evaluate a given boundary integral.
thee | Pointer to an Aprx allocated memory location | |
u | block NODAL solution vector | |
ud | block NODAL dirichlet vector | |
ut | block NODAL analytical vector |
Build boundary Bnode information into two Bmats.
thee | Pointer to an Aprx allocated memory location | |
A | Pointer to stiffness matrix | |
M | Pointer to Mass matrix |
void Aprx_buildNodes | ( | Aprx * | thee | ) |
Determine the node types via a single fast traversal of the elements by looking at the element faces.
thee | Pointer to an Aprx allocated memory location |
void Aprx_buildProlong | ( | Aprx * | thee, | |
int | numRc, | |||
int | numR, | |||
int | pkey | |||
) |
Create storage for and build prolongation operator; the operator prolongates a vector FROM a coarser level TO this level.
thee | Pointer to an Aprx allocated memory location | |
numRc | number of vertices before refinement | |
numR | number of vertices after refinement | |
pkey | Boolean sets the pkey type to prolongate a vector |
Create block version of prolongation operator; also updates the Bnode object.
thee | Pointer to an Aprx allocated memory location | |
p | Pointer to Bnode matrix |
Construct a linear Approximation object.
vmem | Memory management object | |
tgm | Geometry manager source | |
tpde | PDE object |
Approximation object constructor.
vmem | Memory management object | |
tgm | Geometry manager source | |
tpde | PDE object | |
order | Order of approximation to be constructed |
Deform the mesh.
thee | Pointer to an Aprx allocated memory location | |
def | The corresponding deforming block vector |
void Aprx_dtor | ( | Aprx ** | thee | ) |
Approximation object destructor.
thee | Pointer to object to be destroyed |
int Aprx_estRefine | ( | Aprx * | thee, | |
int | key, | |||
int | color, | |||
int | bkey, | |||
double | elevel, | |||
Bvec * | u, | |||
Bvec * | ud, | |||
Bvec * | f, | |||
Bvec * | r | |||
) |
A posteriori error estimation via several different approaches.
Apply one of several a posteriori error estimators: If (key == 2) ==> Nonlinear residual indicator. If (key == 3) ==> Dual-weighted residual indicator. If (key == 4) ==> Local problem indicator. If (key == 5) ==> Recovered gradient indicator. The error tolerance can be used in several ways: If (bkey == 0) ==> mark if: [error/S] > TOL If (bkey == 1) ==> mark if: [error/S] > (TOL^2/numS)^{1/2} If (bkey == 2) ==> mark fixed number/fraction (dep. on ETOL): if elevel < 1, mark fraction of simplices if elevel >=1, mark that many simplices (rounded to the smaller integer)
thee | Pointer to an Aprx allocated memory location | |
key | index of different types of a posteriori error estimators | |
color | chart type of marked simplices | |
bkey | Bisection type | |
elevel | level to determine the fraction of simplices | |
u | block NODAL solution vector | |
ud | block NODAL dirichlet vector | |
f | block NODAL residual vector | |
r | block NODAL temporary vector |
Evaluate error in solution in one of several norms.
thee | Pointer to an Aprx allocated memory location | |
pcolor | simplex chart type | |
key | If (key == 0) ==> L^2 norm of the error. If (key == 1) ==> L^{\infty} norm of the error. If (key == 2) ==> H^1 norm of the error. | |
u | block NODAL solution vector | |
ud | block NODAL dirichlet vector | |
ut | block NODAL analytical vector |
void Aprx_evalFunc | ( | Aprx * | thee, | |
Bvec * | u, | |||
int | block, | |||
int | numPts, | |||
double * | pts, | |||
double * | vals, | |||
int * | marks | |||
) |
Evaluate a finite element solution at a set of arbitrary points.
thee | Pointer to an Aprx allocated memory location | |
u | solution vector | |
block | index for the block | |
numPts | number of all interpolation pts for one simplex | |
pts | coordinates of the pt | |
vals | solution for the interpolated pts | |
marks | index for marked or not |
Evaluate difference of two functions in H1 norm in one element.
thee | Pointer to an Aprx allocated memory location | |
sm | simplex | |
u | block NODAL solution vector | |
ud | block NODAL dirichlet vector | |
ut | block NODAL analytical vector |
Evaluate the trace information.
thee | Pointer to an Aprx allocated memory location | |
Wud | block NODAL dirichlet work vector | |
Wui | block NODAL interior dirichlet work vector | |
Wut | block NODAL analytical work vector |
void Aprx_fanEmat | ( | Aprx * | thee, | |
int | evalKey, | |||
int | energyKey, | |||
int | residKey, | |||
int | tangKey, | |||
int | massKey, | |||
int | bumpKey, | |||
Emat * | em, | |||
Bmat * | A, | |||
Bmat * | M, | |||
Bvec * | F | |||
) |
Fan out an element matrix to the global matrix, and fan out an element vector to the global vector.
thee | Pointer to an Aprx allocated memory location | |
evalKey | == 0 ==> Primal problem: evaluation at z=[0+ud]. == 1 ==> Primal problem: evaluation at z=[u+ud]. == 2 ==> Dual problem: evaluation at z=[0+ud]. == 3 ==> Dual problem: evaluation at z=[u+ud]. == ? ==> Same as 0. | |
energyKey | == 0 ==> DON'T assemble an energy. == 1 ==> Assemble the energy J(z). == 2 ==> Assemble the energy 0. == ? ==> Same as 0. | |
residKey | == 0 ==> DON'T assemble a residual. == 1 ==> Assemble the residual F(z)(v). == 2 ==> Assemble the residual 0. == ? ==> Same as 0. | |
tangKey | == 0 ==> DON'T assemble a tangent matrix. == 1 ==> Assemble the tangent matrix DF(z)(w,v). == 2 ==> Assemble the tangent matrix 0. == ? ==> Same as 0. | |
massKey | == 0 ==> DON'T assemble a mass matrix. == 1 ==> Assemble the mass matrix p(w,v). == 2 ==> Assemble the mass matrix 0. == ? ==> Same as 0. | |
bumpKey | == 0 ==> Do not assemble with bumps. == 1 ==> Assemble bilinear form with bumps. == 2 ==> Assemble residual form with bumps. == ? ==> Same as 1. | |
em | Pointer to element matrix | |
A | Pointer to stiffness matrix | |
M | Pointer to Mass matrix | |
F | Pointer to residual vector |
void Aprx_fe2fd | ( | Aprx * | thee, | |
Bvec * | u, | |||
int | block, | |||
int | nx, | |||
int | ny, | |||
int | nz, | |||
double | x0, | |||
double | y0, | |||
double | z0, | |||
double | dx, | |||
double | dy, | |||
double | dz, | |||
int | derivs, | |||
int | outTypeKey | |||
) |
Interpolates data onto finite difference grid.
* Example: Below is a matlab example courtesy of Doug Arnold. * * %-- BEGIN MATLAB CODE -------------------------------------------------- * * function colormesh * * % This file demonstrates how to go from an unstructured grid to * % a structured one. Each point on the structured grid is colored * % according to which triangle in the unstructured grid it belongs to. * * % Douglas N. Arnold, 3/4/98 * * % set up a mesh * npts=10; * x=[ 1.13 1.66 0.50 0.46 0.26 0.20 0.70 1.93 1.01 0.65 ]; * y=[ 1.81 1.23 0.99 1.30 1.29 1.50 1.10 0.68 0.60 0.07 ]; * tri = delaunay(x,y); * ntri=size(tri,1); * * % grid spacing * h=.03; * k=.04; * * clf * [z,i]=sort(rand(64,1)); * colormap hsv; * cmap = colormap; * cmap=cmap(i,:,:); * * % loop over triangles * for m = 1:ntri * % get vertices of triangle i * v1=[x(tri(m,1)),y(tri(m,1))]; * v2=[x(tri(m,2)),y(tri(m,2))]; * v3=[x(tri(m,3)),y(tri(m,3))]; * * % plot the triangle * ttt=[v1;v2;v3;v1]; * line(ttt(:,1),ttt(:,2),'color','green') * * % find bounding box of triangle * xmin = min([v1(1),v2(1),v3(1)]); * xmax = max([v1(1),v2(1),v3(1)]); * ymin = min([v1(2),v2(2),v3(2)]); * ymax = max([v1(2),v2(2),v3(2)]); * * % loop through grid points in bounding box, plotting those in triangle * for i = ceil(xmin/h):floor(xmax/h) * for j = ceil(ymin/k):floor(ymax/k) * if intriangle(v1,v2,v3,[i*h,j*k]) * % set color according to triangle number * color = cmap(m,:); * line(i*h,j*k,'marker','.','markersize',10,'color',color) * end * end * end * * end * * function i = intriangle(v1,v2,v3,x) * % INTRIANGLE -- check if a point belongs to a triangle * % * % f = intriangle(v1,v2,v3,x) * % * % input: * % v1,v2,v3 2-vectors giving the vertices of the triangle * % x 2-vector giving the point * % * % output: * % i 1 if x belongs to the close triangle, 0 else * * % map point under affine mapping from given triangle to ref triangle * xhat = (x-v1)/[v2-v1;v3-v1]; * % check if mapped point belongs to reference triangle * i = xhat(1) >= 0 & xhat(2) >= 0 & 1-xhat(1)-xhat(2) >= 0; * %-- END MATLAB CODE -------------------------------------------------- * * ***************************************************************************
thee | Pointer to an Aprx allocated memory location | |
u | solution vector | |
block | index of block | |
nx | number of grids in X coordinates | |
ny | number of grids in y coordinates | |
nz | number of grids in z coordinates | |
x0 | x coordinate of lower grid corner | |
y0 | y coordinate of lower grid corner | |
z0 | z coordinate of lower grid corner | |
dx | Grid spacing in x direction | |
dy | Grid spacing in y direction | |
dz | Grid spacing in z direction | |
derivs | index for different nfunc types | |
outTypeKey | ASCII output or Binary output |
Initialize an element matrix.
thee | Pointer to an Aprx allocated memory location | |
bumpKey | == 0 ==> Do not assemble with bumps. == 1 ==> Assemble bilinear form with bumps. == 2 ==> Assemble residual form with bumps. == ? ==> Same as 1. | |
t | Pointer to class T | |
em | Pointer to element matrix |
void Aprx_initNodes | ( | Aprx * | thee | ) |
Make node information.
thee | Pointer to an Aprx allocated memory location |
void Aprx_initSpace | ( | Aprx * | thee, | |
int | i, | |||
int | r, | |||
double | phi[], | |||
double | phix[][3], | |||
double | V[], | |||
double | dV[][3] | |||
) |
Initialize vector test function "i" and its gradient, generated by scalar test function "r" for this element.
thee | Pointer to an Aprx allocated memory location | |
i | vector test function index | |
r | scalar test function index | |
phi | values of local basis funcs at integ pt | |
phix | first derivs of local basis funcs at integ pt | |
V | vector test function | |
dV | first derivs of vector test function |
void Aprx_interact | ( | Aprx * | thee, | |
Re * | re, | |||
Re * | reT, | |||
int * | numR, | |||
int * | numO, | |||
int * | numOYR, | |||
int ** | IJA, | |||
int ** | IJAYR | |||
) |
Count the number of basis (linear/quadratic/etc) functions which interact with a given basis function. This is basically an a priori count of the number of nonzeros per row in the stiffness matrix, and the number of nozeros per row above the diagonal (and thus the number of nonzeros per column below the diagonal if we have symmetric nonzero structure). By doing this count, we can one-time malloc storage for matrices before assembly, and avoid extra copies as well as avoid the use of linked list structures for representing matrices with a priori unknown nonzero structure.
thee | Pointer to an Aprx allocated memory location | |
re | a reference simplex element object. | |
reT | a reference simplex element object. | |
numR | num of rows in the matrix | |
numO | num of nonzeros we are actually storing in the strict upper-triangle of matrix. (DRC only) | |
numOYR | YSMP-row(col) | |
IJA | integer structure [ IA ; JA ] | |
IJAYR | lower-triangle block from upper triangle block |
void Aprx_interactBlock | ( | Aprx * | thee, | |
Re * | re[MAXV], | |||
Re * | reT[MAXV], | |||
MATsym | sym[MAXV][MAXV], | |||
MATmirror | mirror[MAXV][MAXV], | |||
MATformat | frmt[MAXV][MAXV], | |||
int | numR[MAXV], | |||
int | numO[MAXV][MAXV], | |||
int * | IJA[MAXV][MAXV] | |||
) |
Determine all basis function interations.
thee | Pointer to an Aprx allocated memory location | |
re | array of reference simplex element objects | |
reT | array of reference simplex element objects | |
sym | symmetry keys for the matrix | |
mirror | block mirror keys for the matrix | |
frmt | Matrix format | |
numR | num of rows in the matrix | |
numO | num of nonzeros we are actually storing in the strict upper-triangle of matrix. (DRC only) | |
IJA | integer structure [ IA ; JA ] |
int Aprx_markRefine | ( | Aprx * | thee, | |
int | key, | |||
int | color, | |||
int | bkey, | |||
double | elevel, | |||
Bvec * | u, | |||
Bvec * | ud, | |||
Bvec * | f, | |||
Bvec * | r | |||
) |
Mark simplices to be refined.
thee | Pointer to an Aprx allocated memory location | |
key | If ( (key == -1) || (key == 0) || (key == 1)) ==> Let geometry manager deal with it. If (key == 2) ==> Nonlinear residual indicator. If (key == 3) ==> Dual-weighted residual indicator. If (key == 4) ==> Local problem indicator. | |
color | chart type of marked simplices | |
bkey | Bisection type | |
elevel | level to determine the fraction of simplices | |
u | block NODAL solution vector | |
ud | block NODAL dirichlet vector | |
f | block NODAL residual vector | |
r | block NODAL temporary vector |
int Aprx_markRefineFixed | ( | Aprx * | thee, | |
int | num2ref, | |||
int | color | |||
) |
Selects a fixed number/fraction of simplicies to be refined.
This function is called by Aprx_markRefine if fixed number / fraction of simplices <ilevel> has to be refined. Returns number of marked elements.
thee | Pointer to an Aprx allocated memory location | |
num2ref | index for the refined simplex | |
color | chart type of marked simplices |
Count up the total number of basis functions (nodes). This is basically an a priori count of the number of rows in the stiffness matrix.
thee | Pointer to an Aprx allocated memory location | |
re | a reference simplex element object. |
Compute global index of node of dimension dim in the matrix.
thee | Pointer to an Aprx allocated memory location | |
re | Pointer to the reference elemented to be used | |
i | i-simplex global number | |
dim | interactions index | |
q | number of volume quadrature points |
Given an index of node DF within the simplex and its TT structure, compute global index for that node.
thee | Pointer to an Aprx allocated memory location | |
re | Pointer to the reference elemented to be used | |
t | Pointer to class T | |
idf | index of node DF within the simplex |
int Aprx_numB | ( | Aprx * | thee | ) |
Return the number of blocks.
thee | Pointer to an Aprx allocated memory location |
Return a pointer to the current prolongation matrix.
thee | Pointer to an Aprx allocated memory location |
int Aprx_part | ( | Aprx * | thee, | |
int | pkey, | |||
int | pwht, | |||
int | ppow | |||
) |
Do recursive bisection.
thee | Pointer to an Aprx allocated memory location | |
pkey | index for different partition option | |
pwht | index for weighted partitioning | |
ppow | partitioning steps |
Partition the domain using inertial bisection. Partition sets of points in R^d (d=2 or d=3) by viewing them as point masses of a rigid body, and by then employing the classical mechanics ideas of inertia and Euler axes.
thee | Pointer to an Aprx allocated memory location | |
pcolor | simplex chart type | |
numC | number of columns in the matrix | |
evec | eigenvector | |
simH | simHelper class |
int Aprx_partOne | ( | Aprx * | thee, | |
int | pkey, | |||
int | pwht, | |||
int | pcolor, | |||
int | poff | |||
) |
Do a single bisection step.
thee | Pointer to an Aprx allocated memory location | |
pkey | partition method | |
pwht | partitioning weighting | |
pcolor | simplex chart type | |
poff | an index |
int Aprx_partSet | ( | Aprx * | thee, | |
int | pcolor | |||
) |
Clear the partitioning (set all colors to pcolor).
thee | Pointer to an Aprx allocated memory location | |
pcolor | simplex chart type |
int Aprx_partSmooth | ( | Aprx * | thee | ) |
Smooth the partitioning.
thee | Pointer to an Aprx allocated memory location |
int Aprx_partSpect | ( | Aprx * | thee, | |
int | pcolor, | |||
int | numC, | |||
double * | evec, | |||
simHelper * | simH, | |||
int * | ford, | |||
int * | rord, | |||
int | general | |||
) |
Partition the domain using spectral bisection.
* We solve the following generalized eigenvalue problem: * * Ax = lambda Bx * * for second smallest eigenpair. We then return the eigenvector * from the pair. We make this happen by turning it into a * regular eigenvalue problem: * * B^{-1/2} A B^{-1/2} ( B^{1/2} x ) = lambda ( B^{1/2} x ) * * or rather * * C y = lambda y, where C=B^{-1/2}AB^{-1/2}, y=B^{1/2}x. * * The matrix "B" is simply a diagonal matrix with a (positive) * error estimate for the element on the diagonal. Therefore, the * matrix B^{-1/2} is a well-defined positive diagonal matrix. * * We explicitly form the matrix C and then send it to the inverse * Rayleigh-quotient iteration to recover the second smallest * eigenpair. On return, we scale the eigenvector y by B^{-1/2} to * recover the actual eigenvector x = B^{-1/2} y. * * To handle the possibility that an element has zero error, in * which case the B matrix had a zero on the diagonal, we set the * corresponding entry in B^{-1/2} to be 1. *
thee | Pointer to an Aprx allocated memory location | |
pcolor | simplex chart type | |
numC | number of columns in the matrix | |
evec | eigenvector | |
simH | simHelper class | |
ford | temporary storage ford array | |
rord | temporary storage rord array | |
general | index for creating different scaling matrix |
void Aprx_quadEmat | ( | Aprx * | thee, | |
int | evalKey, | |||
int | energyKey, | |||
int | residKey, | |||
int | tangKey, | |||
int | massKey, | |||
int | bumpKey, | |||
int | face, | |||
TT * | t, | |||
Emat * | em, | |||
Bvec * | Wu, | |||
Bvec * | Wud | |||
) |
Build an element matrix and an element load vector for a single given element, using quadrature.
thee | Pointer to an Aprx allocated memory location | |
evalKey | == 0 ==> Primal problem: evaluation at z=[0+ud]. == 1 ==> Primal problem: evaluation at z=[u+ud]. == 2 ==> Dual problem: evaluation at z=[0+ud]. == 3 ==> Dual problem: evaluation at z=[u+ud]. == ? ==> Same as 0. | |
energyKey | == 0 ==> DON'T assemble an energy. == 1 ==> Assemble the energy J(z). == 2 ==> Assemble the energy 0. == ? ==> Same as 0. | |
residKey | == 0 ==> DON'T assemble a residual. == 1 ==> Assemble the residual F(z)(v). == 2 ==> Assemble the residual 0. == ? ==> Same as 0. | |
tangKey | == 0 ==> DON'T assemble a tangent matrix. == 1 ==> Assemble the tangent matrix DF(z)(w,v). == 2 ==> Assemble the tangent matrix 0. == ? ==> Same as 0. | |
massKey | == 0 ==> DON'T assemble a mass matrix. == 1 ==> Assemble the mass matrix p(w,v). == 2 ==> Assemble the mass matrix 0. == ? ==> Same as 0. | |
bumpKey | == 0 ==> Do not assemble with bumps. == 1 ==> Assemble bilinear form with bumps. == 2 ==> Assemble residual form with bumps. == ? ==> Same as 1. | |
face | face < 0 ==> do volume quadrature; we take the node numbers as just the identity mapping. face >= 0 ==> do face quadrature for face "face", where the node numbers for face are taken from t->loc[face][]. | |
t | Pointer to class T | |
em | Pointer to element matrix | |
Wu | block NODAL solution work vector | |
Wud | block NODAL dirichlet work vector |
int Aprx_read | ( | Aprx * | thee, | |
int | key, | |||
Vio * | sock | |||
) |
Read in the user-specified initial mesh given in the "MCSF" or "MCEF" format, and transform into our internal datastructures.
Do a little more than a "Gem_read", in that we also initialize the extrinsic and intrinsic spatial dimensions corresponding to the input mesh, and we also then build the reference elements.
thee | Pointer to an Aprx allocated memory location | |
key | input format type key=0 ==> simplex format key=1 ==> edge format key=2 ==> simplex-nabor format | |
sock | socket for reading the external mesh data (NULL otherwise) |
int Aprx_refine | ( | Aprx * | thee, | |
int | rkey, | |||
int | bkey, | |||
int | pkey | |||
) |
Refine the mesh.
thee | Pointer to an Aprx allocated memory location | |
rkey | Boolean sets the rkey type for refining the mesh. Input: If (rkey==0) Perform recursive simplex bisection until conformity |
IMPORTANT NOTE: In 2D, (rkey==1) WILL generate a conforming mesh. However, in 3D, this procedure will in general produce nonconforming simplices. To produce a conforming mesh in 3D would require an implementation covering all possible face refinement combinations (something like 169 cases). This has been done e.g. by Jurgen Bey in AGM, but we are not that patient; use (rkey==0) above if you want a conforming mesh...
If (rkey==2) As a test of the conformity procedure, perform quadra-[octa-]-section until conformity, which should produce a uniformly regularly refined mesh. (In 2D, each triangle should be divided into four children, and in 3D each tetrahedron should be divided into eight children.)
bkey | Boolean sets the bkey type for bisecting the mesh If (bkey==0) Bisection type: Longest Edge If (bkey==1) Bisection type: Newest Vertex If (bkey==2) Bisection type: Newest Pair | |
pkey | Boolean sets the pkey type to prolongate a vector |
void Aprx_reset | ( | Aprx * | thee | ) |
Reset all of the datastructures.
thee | Pointer to an Aprx allocated memory location |
int Aprx_unRefine | ( | Aprx * | thee, | |
int | rkey, | |||
int | pkey | |||
) |
Unrefine the mesh.
thee | Pointer to an Aprx allocated memory location | |
rkey | Boolean sets the rkey type for un-refining the mesh. | |
pkey | Boolean sets the pkey type to prolongate a vector |
void Aprx_writeGEOM | ( | Aprx * | thee, | |
Vio * | sock, | |||
int | defKey, | |||
int | colKey, | |||
int | chartKey, | |||
double | gluVal, | |||
int | fkey, | |||
Bvec * | u, | |||
char * | format | |||
) |
Write a finite element mesh or mesh function in some format.
thee | Pointer to an Aprx allocated memory location | |
sock | socket for reading the external mesh data (NULL otherwise) | |
defKey | defKey == 0 ==> draw mesh as it is defKey == 1 ==> use "def??" as new vertex coords (deformation) defKey == 2 ==> add "def??" to old vertex coords (displacement) | |
colKey | colKey == 0 ==> color simplices all same default color colKey == 1 ==> color simplices based on their chart colKey == 2 ==> color boundary simplices based on type | |
chartKey | chartKey < 0 ==> draw all simplices chartKey >= 0 ==> draw only simplices with chart chartKey | |
gluVal | gluVal == 1. ==> draw all simplices glued together 0. < gluVal < 1. ==> draw simplices with some separation | |
fkey | fkey == 0 ==> draw simplices fkey == 1 ==> draw only simplex boundary faces fkey == 2 ==> draw only simplices with a boundary face | |
u | block NODAL solution vector | |
format | GV/MATH format |
Write a finite element mesh or mesh function in some format.
thee | Pointer to an Aprx allocated memory location | |
sock | socket for reading the external mesh data (NULL otherwise) | |
u | block NODAL solution vector | |
format | GV/MATH format |