HBvec class

HBvec class. More...


Files

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

Classes

struct  sHBvec
 struct HBvec Definition More...
struct  sBchar
 struct Bchar Definition More...

Typedefs

typedef struct sHBvec HBvec
 Declaration of the HBvec class as the HBvec structure.

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.
HBvecHBvec_ctor (Vmem *vmem, const char *bname, int pnumB)
 The HBvec constructor.
void HBvec_dtor (HBvec **thee)
 The HBvec destructor.
void HBvec_matvec (HBvec *thee, HBmat *Gmat, int key, Bvec *work)
 Special matrix-vector multiplication for HB.

Detailed Description

HBvec class.


Typedef Documentation

typedef struct sHBvec HBvec

Declaration of the HBvec class as the HBvec structure.

Author:
Michael Holst


Function Documentation

HBvec* HBvec_ctor ( Vmem *  vmem,
const char *  bname,
int  pnumB 
)

The HBvec constructor.

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

void HBvec_dtor ( HBvec **  thee  ) 

The HBvec destructor.

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

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).

Authors:
Burak Aksoylu and Stephen Bond 2002/01/09
Note:
Class Whb: Non-Inlineable methods (whb.c)
 * Warning:  This version is uses RESIDUAL/ERROR form!  The initial guess is
 *           always zero, and the righthand side is always the residual from
 *           the previous iteration.
 *
 *   Notes:  By using a residual/error form for the Vcycle, the code inside
 *           the Vcycle is much cleaner.  In addition, it's much simpler to
 *           decide when the residual tolerance has been met.  If a full
 *           approximation scheme is used, or if the residual is calculated
 *           "on the fly", the magnitude of the residual is not known until
 *           the very bottom of the Vcycle.  This makes Vcycle code a bit
 *           more complicated.
 *
 * HB crash course:
 *
 * A typical multigrid algorithm can be implemented in two different ways:
 * either as a full approximation scheme (FAS) or as a coarse grid correction
 * scheme (CGCS). CGCS seems to be more popular than FAS. Here, Alg_hbVcycle
 * is implemented as CGCS.  As far as we can tell, they both require the same
 * number of operations.
 *
 * Here is the MATLAB equivalent of HBvec_hbVcyc.
 *
 * function [u]=vhbmg(rhs,level) ;
 *
 * %%% global variables: prolongation and stiffness matrices
 * global P_12 P_23 P_34 P_45 P_56 P_67 P_78 P_89;
 * global A_1 A_2 A_3 A_4 A_5 A_6 A_7 A_8 A_9;
 *
 *   %%% get the stiffness matrix on this level
 *   A = eval(['A_' num2str(level) ]);
 *
 *  if (level == 1)
 *   u = A \ rhs;
 *  else
 *    %%%%%%% recover the dimensions
 *    P = eval(['P_' num2str(level-1) num2str(level)]);
 *    [N,M] = size(P);
 *    n = N-M;
 *    P2 = [zeros(M,n);eye(n,n)];
 *    myzeros = zeros(n,1);
 *
 *    %%% get the level-2 block of the stiffness matrix
 *    A_22 = A(N-n+1:N, N-n+1:N);
 *    %%% grab the fine part of rhs
 *    rhs_2 = P2'*rhs;
 *    %%% pre-smoothing by block symmetric Gauss-Seidel
 *    u_2 = smooth_block(A_22,myzeros,rhs_2);
 *    %%% transform basis
 *    w = P2*u_2;
 *    %%% restriction of the residual for the coarse grid
 *    res = P'*(rhs - A*w);
 *
 *        u = vhbmg(res,level-1);
 *
 *    %%% Prolongate the error
 *    %%% Be careful u is the error not the solution because of CGCR
 *    u = P*u;
 *    %%% update fine-grid pseudo residual. Pseudo residual because of
 *    %%% CGCS style. In fact, pseudoRes_2 = residual_2 + A_22*u_2.
 *    pseudoRes_2 = P2'*(rhs - A*u);
 *    %%% post-smoothing by block symmetric Gauss-Seidel
 *    %%% Rather than u_2 = smooth_block(A_22,myzeros,residual_2), we use
 *    u_2 = smooth_block(A_22,u_2,pseudoRes_2);
 *
 *    %%% embed the updated error to the solution by inclusion
 *    error = P2*u_2;
 *    u = u + error;
 *  end;
 * 
Returns:
None
Parameters:
dd Pointer to the HBvec object
Ahb nodal stiffness matrix in doubly linked format
rr Pointer to the HBvec object
Ghb change of basis matrices for each level
ww Pointer to the HBvec object
key index for different solution methods: 0=(Au=f), 1=(A'u=f)
csolv index for different HB options csolv==0 --> nonsymmetric HB
csolv==1 --> adjoint of nonsymmetric HB
csolv==2 --> symmetric HB

void HBvec_initStructure ( HBvec thee,
HBmat Ahb 
)

Creates a chain of HBvec pointers to a Bvec.

Authors:
Burak Aksoylu and Stephen Bond
Note:
Class Whb: Non-Inlineable methods (whb.c)
 *   Notes:  This routine takes a null pointer to an HBvec array, and
 *           returns a linked (hierarchical) chain of fully assembled
 *           HBvec's.  Each of the links contains two Bvec pointers,
 *           which correspond to logical subsections of a global Bvec.
 *           The fully assemble HBvec array consist of only pointers,
 *           sharing the data space used by the specified Bvec.
 *
 *           Basic Inputs:
 *
 *           Bmat *P0, *P1, ...  the prologation matrices for each level
 *           Bvec *u             a block vector on the finest level
 *                               (passed using the HBvec bv pointer)
 *
 *           Basic Outputs (inside the HBvec chain):
 *
 *           Bvec *bv, *bv2     pointers to logical subsections of u.
 * 
Returns:
None
Parameters:
thee Pointer to the HBvec object
Ahb nodal stiffness matrix in doubly linked format

void HBvec_killStructure ( HBvec thee  ) 

Kills a chain of HBvec pointers.

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

void HBvec_matvec ( HBvec thee,
HBmat Gmat,
int  key,
Bvec work 
)

Special matrix-vector multiplication for HB.

Author:
Michael Holst
Note:
Class Whb: Non-Inlineable methods (hbtools.c)
Blocks which are entirely zero should be specified using NULL pointers. The work vector should be of the same length as the number of rows in G21. This vector is only used if A12 != 0 or A22 != 0, otherwise it can be safely set to VNULL.
Returns:
None
Parameters:
thee Pointer to the HBvec object
Gmat Pointer to the source HBmat
key 0 standard product, u += G *u
1 transpose product, u += G'*u
work Pointer to the block vector


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