Files | |
file | am.h |
Class AM: a multilevel finite element algebra object. | |
Classes | |
struct | sAM |
Class AM: Definition. More... | |
Typedefs | |
typedef struct sAM | AM |
Declaraction of the AM class as the AM structure. | |
Functions | |
void | AM_tSolve (AM *thee, int meth, double dt, double t0, int numstep, int pfreq, int efreq, int ekeytotal, DYNtype pdetype, double ltol, int lmax, SVio *vsock) |
Solution of Time-Dependent problems by Method of Lines. | |
AM * | AM_ctor (Vmem *vmem, Aprx *taprx) |
The AM constructor. | |
void | AM_dtor (AM **thee) |
The AM destructor. | |
void | AM_create (AM *thee) |
Create the following internal Alg structures: A ==> The tangent matrix (linearization operator) M ==> The mass matrix W[] ==> The node vectors. | |
void | AM_destroy (AM *thee) |
Destroy the following internal Alg structures: A ==> The tangent matrix (linearization operator) M ==> The mass matrix W[] ==> The node vectors. | |
int | AM_markRefine (AM *thee, int key, int color, int bkey, double elevel) |
Mark a given mesh for refinement. | |
int | AM_refine (AM *thee, int rkey, int bkey, int pkey) |
Refine the mesh. | |
int | AM_unRefine (AM *thee, int rkey, int pkey) |
Un-refine the mesh. | |
int | AM_deform (AM *thee) |
Deform the mesh. | |
int | AM_read (AM *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 "Aprx_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. | |
double | AM_assem (AM *thee, int evalKey, int energyKey, int residKey, int tangKey, int massKey, int bumpKey, Bvec *u, Bvec *ud, Bvec *f, int ip[], double rp[]) |
Assemble the linearized problem at a given level. | |
double | AM_evalJ (AM *thee) |
Assemble the energy functional at the current solution. | |
void | AM_evalFunc (AM *thee, int number, int block, int numPts, double *pts, double *vals, int *marks) |
Evaluate a finite element function. | |
void | AM_bndIntegral (AM *thee) |
Perform a boundary integral. | |
double | AM_evalError (AM *thee, int pcolor, int key) |
Evaluate error in the current solution. | |
void | AM_applyDiriZero (AM *thee, Bvec *which) |
Apply zero dirichlet condition at a given level. | |
void | AM_iniGuess (AM *thee, Bvec *which) |
Setup an initial guess at a given level. | |
int | AM_part (AM *thee, int pkey, int pwht, int ppow) |
Partition the mesh using the matching Alg level. | |
int | AM_partSet (AM *thee, int pcolor) |
Set the partition color. | |
int | AM_partSmooth (AM *thee) |
Do a partition smoothing. | |
void | AM_printJ (AM *thee) |
Print the energy. | |
void | AM_printA (AM *thee) |
Print the system matrix. | |
void | AM_printAnoD (AM *thee) |
Print the system matrix with Dirichlet rows/cols zeroed. | |
void | AM_printAsp (AM *thee, char *fname) |
Print the system matrix in MATLAB sparse format. | |
void | AM_printAspNoD (AM *thee, char *fname) |
Print the system matrix in MATLAB sparse format with Dirichlet rows/cols zeroed. | |
void | AM_printP (AM *thee) |
Print the prolongation matrix. | |
void | AM_printPsp (AM *thee, char *fname) |
Print the prolongation matrix in MATLAB sparse format. | |
void | AM_printV (AM *thee, int num) |
Print a specified vector. | |
void | AM_printVsp (AM *thee, int num, char *fname) |
Print a vector in MATLAB sparse format. | |
void | AM_writeGEOM (AM *thee, Vio *sock, int defKey, int colKey, int chartKey, double gluVal, int fkey, int number, char *format) |
Write out a mesh in some format. | |
void | AM_writeSOL (AM *thee, Vio *sock, int number, char *format) |
Write out a solution in some format. | |
void | AM_memChk (AM *thee) |
Print the exact current malloc usage. | |
void | AM_lSolve (AM *thee, int prob, int meth, int itmax, double etol, int prec, int gues, int pjac) |
Linear solver. | |
void | AM_hlSolve (AM *thee, int prob, int meth, int itmax, double etol, int prec, int gues, int pjac) |
Hierarchical linear solver. | |
void | AM_nSolve (AM *thee, int meth, int itmax, double etol, int lmeth, int litmax, double letol, int lprec, int gues, int pjac) |
void | AM_newton (AM *thee, int itmax, double etol, int lmeth, int litmax, double letol, int lprec, int pjac, double loadParm) |
Damped-Inexact Newton iteration. | |
VPUBLIC void | AM_homotopy (AM *thee, int itmax, double etol, int lmeth, int litmax, double letol, int lprec, int pjac) |
Homotopy and the related "incremental loading" iterations. |
Apply zero dirichlet condition at a given level.
thee | Pointer to class AM | |
which | the block vector |
double AM_assem | ( | AM * | thee, | |
int | evalKey, | |||
int | energyKey, | |||
int | residKey, | |||
int | tangKey, | |||
int | massKey, | |||
int | bumpKey, | |||
Bvec * | u, | |||
Bvec * | ud, | |||
Bvec * | f, | |||
int | ip[], | |||
double | rp[] | |||
) |
Assemble the linearized problem at a given level.
thee | Pointer to class AM | |
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. | |
u | block NODAL solution vector | |
ud | block NODAL dirichlet vector | |
f | block NODAL residual vector | |
ip | index for assembled energies | |
rp | parameter for initially assembled PDE |
void AM_bndIntegral | ( | AM * | thee | ) |
Perform a boundary integral.
thee | Pointer to Class AM |
void AM_create | ( | AM * | thee | ) |
Create the following internal Alg structures:
A ==> The tangent matrix (linearization operator)
M ==> The mass matrix
W[] ==> The node vectors.
thee | Pointer to Class AM |
The AM constructor.
vmem | Memory management object | |
taprx | Pointer to linear Approximation object |
int AM_deform | ( | AM * | thee | ) |
Deform the mesh.
thee | Pointer to class AM |
void AM_destroy | ( | AM * | thee | ) |
Destroy the following internal Alg structures:
A ==> The tangent matrix (linearization operator)
M ==> The mass matrix
W[] ==> The node vectors.
thee | Pointer to Class AM |
void AM_dtor | ( | AM ** | thee | ) |
The AM destructor.
thee | Pointer to Class AM |
double AM_evalError | ( | AM * | thee, | |
int | pcolor, | |||
int | key | |||
) |
Evaluate error in the current solution.
thee | Pointer to class AM | |
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. |
void AM_evalFunc | ( | AM * | thee, | |
int | number, | |||
int | block, | |||
int | numPts, | |||
double * | pts, | |||
double * | vals, | |||
int * | marks | |||
) |
Evaluate a finite element function.
thee | Pointer to class AM | |
number | the value of variable in the environment as an integer | |
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 |
double AM_evalJ | ( | AM * | thee | ) |
Assemble the energy functional at the current solution.
thee | Pointer to Class AM |
void AM_hlSolve | ( | AM * | thee, | |
int | prob, | |||
int | meth, | |||
int | itmax, | |||
double | etol, | |||
int | prec, | |||
int | gues, | |||
int | pjac | |||
) |
Hierarchical linear solver.
thee | Pointer to class AM | |
prob | index for primal or dual problem | |
meth | method choice (0=slu,1=mg,2=cg,3=bcg,4=pcg,5=pbcg) | |
itmax | number of iterations to do (the maximum allowed) | |
etol | error tolerance (currently ignored) | |
prec | index for different preconditioners | |
gues | index for initial guess | |
pjac | index for printing the system matrix |
VPUBLIC void AM_homotopy | ( | AM * | thee, | |
int | itmax, | |||
double | etol, | |||
int | lmeth, | |||
int | litmax, | |||
double | letol, | |||
int | lprec, | |||
int | pjac | |||
) |
Homotopy and the related "incremental loading" iterations.
thee | Pointer to class AM | |
itmax | number of iterations to do (the maximum allowed) | |
etol | error tolerance (currently ignored) | |
lmeth | index for specified linear solver | |
litmax | iteration max for linear solver | |
letol | error tolerance for linear solver | |
lprec | preconditioner for linear solver | |
pjac | index for printing the system matrix |
Setup an initial guess at a given level.
thee | Pointer to class AM | |
which | the block vector |
void AM_lSolve | ( | AM * | thee, | |
int | prob, | |||
int | meth, | |||
int | itmax, | |||
double | etol, | |||
int | prec, | |||
int | gues, | |||
int | pjac | |||
) |
Linear solver.
thee | Pointer to class AM | |
prob | index for primal or dual problem | |
meth | method choice (0=slu,1=mg,2=cg,3=bcg,4=pcg,5=pbcg) | |
itmax | number of iterations to do (the maximum allowed) | |
etol | error tolerance (currently ignored) | |
prec | index for different preconditioners | |
gues | index for initial guess | |
pjac | index for printing the system matrix |
int AM_markRefine | ( | AM * | thee, | |
int | key, | |||
int | color, | |||
int | bkey, | |||
double | elevel | |||
) |
Mark a given mesh for refinement.
thee | Pointer to class AM | |
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 |
void AM_memChk | ( | AM * | thee | ) |
Print the exact current malloc usage.
thee | Pointer to Class AM |
void AM_newton | ( | AM * | thee, | |
int | itmax, | |||
double | etol, | |||
int | lmeth, | |||
int | litmax, | |||
double | letol, | |||
int | lprec, | |||
int | pjac, | |||
double | loadParm | |||
) |
Damped-Inexact Newton iteration.
thee | Pointer to class AM | |
itmax | number of iterations to do (the maximum allowed) | |
etol | error tolerance (currently ignored) | |
lmeth | index for specified linear solver | |
litmax | iteration max for linear solver | |
letol | error tolerance for linear solver | |
lprec | preconditioner for linear solver | |
pjac | index for printing the system matrix | |
loadParm | the incoming load parameter |
void AM_nSolve | ( | AM * | thee, | |
int | meth, | |||
int | itmax, | |||
double | etol, | |||
int | lmeth, | |||
int | litmax, | |||
double | letol, | |||
int | lprec, | |||
int | gues, | |||
int | pjac | |||
) |
thee | Pointer to class AM | |
meth | method choice (0=slu,1=mg,2=cg,3=bcg,4=pcg,5=pbcg) | |
itmax | number of iterations to do (the maximum allowed) | |
etol | error tolerance (currently ignored) | |
lmeth | index for specified linear solver | |
litmax | iteration max for linear solver | |
letol | error tolerance for linear solver | |
lprec | preconditioner for linear solver | |
gues | index for initial guess | |
pjac | index for printing the system matrix |
int AM_part | ( | AM * | thee, | |
int | pkey, | |||
int | pwht, | |||
int | ppow | |||
) |
Partition the mesh using the matching Alg level.
thee | Pointer to Class AM | |
pkey | index for different partition option | |
pwht | index for weighted partitioning | |
ppow | partitioning steps |
int AM_partSet | ( | AM * | thee, | |
int | pcolor | |||
) |
Set the partition color.
thee | Pointer to class AM | |
pcolor | simplex chart type |
int AM_partSmooth | ( | AM * | thee | ) |
Do a partition smoothing.
thee | Pointer to Class AM |
void AM_printA | ( | AM * | thee | ) |
Print the system matrix.
thee | Pointer to Class AM |
void AM_printAnoD | ( | AM * | thee | ) |
Print the system matrix with Dirichlet rows/cols zeroed.
thee | Pointer to Class AM |
void AM_printAsp | ( | AM * | thee, | |
char * | fname | |||
) |
Print the system matrix in MATLAB sparse format.
thee | Pointer to class AM | |
fname | the system matrix file name |
void AM_printAspNoD | ( | AM * | thee, | |
char * | fname | |||
) |
Print the system matrix in MATLAB sparse format with Dirichlet rows/cols zeroed.
thee | Pointer to class AM | |
fname | the system matrix file name |
void AM_printJ | ( | AM * | thee | ) |
Print the energy.
thee | Pointer to Class AM |
void AM_printP | ( | AM * | thee | ) |
Print the prolongation matrix.
thee | Pointer to class AM |
void AM_printPsp | ( | AM * | thee, | |
char * | fname | |||
) |
Print the prolongation matrix in MATLAB sparse format.
thee | Pointer to class AM | |
fname | the system matrix file name |
void AM_printV | ( | AM * | thee, | |
int | num | |||
) |
Print a specified vector.
thee | Pointer to Class AM | |
num | the number of block work vector |
void AM_printVsp | ( | AM * | thee, | |
int | num, | |||
char * | fname | |||
) |
Print a vector in MATLAB sparse format.
thee | Pointer to class AM | |
num | the number of block work vector | |
fname | the output file name |
int AM_read | ( | AM * | 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 "Aprx_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 class AM | |
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 AM_refine | ( | AM * | thee, | |
int | rkey, | |||
int | bkey, | |||
int | pkey | |||
) |
Refine the mesh.
thee | Pointer to class AM | |
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 AM_tSolve | ( | AM * | thee, | |
int | meth, | |||
double | dt, | |||
double | t0, | |||
int | numstep, | |||
int | pfreq, | |||
int | efreq, | |||
int | ekeytotal, | |||
DYNtype | pdetype, | |||
double | ltol, | |||
int | lmax, | |||
SVio * | vsock | |||
) |
Solution of Time-Dependent problems by Method of Lines.
thee | Pointer to class AM | |
meth | method choice (0=Forward Euler;1=Backward Euler;2=Trapezoidal Rule) | |
dt | the time step | |
t0 | the initial time | |
numstep | number of the time steps | |
pfreq | frequency of the printing contour | |
efreq | frequency of energy printing | |
ekeytotal | total number of printing energies | |
pdetype | index for different types of PDEs:SIMP_TYPE/TDEP_TYPE/NLIN_TYPE | |
ltol | error tolerance | |
lmax | number of iterations to do (the maximum allowed) | |
vsock | socket for writing a finite element mesh or mesh function |
int AM_unRefine | ( | AM * | thee, | |
int | rkey, | |||
int | pkey | |||
) |
Un-refine the mesh.
thee | Pointer to class AM | |
rkey | Boolean sets the rkey type for un-refining the mesh. | |
pkey | Boolean sets the pkey type to prolongate a vector |
void AM_writeGEOM | ( | AM * | thee, | |
Vio * | sock, | |||
int | defKey, | |||
int | colKey, | |||
int | chartKey, | |||
double | gluVal, | |||
int | fkey, | |||
int | number, | |||
char * | format | |||
) |
Write out a mesh in some format.
thee | Pointer to class AM | |
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 | |
number | the number of block work vector | |
format | GV/MATH format |
void AM_writeSOL | ( | AM * | thee, | |
Vio * | sock, | |||
int | number, | |||
char * | format | |||
) |
Write out a solution in some format.
thee | Pointer to class AM | |
sock | socket for reading the external mesh data (NULL otherwise) | |
number | the number of block work vector | |
format | Pointer to GV/MATH format |