#include <mc/mc_base.h>
#include <mc/bam.h>
Go to the source code of this file.
Classes | |
struct | sHBmat |
struct HBmat Definition More... | |
struct | sHBvec |
struct HBvec Definition More... | |
struct | sBchar |
struct Bchar Definition More... | |
Typedefs | |
typedef struct sHBmat | HBmat |
Declaration of the HBmat class as the HBmat structure. | |
typedef struct sHBvec | HBvec |
Declaration of the HBvec class as the HBvec structure. | |
typedef struct sBchar | Bchar |
Declaration of the Bchar class as the Bchar structure. | |
Enumerations | |
enum | HBMATtype { ZERO_TYPE, AMIN_TYPE, ANOR_TYPE, GHB_TYPE, GWM_TYPE } |
Class Whb: Parameters and datatypes. | |
Functions | |
void | HBvec_hbVcyc (HBvec *dd, HBmat *Ahb, HBvec *rr, HBmat *Ghb, HBvec *ww, int key, int csolv) |
Hierarchical Basis MG Vcycle (Bank, Dupont, and Yserentant). | |
void | HBvec_initStructure (HBvec *thee, HBmat *Ahb) |
Creates a chain of HBvec pointers to a Bvec. | |
void | HBvec_killStructure (HBvec *thee) |
Kills a chain of HBvec pointers. | |
void | HBmat_initG (HBmat *thee, Bmat *Ppro, Bmat *Mlink, int meth) |
Assembles the change of basis matrices for the HB Vcycle. | |
void | HBmat_initA (HBmat *Ahb, HBmat *Ghb, Bmat *Alink) |
Assembles the hierarchical matrix for the HB Vcycle. | |
void | HBmat_killStructure (HBmat *thee) |
Kills a chain of HBmat pointers. | |
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. | |
void | HBmat_initMulti (HBmat **Ahb, HBmat **Ghb, Bmat *A, Bmat *M, Bmat *Ppro, int meth) |
Assembles the hierarchical matrix structures for the HB Vcycle. | |
void | HBmat_killMulti (HBmat **Ahb, HBmat **Ghb) |
Destruct the hierarchical structures. | |
HBmat * | HBmat_ctor (Vmem *vmem, const char *bname, int pnumB) |
The HBmat constructor. | |
void | HBmat_dtor (HBmat **thee) |
The HBmat destructor. | |
HBvec * | HBvec_ctor (Vmem *vmem, const char *bname, int pnumB) |
The HBvec constructor. | |
void | HBvec_dtor (HBvec **thee) |
The HBvec destructor. | |
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. | |
Vec * | Vec_ctorPoint (Vec *v, const char *name, int sep, int numR) |
Create a Vec pointer to the tail of a Vec. | |
void | Vec_dtorPoint (Vec **thee) |
The vec pointer destructor. | |
Bmat * | Bmat_ctorPoint (Bmat *P, const char *name, int *ibase, int *numR) |
Create a Bmat pointer to the tail of a Bmat. | |
void | Bmat_dtorPoint (Bmat **thee) |
The Bmat pointer destructor. | |
Mat * | Mat_ctorPoint (Mat *P, const char *name, int ibase, int numR) |
Create a Mat pointer to the tail of a Pro. | |
void | Mat_dtorPoint (Mat **thee) |
The Mat pointer destructor. | |
void | HBmat_printSp (HBmat *thee, char *fname, int pflag) |
Print an HBmat in MATLAB sparse form. | |
void | Bmat_printDiriSp (Bmat *thee, char *fname, int pflag) |
Print the Dirichlet info of a Bmat in MATLAB sparse form. | |
void | Mat_printDiriSp (Mat *thee, char *fname, int pflag) |
Print the Dirichlet info of a Bmat in MATLAB sparse form. | |
Mat * | Bmat_AD (Bmat *thee, int p, int q) |
Return the Mat structure AD. | |
void | HBvec_matvec (HBvec *thee, HBmat *Gmat, int key, Bvec *work) |
Special matrix-vector multiplication for HB. | |
double | Bmat_lnDet (Bmat *thee) |
Calculate the log of the determinant of a SPD Bmat matrix. | |
void | BXLN_copyBmat (Bmat *thee, Bmat *Amat) |
Create an BXLN copy of a Bmat. | |
void | Bmat_copyBXLN (Bmat *thee, Bmat *Alink) |
Produces a Bmat copy of a BXLN. | |
void | BXLN_hbTriple (Bmat *thee, HBmat *G) |
Produce the HB sparse triple matrix product:. | |
void | BXLN_copyBlocks (Bmat *A, Bmat *A12, Bmat *A21, Bmat *A22) |
Copies the A12, A21, and A22 blocks of a XLN, storing them in Mats of type COL, ROW, and DRC respectively. | |
void | BXLN_shrinkLogical (Bmat *A, Bmat *A12, Bmat *A21) |
Sets the sizes of the logical blocks inside a BXLN using the number of rows and cols in two different Bmats. | |
void | HBmat_GHB2WMHB (HBmat *thee, Bmat *M) |
Constructs G12/G22 for WMHB using Mhb and G21. | |
void | XLN_hbTriple (Mat *thee, Mat *GL21, Mat *GL12, Mat *GL22, Mat *GR21, Mat *GR12, Mat *GR22) |
Produce the HB sparse triple matrix product:. | |
void | XLN_matmatContrib (Mat *thee, Mat *Ablock1, int ibase, int flag1, Mat *Ablock2, int jbase, int flag2) |
Multiplies two Mat's contributing the result to an XLN. | |
void | XLN_copySubblock (Mat *thee, Mat *Amat, int flag) |
Copies the A12, A21, or A22 subblock of a XLN, storing it in a Mat of type COL, ROW, or DRC respectively. | |
void | XLN_copyBlocks (Mat *Ablock, Mat *A12, Mat *A21, Mat *A22) |
Copies the A12, A21, and A22 blocks of a XLN, storing them in Mats of type COL, ROW, and DRC respectively. | |
void | XLN_shrinkLogical (Mat *A, Mat *A12, Mat *A21) |
Sets the sizes of the logical blocks inside a XLN using the number of rows and cols in two different Mats. | |
void | Mat_initGWMHB (Mat *G12, Mat *G22, Mat *G21, Mat *Mblock) |
Constructs G12/G22 for WMHB using Mhb and G21. | |
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. | |
Bchar * | Bchar_ctor (Vmem *vmem, const char *name, int pnumB, int pnumR[MAXV]) |
The block character constructor. | |
void | Bchar_dtor (Bchar **thee) |
The block character destructor. | |
void | Bchar_assem (Bchar *thee, int key, Bmat *Ppro) |
Build a block list of F and C points. | |
void | Bchar_assem2 (Bchar *thee, int key, Bmat *amat, Bmat *Ppro) |
Build a block list of F and C points. |
* * MC = < Manifold Code > * Copyright (C) 1994--2008 Michael Holst * * This library is free software; you can redistribute it and/or * modify it under the terms of the GNU Lesser General Public * License as published by the Free Software Foundation; either * version 2.1 of the License, or (at your option) any later version. * * This library is distributed in the hope that it will be useful, * but WITHOUT ANY WARRANTY; without even the implied warranty of * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU * Lesser General Public License for more details. * * You should have received a copy of the GNU Lesser General Public * License along with this library; if not, write to the Free Software * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA * *
Copies the A12, A21, and A22 blocks of a XLN, storing them in Mats of type COL, ROW, and DRC respectively.
* Notes: The logical position of eacg block is determined from the * sizes of the other blocks. * * [ A11 | A12 ] * [ sepxsep | sepxn ] * A_hb = [----------|---------] where numR = sep + n * [ A21 | A22 ] * [ nxsep | nxn ] * * A12, A21, A22 must be passed w/ NULL_STATE blocks. *
A | Pointer to the Bmat object | |
A12 | Pointer to the Bmat object | |
A21 | Pointer to the Bmat object | |
A22 | Pointer to the Bmat object |
Create an BXLN copy of a Bmat.
thee | Pointer to the Bmat object | |
Amat | Pointer to the Bmat object |
Produce the HB sparse triple matrix product:.
* A(p,q)_hb = (I + G(p)') A(p,q) (I + G(q)) * * = [ I Gp21' ] [ Apq11 Apq12 ] [ I Gq12 ] * [ Gp12' (I+Gp22') ] [ Apq21 Apq22 ] [ Gq21 (I+Gq22) ] *
* Notes: The G21, G22, and G12 blocks must be ROW, DRC, and COL * formats respectively. Any NULL blocks are treated as * zero matrices. If the subblock of the BXLN is of symmetric * type, only the upper triangular portion is used/modified. * * The matrix A is overwritten with the A11 result. * * Method: CASE 1: Heirarchical Basis (HB) * * [A11 A12] += [ A12 Gq21 + Gp21'(A21 + A22 Gq21); Gp21'A22 ] * [A21 A22] [ A22 Gq21 ; 0 ] * * CASE 2: Wavelet Modified Heirarchical Basis (WMHB) * * [A11 A12] += [ Gp21'(A21+A22*Gq21) + A12*Gq21 ; 0 ] * [A21 A22] [ Gp22'(A21+A22*Gq21) + A22*Gq21 ; 0 ] * * + [ 0 ; Gp21'(A22+A22*Gq22+A21*Gq12) ] * [ Gp21'(A11+A12*Gq21) ; Gp22'(A22+A22*Gq22+A21*Gq12) ] * * + [ 0 ; A11*Gq12 + A12*Gq22 ] * [ 0 ; A21*Gq12 + A22*Gq22 + Gp12'(A12+A12*Gq22+A11*Gq12) ] *
thee | Pointer to the Bmat object | |
G | Pointer to the HBmat object |
Sets the sizes of the logical blocks inside a BXLN using the number of rows and cols in two different Bmats.
* Notes: This method is used exclusively by HB methods to "replace" * the matrix A with its A11 subblock. The size of the A11 * sublock can be uniquely determined using the sizes of the A12 * and A21 sublocks: * * [ A11 | A12 ] * [ sepxsep | sepxn ] * A = [----------|---------] where numR = sep + n * [ A21 | ] * [ nxsep | ] * * In the future, this routine may actually shrink the matrix, * but for now it only sets the internal numR and numC arrays * to reflect the size of the logical A11 block. *
A | Pointer to the Bmat object | |
A12 | Pointer to the Bmat object | |
A21 | Pointer to the Bmat object |
Copies the A12, A21, and A22 blocks of a XLN, storing them in Mats of type COL, ROW, and DRC respectively.
* Notes: The logical position of each block is determined from the * sizes of the other blocks. * * [ A11 | A12 ] * [ sepxsep | sepxn ] * A_hb = [----------|---------] where numR = sep + n * [ A21 | A22 ] * [ nxsep | nxn ] * * A12, A21, A22 must be passed w/ NULL_STATE (i.e. after ctor). *
Ablock | Pointer to the Mat object | |
A12 | Pointer to the Mat object | |
A21 | Pointer to the Mat object | |
A22 | Pointer to the Mat object |
Copies the A12, A21, or A22 subblock of a XLN, storing it in a Mat of type COL, ROW, or DRC respectively.
* Notes: The logical position of each block is determined from the * sizes of the other blocks. * * [ A11 | A12 ] * [ sepxsep | sepxn ] * A_hb = [----------|---------] where numR = sep + n * [ A21 | A22 ] * [ nxsep | nxn ] * * The subblock must be passed w/ NULL_STATE (i.e. after ctor). *
thee | Pointer to the Mat object | |
Amat | Pointer to the Mat object | |
flag | index for different copy options |
void XLN_hbTriple | ( | Mat * | thee, | |
Mat * | GL21, | |||
Mat * | GL12, | |||
Mat * | GL22, | |||
Mat * | GR21, | |||
Mat * | GR12, | |||
Mat * | GR22 | |||
) |
Produce the HB sparse triple matrix product:.
* A_hb = (I + GL') A (I + GR) * * = [ I GL21' ] [ A11 A12 ] [ I GR12 ] * [ GL12' (I+GL22') ] [ A21 A22 ] [ GR21 (I+GR22) ] *
* Notes: The G21, G22, and G12 blocks must be ROW, DRC, and COL * formats respectively. Any NULL blocks are treated as * zero matrices. If the XLN is of symmetric type, only the * upper triangular portion is used/modified. * * The matrix A is overwritten with the A11 result. * * Method: CASE 1: Heirarchical Basis (HB) * * [A11 A12] += [ A12 GL21 + GL21'(A21 + A22 GR21); GL21'A22 ] * [A21 A22] [ A22 GR21 ; 0 ] * * CASE 2: Wavelet Modified Heirarchical Basis (WMHB) * * [A11 A12] += [ GL21'(A21+A22*GR21) + A12*GR21 ; 0 ] * [A21 A22] [ GL22'(A21+A22*GR21) + A22*GR21 ; 0 ] * * + [ 0 ; GL21'(A22+A22*GR22+A21*GR12) ] * [ GL21'(A11+A12*GR21) ; GL22'(A22+A22*GR22+A21*GR12) ] * * + [ 0 ; A11*GR12 + A12*GR22 ] * [ 0 ; A21*GR12 + A22*GR22 + GL12'(A12+A12*GR22+A11*GR12) ] *
thee | Pointer to the Mat object | |
GL21 | Pointer to the Mat object | |
GL12 | Pointer to the Mat object | |
GL22 | Pointer to the Mat object | |
GR21 | Pointer to the Mat object | |
GR12 | Pointer to the Mat object | |
GR22 | Pointer to the Mat object |
void XLN_matmatContrib | ( | Mat * | thee, | |
Mat * | Ablock1, | |||
int | ibase, | |||
int | flag1, | |||
Mat * | Ablock2, | |||
int | jbase, | |||
int | flag2 | |||
) |
Multiplies two Mat's contributing the result to an XLN.
thee | Pointer to the Mat object | |
Ablock1 | Pointer to the Mat object | |
ibase | pointer to size array | |
flag1 | index for transpose options | |
Ablock2 | Pointer to the Mat object | |
jbase | pointer to size array | |
flag2 | index for transpose options |
Sets the sizes of the logical blocks inside a XLN using the number of rows and cols in two different Mats.
* Notes: This method is used exclusively by HB methods to "replace" * the matrix A with its A11 subblock. The size of the A11 * sublock can be uniquely determined using the sizes of the A12 * and A21 sublocks: * * [ A11 | A12 ] * [ sepxsep | sepxn ] * A = [----------|---------] where numR = sep + n * [ A21 | ] * [ nxsep | ] * * In the future, this routine may actually shrink the matrix, * but for now it only sets the internal numR and numC values * to reflect the size of the logical A11 block. *
A | Pointer to the Mat object | |
A12 | Pointer to the Mat object | |
A21 | Pointer to the Mat object |