SCHUR_PCG Preconditioned conjugate gradient to solve the Schur complement system Preconditioned conjugate gradient to solve the Schur complement system [x,iter,res]=schur_pcg(x0, b, tol,maxit,param,... AGG,Amm,AGm,LGG,Am,nvl,novg,D,Rgamma,PSH) Input: x0 = column array for intial guess b = r.h.s for PCG tol = tolerance for PCG stopping test maxit = maximum number of iterations for PCG stopping test param = array of parameters AGG, Amm, AGm, LGG, Am = cells produced in schur_local.m. They contain local matrices and lists nvli = column array. nvli(ie) is the number of nodes of \Omega_ie internal to Omega. novg(i,ie)= 0 if node x_i of Omega_ie does not belong to Gamma = j if node x_i of Omega_ie is the node j of Gamma D = diagonal weighting matrix relative to the interface unknowns, set in partition.m Rgamma = matrix of size (ne, ngamma) (defined in CHQZ3, pag. 399) set in cosrgam.m PSH = pseudoinverse of Sigma_H, set in pinv_sigma.m Output: x = column array of PCG solution iter = number of iterations required to satisfy stopping test ||r^(k)||/||b|| < tol res = value of ||r^(k)||/||b|| at last iteration. 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 [x,iter,res]=schur_pcg(x0, b, tol, maxit,param,... 0002 AGG,Amm,AGm,LGG,Am,nvli,novg,D,Rgamma,PSH) 0003 % SCHUR_PCG Preconditioned conjugate gradient to solve the Schur complement system 0004 % Preconditioned conjugate gradient to solve the Schur complement system 0005 % 0006 % [x,iter,res]=schur_pcg(x0, b, tol,maxit,param,... 0007 % AGG,Amm,AGm,LGG,Am,nvl,novg,D,Rgamma,PSH) 0008 % 0009 % Input: x0 = column array for intial guess 0010 % b = r.h.s for PCG 0011 % tol = tolerance for PCG stopping test 0012 % maxit = maximum number of iterations for PCG stopping test 0013 % param = array of parameters 0014 % AGG, Amm, AGm, LGG, Am = cells produced in schur_local.m. 0015 % They contain local matrices and lists 0016 % nvli = column array. nvli(ie) is the number of nodes of \Omega_ie 0017 % internal to Omega. 0018 % novg(i,ie)= 0 if node x_i of Omega_ie does not belong to Gamma 0019 % = j if node x_i of Omega_ie is the node j of Gamma 0020 % D = diagonal weighting matrix relative to the interface 0021 % unknowns, set in partition.m 0022 % Rgamma = matrix of size (ne, ngamma) (defined in CHQZ3, pag. 399) 0023 % set in cosrgam.m 0024 % PSH = pseudoinverse of Sigma_H, set in pinv_sigma.m 0025 % 0026 % Output: x = column array of PCG solution 0027 % iter = number of iterations required to satisfy stopping test 0028 % ||r^(k)||/||b|| < tol 0029 % res = value of ||r^(k)||/||b|| at last iteration. 0030 % 0031 % References: CHQZ2 = C. Canuto, M.Y. Hussaini, A. Quarteroni, T.A. Zang, 0032 % "Spectral Methods. Fundamentals in Single Domains" 0033 % Springer Verlag, Berlin Heidelberg New York, 2006. 0034 % CHQZ3 = C. Canuto, M.Y. Hussaini, A. Quarteroni, T.A. Zang, 0035 % "Spectral Methods. Evolution to Complex Geometries 0036 % and Applications to Fluid DynamicsSpectral Methods" 0037 % Springer Verlag, Berlin Heidelberg New York, 2007. 0038 0039 % Written by Paola Gervasio 0040 % $Date: 2007/04/01$ 0041 0042 ne=length(nvli); 0043 x=x0; 0044 bb=norm(b); 0045 res=[]; 0046 v=schur_mxv(x,AGG,Amm,AGm,LGG,novg,ne); 0047 r=b-v; 0048 if(param(1)==1) 0049 z=r; 0050 elseif(param(1)==2) 0051 z=schur_preconnl(r,ne,nvli,novg,D,LGG,Am); 0052 elseif(param(1)==3) 0053 z=schur_precobnn(r,ne,nvli,novg,D,LGG,Am,Rgamma,PSH,AGG,Amm,AGm); 0054 end 0055 rr0=z'*r; 0056 normar0=sqrt(r'*r); 0057 p=z; 0058 err=normar0/bb; iter=0; 0059 0060 0061 while err > tol & iter < maxit 0062 0063 %v=A*p; 0064 v=schur_mxv(p,AGG,Amm,AGm,LGG,novg,ne); 0065 alpha=(r'*z)/(p'*v); 0066 x=x+alpha*p; 0067 r=r-alpha*v; 0068 if(param(1)==1) 0069 z=r; 0070 elseif(param(1)==2) 0071 z=schur_preconnl(r,ne,nvli,novg,D,LGG,Am); 0072 elseif(param(1)==3) 0073 z=schur_precobnn(r,ne,nvli,novg,D,LGG,Am,Rgamma,PSH,AGG,Amm,AGm); 0074 end 0075 rr=z'*r; 0076 beta=rr/rr0; 0077 p=z+beta*p; 0078 err=sqrt(r'*r)/normar0; 0079 %res=[res;err]; 0080 iter=iter+1; 0081 rr0=rr; 0082 0083 0084 %fprintf('It. GC %d, err %13.6e \n',iter,err) 0085 end 0086 res=err; 0087 if err>tol & iter>=maxit 0088 fprintf('No convergence is reached in %d iterations \n',iter) 0089 end 0090 0091 return