HBmat class

HBmat class. More...


Files

file  whb.h
 Class Whb: stabilized hierarchical basis library.

Classes

struct  sHBmat
 struct HBmat Definition More...

Typedefs

typedef struct sHBmat HBmat
 Declaration of the HBmat class as the HBmat structure.

Functions

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 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.
HBmatHBmat_ctor (Vmem *vmem, const char *bname, int pnumB)
 The HBmat constructor.
void HBmat_dtor (HBmat **thee)
 The HBmat destructor.
void HBmat_printSp (HBmat *thee, char *fname, int pflag)
 Print an HBmat in MATLAB sparse form.
void HBmat_GHB2WMHB (HBmat *thee, Bmat *M)
 Constructs G12/G22 for WMHB using Mhb and G21.

Detailed Description

HBmat class.


Typedef Documentation

typedef struct sHBmat HBmat

Declaration of the HBmat class as the HBmat structure.

Author:
Michael Holst


Function Documentation

HBmat* HBmat_ctor ( Vmem *  vmem,
const char *  bname,
int  pnumB 
)

The HBmat constructor.

Author:
Stephen Bond 2002/01/07
Note:
Class Whb: Non-Inlineable methods (hbtools.c)
Returns:
Pointer to the HBmat object
Parameters:
vmem Memory management object
bname character string name for the matrix
pnumB num vector blocks

void HBmat_dtor ( HBmat **  thee  ) 

The HBmat destructor.

Author:
Stephen Bond 2002/01/07
Note:
Class Whb: Non-Inlineable methods (hbtools.c)
Returns:
None
Parameters:
thee Pointer to the HBmat object

void HBmat_GHB2WMHB ( HBmat thee,
Bmat M 
)

Constructs G12/G22 for WMHB using Mhb and G21.

