


SCHUR_PRECONNL: Solves the sistem (P^{NN})^{-1} z=r where P^{NN} is the Neumann-Neumann preconditioner for Schur compl.
Neumann-Neumann preconditioner for Schur
Complement (Algorithm 6.4.3., pag. 398 CHQZ3)
Solves the sistem (P^{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)
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.

0001 function [z]=preconnl(r,ne,nvli,novg,D,LGG,Am); 0002 % SCHUR_PRECONNL: Solves the sistem (P^{NN})^{-1} z=r where P^{NN} is the Neumann-Neumann preconditioner for Schur compl. 0003 % 0004 % Neumann-Neumann preconditioner for Schur 0005 % Complement (Algorithm 6.4.3., pag. 398 CHQZ3) 0006 % Solves the sistem (P^{NN})^{-1} z=r 0007 % 0008 % Input: r= r.h.s, column array of length ngamma 0009 % ne = number of subdomains 0010 % nvli = column array with number of internal nodes of each 0011 % subdomain 0012 % novg = restriction map from (\cup_i \Omega_i) to \Gamma 0013 % set in cosnovg.m 0014 % D = diagonal weighting matrix (column array of lenght ngamam) 0015 % LGG = structure with maps from \Gamma_i to \Gamma 0016 % set in schur_local.m 0017 % Am = structure of choleski factors of local matrices 0018 % \tilde A_m= A_mm+ M_m/H_m^2 (H_m is the size of the 0019 % subdomain of index m, while M_m si the local mass matrix) 0020 % 0021 % Ouput: z = solution, column array of length ngamma 0022 % 0023 % References: CHQZ2 = C. Canuto, M.Y. Hussaini, A. Quarteroni, T.A. Zang, 0024 % "Spectral Methods. Fundamentals in Single Domains" 0025 % Springer Verlag, Berlin Heidelberg New York, 2006. 0026 % CHQZ3 = C. Canuto, M.Y. Hussaini, A. Quarteroni, T.A. Zang, 0027 % "Spectral Methods. Evolution to Complex Geometries 0028 % and Applications to Fluid DynamicsSpectral Methods" 0029 % Springer Verlag, Berlin Heidelberg New York, 2007. 0030 0031 % Written by Paola Gervasio 0032 % $Date: 2007/04/01$ 0033 0034 0035 ngamma=length(r); 0036 0037 z=zeros(ngamma,1); 0038 for ie=1:ne 0039 % 0040 lgamma=LGG{ie}; 0041 mn=nvli(ie); 0042 % compute D_ie (restriction of global weighting matrix D to Gamma_ie) 0043 Die=D(novg(lgamma,ie)); 0044 f=zeros(mn,1); 0045 f(lgamma)=Die.*r(novg(lgamma,ie)); 0046 R=Am{ie}; 0047 zloc=R\(R'\f); 0048 z(novg(lgamma,ie))=z(novg(lgamma,ie))+Die.*zloc(lgamma); 0049 end 0050 return