Home > Src > Elliptic_2d > Schur > schur_precobnn.m

schur_precobnn

PURPOSE ^

SCHUR_PRECOBNN: Solves the sistem (P_B^{NN})^{-1} z=r where P_B^{NN} is the Balancing Neumann-Neumann preconditioner for Schur compl.

SYNOPSIS ^

function [z]=schur_precobnn(r,ne,nvli,novg,D,LGG,Am,Rgamma,PSH,AGG,Amm,AGm);

DESCRIPTION ^

 SCHUR_PRECOBNN:   Solves the sistem (P_B^{NN})^{-1} z=r where P_B^{NN} is the Balancing Neumann-Neumann preconditioner for Schur compl.
    (Algorithm 6.4.4., pag. 399 CHQZ3)
                   Solves the sistem (P_B^{NN})^{-1} z=r

 Input: r= r.h.s, column array of length ngamma
        ne = number of subdomains
        nvli = column array with number of internal nodes of each
        subdomain
        novg = restriction map from (\cup_i \Omega_i) to \Gamma
                 set in cosnovg.m
        D = diagonal weighting matrix (column array of lenght ngamam)
        LGG = structure with maps from \Gamma_i to \Gamma
                 set in schur_local.m
        Am  = structure of choleski factors of local matrices 
              \tilde A_m= A_mm+ M_m/H_m^2 (H_m is the size of the
              subdomain of index m, while M_m si the local mass matrix)
        Rgamma = matrix of size (ne, ngamma) (defined in CHQZ3, pag. 399)
                  (R_\Gamma)_ij: =  1/n_j     if  x_j \in \Gamma_i
                                 =  0         otherwise
        PSH  = pseudoinverse of Sigma_H, i.e. Rgamma^T A_H^+ Rgamma
        AGG = structure of local matrices A_{\Gamma_m,\Gamma_m}
        Amm = structure of local matrices (A_{m,m})^{-1}
        AGm = structure of local matrices A_{\Gamma_m,m}

 Ouput: z = solution, column array of length ngamma


 References: CHQZ2 = C. Canuto, M.Y. Hussaini, A. Quarteroni, T.A. Zang,
                    "Spectral Methods. Fundamentals in Single Domains"
                    Springer Verlag, Berlin Heidelberg New York, 2006.
             CHQZ3 = C. Canuto, M.Y. Hussaini, A. Quarteroni, T.A. Zang,
                    "Spectral Methods. Evolution to Complex Geometries 
                     and Applications to Fluid DynamicsSpectral Methods"
                    Springer Verlag, Berlin Heidelberg New York, 2007.

CROSS-REFERENCE INFORMATION ^

This function calls: This function is called by:

SOURCE CODE ^

0001 function [z]=schur_precobnn(r,ne,nvli,novg,D,LGG,Am,Rgamma,PSH,AGG,Amm,AGm);
0002 % SCHUR_PRECOBNN:   Solves the sistem (P_B^{NN})^{-1} z=r where P_B^{NN} is the Balancing Neumann-Neumann preconditioner for Schur compl.
0003 %    (Algorithm 6.4.4., pag. 399 CHQZ3)
0004 %                   Solves the sistem (P_B^{NN})^{-1} z=r
0005 %
0006 % Input: r= r.h.s, column array of length ngamma
0007 %        ne = number of subdomains
0008 %        nvli = column array with number of internal nodes of each
0009 %        subdomain
0010 %        novg = restriction map from (\cup_i \Omega_i) to \Gamma
0011 %                 set in cosnovg.m
0012 %        D = diagonal weighting matrix (column array of lenght ngamam)
0013 %        LGG = structure with maps from \Gamma_i to \Gamma
0014 %                 set in schur_local.m
0015 %        Am  = structure of choleski factors of local matrices
0016 %              \tilde A_m= A_mm+ M_m/H_m^2 (H_m is the size of the
0017 %              subdomain of index m, while M_m si the local mass matrix)
0018 %        Rgamma = matrix of size (ne, ngamma) (defined in CHQZ3, pag. 399)
0019 %                  (R_\Gamma)_ij: =  1/n_j     if  x_j \in \Gamma_i
0020 %                                 =  0         otherwise
0021 %        PSH  = pseudoinverse of Sigma_H, i.e. Rgamma^T A_H^+ Rgamma
0022 %        AGG = structure of local matrices A_{\Gamma_m,\Gamma_m}
0023 %        Amm = structure of local matrices (A_{m,m})^{-1}
0024 %        AGm = structure of local matrices A_{\Gamma_m,m}
0025 %
0026 % Ouput: z = solution, column array of length ngamma
0027 %
0028 %
0029 % References: CHQZ2 = C. Canuto, M.Y. Hussaini, A. Quarteroni, T.A. Zang,
0030 %                    "Spectral Methods. Fundamentals in Single Domains"
0031 %                    Springer Verlag, Berlin Heidelberg New York, 2006.
0032 %             CHQZ3 = C. Canuto, M.Y. Hussaini, A. Quarteroni, T.A. Zang,
0033 %                    "Spectral Methods. Evolution to Complex Geometries
0034 %                     and Applications to Fluid DynamicsSpectral Methods"
0035 %                    Springer Verlag, Berlin Heidelberg New York, 2007.
0036 
0037 %   Written by Paola Gervasio
0038 %   $Date: 2007/04/01$
0039 
0040 
0041 ra=PSH*r;
0042 
0043 y=r-schur_mxv(ra,AGG,Amm,AGm,LGG,novg,ne);
0044 z=schur_preconnl(y,ne,nvli,novg,D,LGG,Am);
0045 y=schur_mxv(z,AGG,Amm,AGm,LGG,novg,ne);
0046 z=z-PSH*y+ra;
0047 return

Generated on Fri 21-Sep-2007 10:07:00 by m2html © 2003