Author:
Stephen Bond 2002/02/03
Note:
Class Whb: Non-Inlineable methods (brcmat.c)
 *   Notes:  This routine is the second step in the construction of the
 *           change of basis matrix for WMHB.  It constructs the G12 and
 *           G22 blocks using the mass matrix in the HB basis:
 *
 *           [ G11 ] = [ 0 ]  ;   [ G12 ] =  [ I ] inv(-Mhb_11) Mhb_12
 *           [ G21 ]   [ R ]      [ G22 ]    [ R ]
 *
 *           where Mhb is the previously formed mass matrix in the HB basis:
 *
 *           Mhb = [ Mhb_11 Mhb_12 ] = [ I R'] [ M_11 M_12 ] [ I 0 ]
 *                 [ Mhb_21 Mhb_22 ]   [ 0 I ] [ M_21 M_22 ] [ R I ]
 *
 *           Here inv( . ) represents any approximation to the inverse.
 *           To keep optimal complexity this should be chosen very carefully.
 *           Some possible choices include special matrix polynomials,
 *           mass lumping, or simply the inverse of the diagonal.
 *
 *           The G12 and G22 blocks should be passed in with NULL_STATE
 *           subblocks, with only their shape determined (i.e. just after
 *           ctor).  The G21 block must already exist, and contain the tail
 *           of the prolongation matrix, i.e. P = [I ; R].
 * 
Returns:
None
Parameters:
thee Pointer to the HBmat object
M Pointer to the Bmat object

void HBmat_initA ( HBmat Ahb,
HBmat Ghb,
Bmat Alink 
)

Assembles the hierarchical matrix for the HB Vcycle.

Author:
Burak Aksoylu and Stephen Bond
Note:
Class Whb: Non-Inlineable methods (whb.c)
 *   Notes:  This routine takes a pointer to an HBmat array, and
 *           returns a linked (hierarchical) chain of fully assembled
 *           HBmat's.  Each of the links contains several block matrices,
 *           which correspond to the stabilized subsections of the global
 *           block (stiffness) matrix.  All together, these stabilized
 *           blocks form the entire stiffness matrix, with respect to the
 *           recursively defined "hierarchical basis".
 *
 *           Basic Inputs:
 *
 *           HBmat *G    change of basis matrices for each level
 *           BXLN *A   nodal stiffness matrix in doubly linked format
 *           int minLev  coarsest level in the hierarchy
 *
 *           Basic Outputs (inside the HBmat chain):
 *
 *           Bmat *A21, *A22, *A12   stabilized Mat blocks on each level
 * 
Returns:
None
Parameters:
Ahb nodal stiffness matrix in doubly linked format
Ghb change of basis matrices for each level
Alink coarsest level in the hierarchy

void HBmat_initG ( HBmat thee,
Bmat Ppro,
Bmat Mlink,
int  meth 
)

Assembles the change of basis matrices for the HB Vcycle.

Author:
Burak Aksoylu and Stephen Bond
Note:
Class Whb: Non-Inlineable methods (whb.c)
 *   Notes:  This routine takes a pointer to an HBmat array, and
 *           returns a linked (hierarchical) chain of fully assembled
 *           HBmat's.  Each of the links contains several block matrices,
 *           which correspond to the logical subsections of the change
 *           of basis matrix:
 *
 *           [ S11 S12 ] = [  I  0  ] + [ G11 G12 ]
 *           [ S21 S22 ]   [  0  I  ]   [ G21 G22 ]
 *
 *           To keep optimal storage complexity, only G is stored, and
 *           any blocks which are entirely zero have VNULL pointers.
 *
 *           Basic Inputs:
 *
 *           BXLN *M   nodal mass matrix in doubly linked format
 *           Bmat *R     tails of the prolongation matrices (inside Algs)
 *           int minLev  coarsest level in the hierarchy
 *
 *           Basic Outputs (inside the HBmat chain):
 *
 *           Bmat *G21, *G22, *G12   logical blocks of G = S - I.
 *
 *           The mass matrix can be safely passed as VNULL in the HB
 *           case, since it is only touched in the WMHB case.
 * 
Returns:
None
Parameters:
thee Pointer to the HBmat object
Ppro tails of the prolongation matrices (inside Algs) 3 --> Standard HB
5 --> Wavelet Modified HB
Mlink nodal mass matrix in doubly linked format
meth method choice (0=Standard HB,1=Wavelet Modified HB)

void HBmat_initMulti ( HBmat **  Ahb,
HBmat **  Ghb,
Bmat A,
Bmat M,
Bmat Ppro,
int  meth 
)

Assembles the hierarchical matrix structures for the HB Vcycle.

Authors:
Burak Aksoylu and Stephen Bond
Note:
Class Whb: Non-Inlineable methods (bvechb.c)
 *   Notes:  This routine takes a null pointer to two HBmat pointers and
 *           returns a linked (hierarchical) chain of fully assembled
 *           HBmat's.
 *
 *           Basic Inputs:
 *
 *           Bmat *A             the finest level stiffness matrix
 *           Bmat *P0, *P1, ...  the prologation matrices for each level
 *
 *           Basic Outputs (inside the HBmat chain):
 *
 *           Bmat *A21, *A22, *A12   stabilized Mat blocks on each level
 *           Bmat *A11               stabilized A11 block on coarse level
 *           Bvec *d1, *d2           defect section pointers
 *           Bvec *r1, *r2           residual section pointers
 *
 * Example:  Suppose we have a two-level system, resulting from a
 *           mesh refinement.  Let P be the prologation (restriction)
 *           matrix which prolongates (restricts) a trial solution
 *           from level 0 to 1 (1 to 0).  The stiffness matrix, A,
 *           can be naturally decomposed into four blocks, if the nodes
 *           added during refinement (fine nodes) are logically
 *           "appended" to the end of the solution vector:
 *
 *           A = [ A11 A12 ] ; P = [ I ] ; d = [ d1 ] ; r = [ r1 ]
 *               [ A21 A22 ]       [ R ]       [ d2 ]       [ r2 ]
 *
 *           Hierarchical Basis Multigrid (HB), doesn't use the standard
 *           nodal basis, using instead the "hierarchical basis".  The
 *           HB-basis stiffness matrix is related to the nodal basis
 *           stiffness matrix through a "change-of-basis" transform:
 *
 *           A_hb = S'*A*S ; d = S*d_hb ; r_hb = S'*r ; S = [ I 0 ]
 *                                                          [ R I ]
 *           Written out in blocks, A_hb has the form:
 *
 *           A_hb = [ A11 + A12*R + R'*A21 + R'*A22*R ; A12 + R'*A22 ]
 *                  [ A21 + A22*R                     ; A22          ]
 *
 *           The resulting HBmat and HBvec chains will have two levels,
 *           0 and 1, with a pointer to the chains returned as output.
 *           Level 0 contains the only non-null pointer to A11.
 * 
Returns:
None
Parameters:
Ahb nodal stiffness matrix in doubly linked format
Ghb change of basis matrices for each level
A system matrix
M mass matrix
Ppro tails of the prolongation matrices (inside Algs) 3 --> Standard HB
5 --> Wavelet Modified HB
meth method choice (0=Standard HB,1=Wavelet Modified HB)

void HBmat_killMulti ( HBmat **  Ahb,
HBmat **  Ghb 
)

Destruct the hierarchical structures.

Authors:
Burak Aksoylu and Stephen Bond
Note:
Class Whb: Non-Inlineable methods (bvechb.c)
Returns:
None
Parameters:
Ahb nodal stiffness matrix in doubly linked format
Ghb change of basis matrices for each level

void HBmat_killStructure ( HBmat thee  ) 

Kills a chain of HBmat pointers.

Author:
Burak Aksoylu and Stephen Bond
Note:
Class Whb: Non-Inlineable methods (whb.c)
Returns:
None
Parameters:
thee Pointer to the HBmat object

void HBmat_printSp ( HBmat thee,
char *  fname,
int  pflag 
)

Print an HBmat in MATLAB sparse form.

Author:
Stephen Bond 2002/08/11
Note:
Class Whb: Non-Inlineable methods (hbtools.c)
Returns:
None
Parameters:
thee Pointer to the HBmat object
fname character string name for the HBmatrix
pflag index for write/append


Generated on Mon Aug 9 11:13:45 2010 for MC by  doxygen 1.5.6