Home > Src > Elliptic_2d > Schur > schur_mxv.m

schur_mxv

PURPOSE ^

SCHUR_MXV Computes matrix vector product, where the matrix is the Schur compl.

SYNOPSIS ^

function [v]=schur_mxv(x,AGG,Amm,AGm,LGG,novg,ne);

DESCRIPTION ^

 SCHUR_MXV Computes matrix vector product, where the matrix is the Schur compl.

 Computes matrix vector product v=A*x, where matrix A is the Schur complement 
 matrix

   [v]=schur_mxv(x,AGG,Amm,AGm,LGG,novg,ne);

 Input : x = column array, to perform v=A*x
         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}
         LGG = structure with maps from \Gamma_i to \Gamma
                 set in schur_local.m
         novg = restriction map from (\cup_i \Omega_i) to \Gamma
                 set in cosnovg.m
         ne = number of spectral elements

 Output : v= column array v=A*x, where A is the Schur complement matrix

 References: 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 [v]=schur_mxv(x,AGG,Amm,AGm,LGG,novg,ne);
0002 % SCHUR_MXV Computes matrix vector product, where the matrix is the Schur compl.
0003 %
0004 % Computes matrix vector product v=A*x, where matrix A is the Schur complement
0005 % matrix
0006 %
0007 %   [v]=schur_mxv(x,AGG,Amm,AGm,LGG,novg,ne);
0008 %
0009 % Input : x = column array, to perform v=A*x
0010 %         AGG = structure of local matrices A_{\Gamma_m,\Gamma_m}
0011 %         Amm = structure of local matrices (A_{m,m})^{-1}
0012 %         AGm = structure of local matrices A_{\Gamma_m,m}
0013 %         LGG = structure with maps from \Gamma_i to \Gamma
0014 %                 set in schur_local.m
0015 %         novg = restriction map from (\cup_i \Omega_i) to \Gamma
0016 %                 set in cosnovg.m
0017 %         ne = number of spectral elements
0018 %
0019 % Output : v= column array v=A*x, where A is the Schur complement matrix
0020 %
0021 % References: CHQZ3 = C. Canuto, M.Y. Hussaini, A. Quarteroni, T.A. Zang,
0022 %                    "Spectral Methods. Evolution to Complex Geometries
0023 %                     and Applications to Fluid DynamicsSpectral Methods"
0024 %                    Springer Verlag, Berlin Heidelberg New York, 2007.
0025 
0026 %   Written by Paola Gervasio
0027 %   $Date: 2007/04/01$
0028 
0029 
0030  
0031 n=length(x);
0032 v=zeros(n,1);
0033 
0034 for ie=1:ne
0035 lg=LGG{ie};
0036 ngammam=length(lg);
0037 agg=AGG{ie};agm=AGm{ie};amm=Amm{ie};
0038 xloc=x(novg(lg,ie));
0039 vloc=agg*xloc-agm*(amm*(agm'*xloc));
0040 v(novg(lg,ie))=v(novg(lg,ie))+vloc;
0041 end
0042 
0043 
0044

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