Home > Src > Elliptic_2d > Schur > schur_assemb.m

schur_assemb

PURPOSE ^

SCHUR_ASSEMB Assembles global Schur complement matrix

SYNOPSIS ^

function [Sigma]=schur_assemb(AGG,Amm,AGm,LGG,novg,lint,lgamma,param);

DESCRIPTION ^

 SCHUR_ASSEMB Assembles global Schur complement matrix

 Input: 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}
        Lmm = structure containing maps from internal{\Omega_m to \Omega_m
                set in schur_local.m
        LGG = structure with maps from \Gamma_i to \Gamma
                set in schur_local.m
        novi = 2-indexes array of size (max(nvli),ne), computed in cosnovi
        nvli = column array. nvli(ie) is the number of nodes of \Omega_ie
        internal to Omega.
        nov = 2-index array of local to global map, 
                size(nov)=[max(npdx*npdy),ne]
        novg = restriction map from (\cup_i \Omega_i) to \Gamma
                set in cosnovg.m
        lint =  list of those nodes of Omega (close) which
                are internal to Omega.
        lgamma = list of internal nodes which are on the interface between
                spectral elements
        ugamma = u on the interface  Gamma
        f = r.h.s.
        ub = solution on \partial\Omega

 Output: un = solution in Omega


 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 [Sigma]=schur_assemb(AGG,Amm,AGm,LGG,novg,...
0002     lint,lgamma,param);
0003 % SCHUR_ASSEMB Assembles global Schur complement matrix
0004 %
0005 % Input: AGG =structure of local matrices A_{\Gamma_m,\Gamma_m}
0006 %        Amm = structure of local matrices (A_{m,m})^{-1}
0007 %        AGm = structure of local matrices A_{\Gamma_m,m}
0008 %        Lmm = structure containing maps from internal{\Omega_m to \Omega_m
0009 %                set in schur_local.m
0010 %        LGG = structure with maps from \Gamma_i to \Gamma
0011 %                set in schur_local.m
0012 %        novi = 2-indexes array of size (max(nvli),ne), computed in cosnovi
0013 %        nvli = column array. nvli(ie) is the number of nodes of \Omega_ie
0014 %        internal to Omega.
0015 %        nov = 2-index array of local to global map,
0016 %                size(nov)=[max(npdx*npdy),ne]
0017 %        novg = restriction map from (\cup_i \Omega_i) to \Gamma
0018 %                set in cosnovg.m
0019 %        lint =  list of those nodes of Omega (close) which
0020 %                are internal to Omega.
0021 %        lgamma = list of internal nodes which are on the interface between
0022 %                spectral elements
0023 %        ugamma = u on the interface  Gamma
0024 %        f = r.h.s.
0025 %        ub = solution on \partial\Omega
0026 %
0027 % Output: un = solution in Omega
0028 %
0029 %
0030 % References: CHQZ2 = C. Canuto, M.Y. Hussaini, A. Quarteroni, T.A. Zang,
0031 %                    "Spectral Methods. Fundamentals in Single Domains"
0032 %                    Springer Verlag, Berlin Heidelberg New York, 2006.
0033 %             CHQZ3 = C. Canuto, M.Y. Hussaini, A. Quarteroni, T.A. Zang,
0034 %                    "Spectral Methods. Evolution to Complex Geometries
0035 %                     and Applications to Fluid DynamicsSpectral Methods"
0036 %                    Springer Verlag, Berlin Heidelberg New York, 2007.
0037 
0038 %   Written by Paola Gervasio
0039 %   $Date: 2007/04/01$
0040 
0041 
0042 n=length(lgamma);
0043 [ldnov,ne]=size(novg);
0044 Sigma=zeros(n,n);
0045 
0046     % Schur complement
0047 for ie=1:ne
0048 lg=LGG{ie};
0049 ngammam=length(lg);
0050 agg=AGG{ie};agm=AGm{ie};amm=Amm{ie};
0051 Sigma_ie=agg-agm*amm*agm';
0052 Sigma(novg(lg,ie),novg(lg,ie))=Sigma(novg(lg,ie),novg(lg,ie))+Sigma_ie;
0053 end
0054 
0055 if(param(1)==2)
0056     % (Neumann Neumann preconditioner)
0057 for ie=1:ne
0058 lg=LGG{ie};
0059 ngammam=length(lg);
0060 agg=AGG{ie};agm=AGm{ie};amm=Amm{ie};
0061 Sigma_ie=agg-agm*amm*agm';
0062 Sigma(novg(lg,ie),novg(lg,ie))=Sigma(novg(lg,ie),novg(lg,ie))+Sigma_ie;
0063 end
0064 end
0065 
0066 return

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