PINV_SIGMA Computes the pseudoinverse of Sigma_H Computes the pseudoinverse of Sigma_H (Algorithm 6.4.4, pag. 399, CHQZ3) 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} 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 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 Output: PSH = pseudoinverse of Sigma_H, i.e. Rgamma^T A_H^+ Rgamma 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 [PSH]=pinv_sigma(AGG,Amm,AGm,LGG,novg,Rgamma); 0002 % PINV_SIGMA Computes the pseudoinverse of Sigma_H 0003 % 0004 % Computes the pseudoinverse of Sigma_H (Algorithm 6.4.4, pag. 399, CHQZ3) 0005 % 0006 % Input : AGG = structure of local matrices A_{\Gamma_m,\Gamma_m} 0007 % Amm = structure of local matrices (A_{m,m})^{-1} 0008 % AGm = structure of local matrices A_{\Gamma_m,m} 0009 % LGG = structure with maps from \Gamma_i to \Gamma 0010 % set in schur_local.m 0011 % novg = restriction map from (\cup_i \Omega_i) to \Gamma 0012 % set in cosnovg.m 0013 % Rgamma = matrix of size (ne, ngamma) (defined in CHQZ3, pag. 399) 0014 % (R_\Gamma)_ij: = 1/n_j if x_j \in \Gamma_i 0015 % = 0 otherwise 0016 % 0017 % Output: PSH = pseudoinverse of Sigma_H, i.e. Rgamma^T A_H^+ Rgamma 0018 % 0019 % 0020 % References: CHQZ2 = C. Canuto, M.Y. Hussaini, A. Quarteroni, T.A. Zang, 0021 % "Spectral Methods. Fundamentals in Single Domains" 0022 % Springer Verlag, Berlin Heidelberg New York, 2006. 0023 % CHQZ3 = C. Canuto, M.Y. Hussaini, A. Quarteroni, T.A. Zang, 0024 % "Spectral Methods. Evolution to Complex Geometries 0025 % and Applications to Fluid DynamicsSpectral Methods" 0026 % Springer Verlag, Berlin Heidelberg New York, 2007. 0027 0028 % Written by Paola Gervasio 0029 % $Date: 2007/04/01$ 0030 0031 [ne,ngamma]=size(Rgamma); 0032 Sigma=zeros(ngamma,ngamma); 0033 for ie=1:ne 0034 lg=LGG{ie}; 0035 ngammam=length(lg); nl=novg(lg,ie); 0036 agg=AGG{ie};agm=AGm{ie};amm=Amm{ie}; 0037 Sigma_loc=agg-agm*amm*agm'; 0038 Sigma(nl,nl)=Sigma(nl,nl)+Sigma_loc; 0039 end 0040 0041 % A_H=Rgamma*Sigma*Rgamma' 0042 0043 PSH=pinv(full(Rgamma*Sigma*Rgamma')); 0044 clear Sigma; 0045 % Sigma_H^+ 0046 0047 PSH=Rgamma'*PSH*Rgamma;