Files | |
file | bvec.h |
Class Bvec: a block vector object. | |
Classes | |
struct | sBvec |
Contains public data members for Bvec class. More... | |
Typedefs | |
typedef struct sBvec | Bvec |
Declaration of the Bvec class as the Bvec structure. | |
Functions | |
int | Bvec_numB (Bvec *thee) |
Return the number of blocks. | |
Bvec * | Bvec_ctor (Vmem *vmem, const char *name, int pnumB, int pnumR[MAXV]) |
The block vector constructor (data array malloc'd by us). | |
Bvec * | Bvec_ctor2 (Vmem *vmem, const char *name, int pnumB, int pnumR[MAXV], double *data) |
The block vector constructor (data array passed in). | |
Bvec * | Bvec_ctor3 (Vmem *vmem, const char *name, int length) |
The vector constructor (data array malloc'd by us). | |
Bvec * | Bvec_ctor4 (Vmem *vmem, const char *name, int length, double *data) |
The vector constructor (data array malloc'd by us). | |
void | Bvec_dtor (Bvec **thee) |
The block vector destructor. | |
void | Bvec_createVectors (Bvec *thee, Bvec *vecs[], int num) |
Create a set of vectors with duplicate attributes. | |
void | Bvec_destroyVectors (Bvec *thee, Bvec *vecs[], int num) |
Destroy set of vectors previously created by Bvec_createVectors. | |
void | Bvec_createVecMat (Bvec *thee, Bvec *vecs[], int num, Mat **mat) |
Create a matrix from columns of vectors. | |
void | Bvec_destroyVecMat (Bvec *thee, Bvec *vecs[], int num, Mat **mat) |
Destroy matrix previously created by Bvec_createVecMat. | |
int | Bvec_len (Bvec *thee) |
Return the total number of rows. | |
double * | Bvec_addr (Bvec *thee) |
Return a pointer to the total data array. | |
double | Bvec_val (Bvec *thee, int i) |
Return value of component of vector. | |
void | Bvec_set (Bvec *thee, int i, double val) |
Set value of component of vector. | |
double | Bvec_nrm1 (Bvec *thee) |
1-norm of a block vector. | |
double | Bvec_nrm2 (Bvec *thee) |
2-norm of a block vector. | |
double | Bvec_nrm8 (Bvec *thee) |
oo-norm of a block vector. | |
double | Bvec_dif1 (Bvec *thee, Bvec *v) |
1-norm of the difference of two block vectors. | |
double | Bvec_dif2 (Bvec *thee, Bvec *v) |
2-norm of the difference of two block vectors. | |
double | Bvec_dif8 (Bvec *thee, Bvec *v) |
oo-norm of the difference of two block vectors. | |
double | Bvec_dot (Bvec *thee, Bvec *v) |
Dot product of two block vectors. | |
void | Bvec_init (Bvec *thee, double val) |
Initialize a block vector. | |
void | Bvec_scal (Bvec *thee, double val) |
Scale a block vector. | |
void | Bvec_copy (Bvec *thee, Bvec *v) |
Copy a block vector. | |
void | Bvec_axpy (Bvec *thee, Bvec *v, double val) |
Saxpy for a block vector. | |
void | Bvec_print (Bvec *thee) |
Print the block vector. | |
void | Bvec_printSp (Bvec *thee, char *fname) |
Print the block vector. | |
void | Bvec_diagScale (Bvec *thee, Bmat *A, Bvec *f) |
Apply inverse of diagonal to a vector. | |
void | Bvec_matvec (Bvec *thee, Bmat *A, Bvec *v, int key, int part) |
Block matrix times vector. | |
void | Bvec_smooth (Bvec *thee, Bmat *A, Bvec *f, Bvec *w, int key, int ioflag, int meth, int adj, int itmax, double etol, double omega) |
Generic smoothing operator (for a general YSMP-based matrix). | |
void | Bvec_bnd (Bvec *thee, Bmat *bmat, int key) |
Apply boundary condition. | |
void | Bvec_memChk (Bvec *thee) |
Print the exact current malloc usage. | |
int | Bvec_numRT (Bvec *thee) |
Return the total number of rows. | |
int | Bvec_numRB (Bvec *thee, int i) |
Return the size of a block. | |
double * | Bvec_addrB (Bvec *thee, int i) |
Return a pointer to the data in a block. | |
double | Bvec_valB (Bvec *thee, int i, int which) |
Return value of component of a vector. | |
void | Bvec_setB (Bvec *thee, int i, int which, double val) |
Set value of component of a vector. | |
void | Bvec_addToB (Bvec *thee, int i, int which, double val) |
Add value to component of a vector. | |
void | Bvec_initB (Bvec *thee, int i, double val) |
Initialize a particular block of a block vector. | |
void | Bvec_absLog (Bvec *thee, double val) |
Log of abs value of entries of block vector, shifted by < val >. | |
void | Bvec_lmethod (Bvec *thee, Bmat *A, Bvec *f, Bvec *r, Bvec *ut, int key, int flag, int itmax, double etol, int prec, int cycle, Bmat *P, int meth) |
Executes one of a number of linear solvers. | |
void | Bvec_slu (Bvec *thee, Bmat *A, Bvec *f, Bvec *r, Bvec *ut, int key, int flag) |
Sparse LU direct solver (requires only a nonsingular A). | |
void | Bvec_mg (Bvec *thee, Bmat *A, Bvec *f, Bvec *r, Bvec *ut, int key, int flag, int itmax, double etol, int prec, int cycle, Bmat *P) |
MG algorithm front end (NOT recursively callable). | |
void | Bvec_mgInit (Bvec *thee, Bmat *A, Bvec *f, Bvec *r, Bmat *P) |
Initialize the algebraic multilevel hierarchy for Bvec_mg. | |
void | Bvec_mgDestroy (Bvec *thee, Bmat *A, Bvec *f, Bvec *r) |
Destroy algebraic multilevel hierarchy created by Bvec_mgInit. | |
void | Bvec_cg (Bvec *thee, Bmat *A, Bvec *f, Bvec *r, Bvec *ut, int key, int flag, int itmax, double etol, int prec, int cycle, Bmat *P, Bvec *p, Bvec *ap, Bvec *bap, Bvec *po, Bvec *apo, Bvec *tp) |
CG method (ODIR algorithm, requires SPD operators A, B, and M). | |
void | Bvec_bcg (Bvec *thee, Bmat *A, Bvec *f, Bvec *r, Bvec *ut, int key, int flag, int itmax, double etol, int prec, int cycle, Bmat *P, Bvec *s, Bvec *p, Bvec *ap, Bvec *bap, Bvec *po, Bvec *apo, Bvec *q, Bvec *atq, Bvec *btatq, Bvec *qo, Bvec *atqo) |
BCG method (pure bi-orthogonal extension of CG). | |
void | Bvec_eig (Bvec *thee, Bmat *A, int litmax, double letol, double *lambda, int key, int flag, int itmax, double etol) |
Calculate the eigenvector corresponding to the second smallest eigenvalue of the given (symmetric adjacency) matrix. | |
void | Bvec_hlmethod (Bvec *thee, Bmat *A, Bvec *f, Bvec *r, Bvec *ut, int key, int flag, int itmax, double etol, int prec, int cycle, Bmat *P, Bmat *M, int meth) |
Executes one of a number of hierarchical linear solvers. | |
void | Bvec_hb (Bvec *thee, Bmat *A, Bvec *f, Bvec *r, Bvec *ut, int key, int flag, int itmax, double etol, int meth, Bmat *Ppro, Bmat *M) |
HBMG algorithm due to Bank, Dupont, and Yserentant. | |
void | Bvec_hcg (Bvec *thee, Bmat *A, Bvec *f, Bvec *r, Bvec *ut, int key, int flag, int itmax, double etol, int prec, Bmat *P, Bmat *M, Bvec *p, Bvec *ap, Bvec *bap, Bvec *po, Bvec *apo, Bvec *tp) |
CG method with Bvec_hb as an optional preconditioner. | |
void | Bvec_hbcg (Bvec *thee, Bmat *A, Bvec *f, Bvec *r, Bvec *ut, int key, int flag, int itmax, double etol, int prec, Bmat *P, Bmat *M, Bvec *s, Bvec *p, Bvec *ap, Bvec *bap, Bvec *po, Bvec *apo, Bvec *q, Bvec *atq, Bvec *btatq, Bvec *qo, Bvec *atqo) |
Bvec_hbcg with Bvec_hb as an optional preconditioner. | |
Bvec * | Bvec_ctorPoint (Bvec *v, const char *name, int *ibase, int *numR) |
Create a Bvec pointer to the tail of a Bvec. | |
void | Bvec_dtorPoint (Bvec **thee) |
The block vec pointer destructor. | |
void | Bvec_submethod (Bvec *thee, Bmat *A, Bvec *f, Bvec *r, Bvec *ut, int key, int flag, int itmax, double etol, int prec, int cycle, Bmat *P, int meth) |
Executes one of a number of linear solvers. | |
void | Bvec_fsmooth (Bvec *thee, Bmat *amat, Bvec *f, Bvec *w, int key, int ioflag, int meth, int adj, int itmax, Bchar *fc) |
F-point smoothing operator. |
void Bvec_absLog | ( | Bvec * | thee, | |
double | val | |||
) |
Log of abs value of entries of block vector, shifted by < val >.
thee | pointer to the block vector | |
val | value of the i-th component of the block vector |
double* Bvec_addr | ( | Bvec * | thee | ) |
Return a pointer to the total data array.
thee | pointer to the block vector |
double* Bvec_addrB | ( | Bvec * | thee, | |
int | i | |||
) |
Return a pointer to the data in a block.
thee | pointer to the block vector | |
i | index of the block vector |
void Bvec_addToB | ( | Bvec * | thee, | |
int | i, | |||
int | which, | |||
double | val | |||
) |
Add value to component of a vector.
thee | pointer to the block vector | |
i | index of the block vector | |
which | index of the component in the vector | |
val | value of the i-th component of the block vector |
Saxpy for a block vector.
thee | pointer to the block vector | |
v | the source block vector | |
val | coeficient for scaling the source block vector |
void Bvec_bcg | ( | Bvec * | thee, | |
Bmat * | A, | |||
Bvec * | f, | |||
Bvec * | r, | |||
Bvec * | ut, | |||
int | key, | |||
int | flag, | |||
int | itmax, | |||
double | etol, | |||
int | prec, | |||
int | cycle, | |||
Bmat * | P, | |||
Bvec * | s, | |||
Bvec * | p, | |||
Bvec * | ap, | |||
Bvec * | bap, | |||
Bvec * | po, | |||
Bvec * | apo, | |||
Bvec * | q, | |||
Bvec * | atq, | |||
Bvec * | btatq, | |||
Bvec * | qo, | |||
Bvec * | atqo | |||
) |
BCG method (pure bi-orthogonal extension of CG).
* PG-ODIR(M,B,A) (for Au=f, preconditioner B, inner-product M) * ------------------------------------------------------------ * * Given u^0 * r^0 = f - A u^0 * s^0 = f - A' u^0 * p^0 = B r^0 / ||B r^0|| * q^0 = B's^0 / ||B's^0|| * for i=0,1,2,.... * alpha_i = <M e^i, q^i> / <M p^i, q^i> * u^{i+1} = u^i + alpha_i p^i * r^{i+1} = r^i - alpha_i A p^i * gamma_i = <M B A p^i, q^i> / <M p^i, q^i> * sigma_i = <M B A p^i, q^{i-1}> / <M p^{i-1}, q^{i-1}> * delta_i = <M B A p^{i-1}, q^i> / <M p^{i-1}, q^{i-1}> * v = B A p^i - gamma_i p^i - sigma_i p^{i-1} * w = B'A' q^i - gamma_i q^i - delta_i q^{i-1} * p^{i+1} = v / ||v|| * q^{i+1} = w / ||w|| * end for * * BCG=PG-ODIR(M=A,B=B,A=A) * ------------------------ * Given u^0 * r^0 = f - A u^0 * s^0 = f - A' u^0 * p^0 = B r^0 / ||B r^0|| * q^0 = B's^0 / ||B's^0|| * for i=0,1,2,.... * alpha_i = <r^i, q^i> / <A p^i, q^i> * u^{i+1} = u^i + alpha_i p^i * r^{i+1} = r^i - alpha_i A p^i * gamma_i = <B A p^i, A' q^i> / <A p^i, q^i> * sigma_i = <B A p^i, A' q^{i-1}> / <A p^{i-1}, q^{i-1}> * delta_i = <A p^{i-1}, B'A' q^i> / <A p^{i-1}, q^{i-1}> * v = B A p^i - gamma_i p^i - sigma_i p^{i-1} * w = B'A' q^i - gamma_i q^i - delta_i q^{i-1} * p^{i+1} = v / ||v|| * q^{i+1} = w / ||w|| * end for *
thee | pointer to the block vector | |
A | system matrix | |
f | source vector slot | |
r | residual slot | |
ut | block NODAL analytical solution | |
key | smooth with A or A' (0=A, 1=A') | |
flag | determines which "mode" we run in: 0 --> Normal: check itmax, error tolerance; normal i/o 1 --> Silent: check only itmax; no i/o 2 --> Subcycle: check itmax, error tolerance; subcycle i/o 3 --> Subcycle: check itmax, error tolerance; no i/o | |
itmax | number of iterations to do (the maximum allowed) | |
etol | error tolerance | |
prec | index for different preconditioners | |
cycle | index for symmetric/nonsymmetric multigrids | |
P | prolongation matrix maintained by aprx | |
s | residual block vector s: f - A'u | |
p | the block vector p | |
ap | another block vector ap | |
bap | third block vector bap | |
po | the block vector po | |
apo | the block vector apo | |
q | the block vector q | |
atq | the block vector atq | |
btatq | the block vector btatq | |
qo | the block vector qo | |
atqo | the block vector atqo |
Apply boundary condition.
thee | pointer to the block vector | |
bmat | the block sparse matrix | |
key | index for setting zeroes in the row/column components. |
void Bvec_cg | ( | Bvec * | thee, | |
Bmat * | A, | |||
Bvec * | f, | |||
Bvec * | r, | |||
Bvec * | ut, | |||
int | key, | |||
int | flag, | |||
int | itmax, | |||
double | etol, | |||
int | prec, | |||
int | cycle, | |||
Bmat * | P, | |||
Bvec * | p, | |||
Bvec * | ap, | |||
Bvec * | bap, | |||
Bvec * | po, | |||
Bvec * | apo, | |||
Bvec * | tp | |||
) |
CG method (ODIR algorithm, requires SPD operators A, B, and M).
* ODIR(M,B,A) (for Au=f, preconditioner B, inner-product M) * --------------------------------------------------------- * * Given u^0 * r^0 = f - A u^0 * p^0 = B r^0 / ||B r^0|| * for i=0,1,2,.... * alpha_i = <M e^i, p^i> / <M p^i, p^i> * u^{i+1} = u^i + alpha_i p^i * r^{i+1} = r^i - alpha_i A p^i * gamma_i = <M B A p^i, p^i> / <M p^i, p^i> * sigma_i = <M B A p^i, p^{i-1}> / <M p^{i-1}, p^{i-1}> * v = B A p^i - gamma_i p^i - sigma_i p^{i-1} * p^{i+1} = v / ||v|| * end for * * CG=ODIR(M=A,B=B,A=A) * -------------------- * * Given u^0 * r^0 = f - A u^0 * p^0 = B r^0 / ||B r^0|| * for i=0,1,2,.... * alpha_i = <r^i, p^i> / <A p^i, p^i> * u^{i+1} = u^i + alpha_i p^i * r^{i+1} = r^i - alpha_i A p^i * gamma_i = <B A p^i, A p^i> / <A p^i, p^i> * sigma_i = <B A p^i, A p^{i-1}> / <A p^{i-1}, p^{i-1}> * v = B A p^i - gamma_i p^i - sigma_i p^{i-1} * p^{i+1} = v / ||v|| * end for *
thee | pointer to the block vector | |
A | system matrix | |
f | source vector slot | |
r | residual slot | |
ut | block NODAL analytical solution | |
key | index for different solution methods: 0=(Au=f), 1=(A'u=f) | |
flag | 0 --> Normal: check itmax, error tolerance; normal i/o 1 --> Silent: check only itmax; no i/o 2 --> Subcycle: check itmax, error tolerance; subcycle i/o 3 --> Subcycle: check itmax, error tolerance; no i/o | |
itmax | number of iterations to do (the maximum allowed) | |
etol | error tolerance | |
prec | index for different preconditioners | |
cycle | index for symmetric/nonsymmetric multigrids | |
P | prolongation matrix maintained by aprx | |
p | the block vector p | |
ap | another block vector ap | |
bap | third block vector bap | |
po | the block vector po | |
apo | the block vector apo | |
tp | temporary block vector tp |
Copy a block vector.
thee | pointer to the block vector | |
v | the source block vector |
Create a matrix from columns of vectors.
thee | pointer to the block vector | |
vecs | the created block vector | |
num | the number of created vectors | |
mat | Pointer to the sparse matrix |
Create a set of vectors with duplicate attributes.
thee | pointer to the block vector | |
vecs | the created block vector | |
num | the number of created vectors |
Bvec* Bvec_ctor | ( | Vmem * | vmem, | |
const char * | name, | |||
int | pnumB, | |||
int | pnumR[MAXV] | |||
) |
The block vector constructor (data array malloc'd by us).
vmem | Memory management object | |
name | character string name for this vector | |
pnumB | num vector blocks | |
pnumR | num of rows in each block vector |
Bvec* Bvec_ctor2 | ( | Vmem * | vmem, | |
const char * | name, | |||
int | pnumB, | |||
int | pnumR[MAXV], | |||
double * | data | |||
) |
The block vector constructor (data array passed in).
vmem | Memory management object | |
name | character string name for this vector | |
pnumB | num vector blocks | |
pnumR | num of rows in each block vector | |
data | the vector components |
Bvec* Bvec_ctor3 | ( | Vmem * | vmem, | |
const char * | name, | |||
int | length | |||
) |
The vector constructor (data array malloc'd by us).
vmem | Memory management object | |
name | character string name for this vector | |
length | number of components in the vector |
Bvec* Bvec_ctor4 | ( | Vmem * | vmem, | |
const char * | name, | |||
int | length, | |||
double * | data | |||
) |
The vector constructor (data array malloc'd by us).
vmem | Memory management object | |
name | character string name for this vector | |
length | number of components in the vector | |
data | the vector components |
Create a Bvec pointer to the tail of a Bvec.
v | Pointer to the Bvec object | |
name | character string name for the block vector | |
ibase | pointer to size array | |
numR | num of rows in each block vector |
Destroy matrix previously created by Bvec_createVecMat.
thee | pointer to the block vector | |
vecs | the created block vector | |
num | the number of created vectors | |
mat | Pointer to the sparse matrix |
Destroy set of vectors previously created by Bvec_createVectors.
thee | pointer to the block vector | |
vecs | the created block vector | |
num | the number of created vectors |
Apply inverse of diagonal to a vector.
thee | pointer to the block vector | |
A | system matrix | |
f | source vector slot |
1-norm of the difference of two block vectors.
thee | pointer to the block vector | |
v | the source block vector |
2-norm of the difference of two block vectors.
thee | pointer to the block vector | |
v | the source block vector |
oo-norm of the difference of two block vectors.
thee | pointer to the block vector | |
v | the source block vector |
Dot product of two block vectors.
thee | pointer to the block vector | |
v | the source block vector |
void Bvec_dtor | ( | Bvec ** | thee | ) |
The block vector destructor.
thee | pointer to the block vector |
void Bvec_dtorPoint | ( | Bvec ** | thee | ) |
The block vec pointer destructor.
thee | Pointer to Pointer of the Bvec object |
void Bvec_eig | ( | Bvec * | thee, | |
Bmat * | A, | |||
int | litmax, | |||
double | letol, | |||
double * | lambda, | |||
int | key, | |||
int | flag, | |||
int | itmax, | |||
double | etol | |||
) |
Calculate the eigenvector corresponding to the second smallest eigenvalue of the given (symmetric adjacency) matrix.
thee | pointer to the block vector | |
A | system matrix | |
litmax | number of iterations to do (the maximum allowed) | |
letol | error tolerance | |
lambda | Point to the eigenvalue | |
key | smooth with A or A' (0=A, 1=A') | |
flag | determines which "mode" we run in: 0 --> Normal: check itmax, error tolerance; normal i/o 1 --> Silent: check only itmax; no i/o 2 --> Subcycle: check itmax, error tolerance; subcycle i/o 3 --> Subcycle: check itmax, error tolerance; no i/o | |
itmax | number of iterations to do (the maximum allowed) | |
etol | error tolerance |
void Bvec_fsmooth | ( | Bvec * | thee, | |
Bmat * | amat, | |||
Bvec * | f, | |||
Bvec * | w, | |||
int | key, | |||
int | ioflag, | |||
int | meth, | |||
int | adj, | |||
int | itmax, | |||
Bchar * | fc | |||
) |
F-point smoothing operator.
* Notes: Consider the following partitioned linear system * * [ Acc Acf ] [ xc ] = [ bc ] * [ Afc Aff ] [ xf ] [ bf ] * * where xc is a vector of "coarse" or c-points and * xf is a vector of "fine" or f-points. This partitioning * induces an analogous partitioning of the matrix A and the * vector b as illustrated above. Given such a partitioning, * an f-point smoother applies a Gauss-Seidel (or Jacobi) * smoother to the f-points, leaving the c-points unchanged. * * The forward Gauss-Seidel f-point smoother can be written as * xc = xc * xf = ( Dff + Lff )^(-1) ( bf - Afc xc - Uff xf ) * * The adjoint Gauss-Seidel f-point smoother can be written as * xc = xc * xf = ( Dff + Uff )^(-1) ( bf - Afc xc - Lff xf ) * * where, Aff = Uff + Lff + Dff is the standard gs splitting. * * This method is invoked in the same manner as Bvec_smooth * with one additional input. This additional input is a char * array of point types, with coarse and fine points indicated * by 'c' and 'f' respectively. * * For example, suppose we have 3 blocks with fpoints (1,5) * in the first block, (0,3) in the second block, and (2,3,4) * in the third block. In this case fc should be * * fc = {{'c', 'f', 'c', 'c', 'c', 'f'}, * {'f', 'c', 'c', 'f', 'c', 'c'}, * {'c', 'c', 'f', 'f', 'f', 'c'}} * * This algorithm is not fast and its only practical purpose is * for understanding and debugging HB and BPX style preconditioners. * The HBMG algorithm is output-equivalent to standard * multigrid with f-point smoothing where the f-points are the * degrees of freedom introduced on the current level. For * multiplicative BPX, the list of f-points is expanded to include * the one-ring of the f-points used by HB. * * The optimal versions of multiplicative HB and BPX require * the use of a change of basis so the f-point smoother does * not need to touch the c-points explicitly. However, the * optimal version is mathematically "output-equivalent" to * the suboptimal version described above. *
thee | pointer to the block vector | |
amat | system matrix | |
f | source vector slot | |
w | the work block vector | |
key | smooth with A or A' (0=A, 1=A') | |
ioflag | debug output level (0=normal, 1=none, 2=lots, ... ) | |
meth | smoother choice (0=jac, 1=gs, ... ) | |
adj | adjoint choice (0=normal, 1=adjoint) | |
itmax | number of iterations to do (the maximum allowed) | |
fc | Pointer to the block character |
void Bvec_hb | ( | Bvec * | thee, | |
Bmat * | A, | |||
Bvec * | f, | |||
Bvec * | r, | |||
Bvec * | ut, | |||
int | key, | |||
int | flag, | |||
int | itmax, | |||
double | etol, | |||
int | meth, | |||
Bmat * | Ppro, | |||
Bmat * | M | |||
) |
HBMG algorithm due to Bank, Dupont, and Yserentant.
thee | pointer to the block vector | |
A | system matrix | |
f | source vector slot | |
r | residual slot | |
ut | block NODAL analytical solution | |
key | 0 --> solve Au=f 1 --> solve A'u=f | |
flag | 0 --> Normal: check itmax, error tolerance; normal i/o 1 --> Normal: check itmax; no error tolerance; no i/o | |
itmax | number of iterations to do (the maximum allowed) | |
etol | error tolerance | |
meth | method choice (0=Standard HB,1=Wavelet Modified HB) | |
Ppro | tails of the prolongation matrices (inside Algs) 3 --> Standard HB 5 --> Wavelet Modified HB | |
M | the mass matrix |
void Bvec_hbcg | ( | Bvec * | thee, | |
Bmat * | A, | |||
Bvec * | f, | |||
Bvec * | r, | |||
Bvec * | ut, | |||
int | key, | |||
int | flag, | |||
int | itmax, | |||
double | etol, | |||
int | prec, | |||
Bmat * | P, | |||
Bmat * | M, | |||
Bvec * | s, | |||
Bvec * | p, | |||
Bvec * | ap, | |||
Bvec * | bap, | |||
Bvec * | po, | |||
Bvec * | apo, | |||
Bvec * | q, | |||
Bvec * | atq, | |||
Bvec * | btatq, | |||
Bvec * | qo, | |||
Bvec * | atqo | |||
) |
Bvec_hbcg with Bvec_hb as an optional preconditioner.
thee | pointer to the block vector | |
A | system matrix | |
f | source vector slot | |
r | residual slot | |
ut | block NODAL analytical solution | |
key | smooth with A or A' (0=A, 1=A') | |
flag | Determines which "mode" we run in: 0 --> Normal: check itmax, error tolerance; normal i/o 1 --> Silent: check only itmax; no i/o 2 --> Subcycle: check itmax, error tolerance; subcycle i/o 3 --> Subcycle: check itmax, error tolerance; no i/o | |
itmax | number of iterations to do (the maximum allowed) | |
etol | error tolerance | |
prec | index for different preconditioners | |
P | prolongation matrix maintained by aprx | |
M | the mass matrix | |
s | residual block vector s: f - A'u | |
p | the block vector p | |
ap | another block vector ap | |
bap | third block vector bap | |
po | the block vector po | |
apo | the block vector apo | |
q | the block vector q | |
atq | the block vector atq | |
btatq | the block vector btatq | |
qo | the block vector qo | |
atqo | the block vector atqo |
void Bvec_hcg | ( | Bvec * | thee, | |
Bmat * | A, | |||
Bvec * | f, | |||
Bvec * | r, | |||
Bvec * | ut, | |||
int | key, | |||
int | flag, | |||
int | itmax, | |||
double | etol, | |||
int | prec, | |||
Bmat * | P, | |||
Bmat * | M, | |||
Bvec * | p, | |||
Bvec * | ap, | |||
Bvec * | bap, | |||
Bvec * | po, | |||
Bvec * | apo, | |||
Bvec * | tp | |||
) |
CG method with Bvec_hb as an optional preconditioner.
thee | pointer to the block vector | |
A | system matrix | |
f | source vector slot | |
r | residual slot | |
ut | block NODAL analytical solution | |
key | index for different solution methods: 0=(Au=f), 1=(A'u=f) | |
flag | Determines which "mode" we run in: 0 --> Normal: check itmax, error tolerance; normal i/o 1 --> Silent: check only itmax; no i/o 2 --> Subcycle: check itmax, error tolerance; subcycle i/o 3 --> Subcycle: check itmax, error tolerance; no i/o | |
itmax | number of iterations to do (the maximum allowed) | |
etol | error tolerance | |
prec | index for different preconditioners | |
P | prolongation matrix maintained by aprx | |
M | the mass matrix | |
p | the block vector p | |
ap | another block vector ap | |
bap | third block vector bap | |
po | the block vector po | |
apo | the block vector apo | |
tp | temporary block vector tp |
void Bvec_hlmethod | ( | Bvec * | thee, | |
Bmat * | A, | |||
Bvec * | f, | |||
Bvec * | r, | |||
Bvec * | ut, | |||
int | key, | |||
int | flag, | |||
int | itmax, | |||
double | etol, | |||
int | prec, | |||
int | cycle, | |||
Bmat * | P, | |||
Bmat * | M, | |||
int | meth | |||
) |
Executes one of a number of hierarchical linear solvers.
thee | pointer to the block vector | |
A | system matrix | |
f | source vector slot | |
r | residual slot | |
ut | block NODAL analytical solution | |
key | index for different solution methods: 0=(Au=f), 1=(A'u=f) | |
flag | index for i/o control | |
itmax | number of iterations to do (the maximum allowed) | |
etol | error tolerance | |
prec | index for different preconditioners | |
cycle | index for symmetric/nonsymmetric multigrids | |
P | prolongation matrix maintained by aprx | |
M | the mass matrix | |
meth | method choice (0=Standard HB,1=Wavelet Modified HB) |
void Bvec_init | ( | Bvec * | thee, | |
double | val | |||
) |
Initialize a block vector.
thee | pointer to the block vector | |
val | value to be initialized for the block vector |
void Bvec_initB | ( | Bvec * | thee, | |
int | i, | |||
double | val | |||
) |
Initialize a particular block of a block vector.
thee | pointer to the block vector | |
i | index of the block vector | |
val | value of the i-th component of the block vector |
int Bvec_len | ( | Bvec * | thee | ) |
Return the total number of rows.
thee | pointer to the block vector |
void Bvec_lmethod | ( | Bvec * | thee, | |
Bmat * | A, | |||
Bvec * | f, | |||
Bvec * | r, | |||
Bvec * | ut, | |||
int | key, | |||
int | flag, | |||
int | itmax, | |||
double | etol, | |||
int | prec, | |||
int | cycle, | |||
Bmat * | P, | |||
int | meth | |||
) |
Executes one of a number of linear solvers.
thee | pointer to the block vector | |
A | system matrix | |
f | source vector slot | |
r | residual slot | |
ut | block NODAL analytical solution | |
key | index for different solution methods: 0=(Au=f), 1=(A'u=f) | |
flag | index for i/o control | |
itmax | number of iterations to do (the maximum allowed) | |
etol | error tolerance | |
prec | index for different preconditioners | |
cycle | index for symmetric/nonsymmetric multigrids | |
P | prolongation matrix maintained by aprx | |
meth | method choice (0=slu,1=mg,2=cg,3=bcg,4=pcg,5=pbcg) |
Block matrix times vector.
thee | pointer to the block vector | |
A | system matrix | |
v | the source block vector | |
key | index for different accumulation options | |
part | index for future use |
void Bvec_memChk | ( | Bvec * | thee | ) |
Print the exact current malloc usage.
thee | pointer to the block vector |
void Bvec_mg | ( | Bvec * | thee, | |
Bmat * | A, | |||
Bvec * | f, | |||
Bvec * | r, | |||
Bvec * | ut, | |||
int | key, | |||
int | flag, | |||
int | itmax, | |||
double | etol, | |||
int | prec, | |||
int | cycle, | |||
Bmat * | P | |||
) |
MG algorithm front end (NOT recursively callable).
thee | solution slot | |
A | system matrix | |
f | source vector slot | |
r | residual slot | |
ut | true solution slot | |
key | just passed along to Bvec_mgCore() | |
flag | just passed along to Bvec_mgCore() | |
itmax | just passed along to Bvec_mgCore() | |
etol | just passed along to Bvec_mgCore() | |
prec | just passed along to Bvec_mgCore() | |
cycle | just passed along to Bvec_mgCore() | |
P | prolongation operators (linked list) |
Destroy algebraic multilevel hierarchy created by Bvec_mgInit.
thee | solution slot | |
A | system matrix | |
f | source vector slot | |
r | residual slot |
Initialize the algebraic multilevel hierarchy for Bvec_mg.
thee | solution slot | |
A | system matrix | |
f | source vector slot | |
r | residual slot | |
P | prolongation operators (linked list) |
double Bvec_nrm1 | ( | Bvec * | thee | ) |
1-norm of a block vector.
thee | pointer to the block vector |
double Bvec_nrm2 | ( | Bvec * | thee | ) |
2-norm of a block vector.
thee | pointer to the block vector |
double Bvec_nrm8 | ( | Bvec * | thee | ) |
oo-norm of a block vector.
thee | pointer to the block vector |
int Bvec_numB | ( | Bvec * | thee | ) |
Return the number of blocks.
thee | pointer to the block vector |
int Bvec_numRB | ( | Bvec * | thee, | |
int | i | |||
) |
Return the size of a block.
thee | pointer to the block vector | |
i | index of the block vector |
int Bvec_numRT | ( | Bvec * | thee | ) |
Return the total number of rows.
thee | pointer to the block vector |
void Bvec_print | ( | Bvec * | thee | ) |
Print the block vector.
thee | pointer to the block vector |
void Bvec_printSp | ( | Bvec * | thee, | |
char * | fname | |||
) |
Print the block vector.
thee | pointer to the block vector | |
fname | the file name of the block vector |
void Bvec_scal | ( | Bvec * | thee, | |
double | val | |||
) |
Scale a block vector.
thee | pointer to the block vector | |
val | value of the i-th component of the block vector |
void Bvec_set | ( | Bvec * | thee, | |
int | i, | |||
double | val | |||
) |
Set value of component of vector.
thee | pointer to the block vector | |
i | index of the block vector | |
val | value of the i-th component of the block vector |
void Bvec_setB | ( | Bvec * | thee, | |
int | i, | |||
int | which, | |||
double | val | |||
) |
Set value of component of a vector.
thee | pointer to the block vector | |
i | index of the block vector | |
which | index of the component in the vector | |
val | value of the i-th component of the block vector |
Sparse LU direct solver (requires only a nonsingular A).
thee | pointer to the block vector | |
A | system matrix | |
f | source vector slot | |
r | residual slot | |
ut | block NODAL analytical solution | |
key | key ==0 --> Solve: A u = f key ==1 --> Solve: A' u = f | |
flag | Determines which "mode" we run in: flag==0,2 --> Normal: normal i/o flag==1,3 --> Silent: no i/o |
void Bvec_smooth | ( | Bvec * | thee, | |
Bmat * | A, | |||
Bvec * | f, | |||
Bvec * | w, | |||
int | key, | |||
int | ioflag, | |||
int | meth, | |||
int | adj, | |||
int | itmax, | |||
double | etol, | |||
double | omega | |||
) |
Generic smoothing operator (for a general YSMP-based matrix).
thee | pointer to the block vector | |
A | system matrix | |
f | source vector slot | |
w | the work block vector | |
key | smooth with A or A' (0=A, 1=A') | |
ioflag | debug output level (0=normal, 1=none, 2=lots, ... ) | |
meth | smoother choice (0=jac, 1=gs, ... ) | |
adj | adjoint choice (0=normal, 1=adjoint) | |
itmax | number of iterations to do (the maximum allowed) | |
etol | error tolerance (currently ignored) | |
omega | SOR parameter (currently ignored) |
void Bvec_submethod | ( | Bvec * | thee, | |
Bmat * | A, | |||
Bvec * | f, | |||
Bvec * | r, | |||
Bvec * | ut, | |||
int | key, | |||
int | flag, | |||
int | itmax, | |||
double | etol, | |||
int | prec, | |||
int | cycle, | |||
Bmat * | P, | |||
int | meth | |||
) |
Executes one of a number of linear solvers.
thee | pointer to the block vector | |
A | system matrix | |
f | source vector slot | |
r | residual slot | |
ut | block NODAL analytical solution | |
key | index for different solution methods: 0=(Au=f), 1=(A'u=f) | |
flag | index for i/o control | |
itmax | number of iterations to do (the maximum allowed) | |
etol | error tolerance | |
prec | index for different preconditioners | |
cycle | index for symmetric/nonsymmetric multigrids | |
P | prolongation matrix maintained by aprx | |
meth | method choice (0=Standard HB,1=Wavelet Modified HB) |
double Bvec_val | ( | Bvec * | thee, | |
int | i | |||
) |
Return value of component of vector.
thee | pointer to the block vector | |
i | index of the block vector |
double Bvec_valB | ( | Bvec * | thee, | |
int | i, | |||
int | which | |||
) |
Return value of component of a vector.
thee | pointer to the block vector | |
i | index of the block vector | |
which | index of the component in the vector |