STIFFQ1H Construction of stiffness Q1 matrix on the coarse grid. [Ac,Acb,wwc,lista_coarse,noec,novc,lintc,ldirc]=... stiffq1H(nx,nex,ny,ney,xy,nov,ifro); called by schwarz_2d.m and b eig_schwarz_2d.m Input: nx = polynomial degree in each element (the same in each element) along x-direction nex = number of elements (equally spaced) along x-direction ny = polynomial degree in each element (the same in each element) along y-direction% gammax = column or row array of length nex-1. ney = number of elements (equally spaced) along y-direction xy = column array with global mesh, length: noe=nov(npdx*npdy,ne) nov = local -global map, previously generated by cosnov_2d ifro = column array of length noe=nov(npdx*npdy,ne): if (x_i,y_i) is internal to Omega then ifro(i)=0, if (x_i,y_i) is on \partial\Omega then ifro(i)=1, Output: Ac = Cholesky factor of the Q1 stiffness matrix (internal/internal nodes) on the coarse mesh Acb = Q1 stiffness matrix (internal/boundary nodes) on the coarse mesh wwc = Q1 mass matrix (internal nodes) on the coarse mesh lista_coarse = list of nodes of Omega wich are nodes of teh coarse mesh noec = number of nodes of the coarse mesh novc(4,ne) = it maps each macro-Q1 element on the coarse mesh Q1 lintc = lista of internal nodes of the coarse mesh ldirc = lista of boudanry nodes of the coarse mesh 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 [Ac,Acb,wwc,lista_coarse,noec,novc,lintc,ldirc]=... 0002 stiffq1H(nx,nex,ny,ney,xy,nov,ifro); 0003 % STIFFQ1H Construction of stiffness Q1 matrix on the coarse grid. 0004 % 0005 % [Ac,Acb,wwc,lista_coarse,noec,novc,lintc,ldirc]=... 0006 % stiffq1H(nx,nex,ny,ney,xy,nov,ifro); 0007 % 0008 % called by schwarz_2d.m and b eig_schwarz_2d.m 0009 % 0010 % Input: 0011 % nx = polynomial degree in each element (the same in each element) 0012 % along x-direction 0013 % nex = number of elements (equally spaced) along x-direction 0014 % ny = polynomial degree in each element (the same in each element) 0015 % along y-direction% gammax = column or row array of length nex-1. 0016 % ney = number of elements (equally spaced) along y-direction 0017 % xy = column array with global mesh, length: noe=nov(npdx*npdy,ne) 0018 % nov = local -global map, previously generated by cosnov_2d 0019 % ifro = column array of length noe=nov(npdx*npdy,ne): 0020 % if (x_i,y_i) is internal to Omega then ifro(i)=0, 0021 % if (x_i,y_i) is on \partial\Omega then ifro(i)=1, 0022 % 0023 % Output: 0024 % Ac = Cholesky factor of the Q1 stiffness matrix (internal/internal nodes) 0025 % on the coarse mesh 0026 % Acb = Q1 stiffness matrix (internal/boundary nodes) 0027 % on the coarse mesh 0028 % wwc = Q1 mass matrix (internal nodes) 0029 % on the coarse mesh 0030 % lista_coarse = list of nodes of Omega wich are nodes of teh coarse mesh 0031 % noec = number of nodes of the coarse mesh 0032 % novc(4,ne) = it maps each macro-Q1 element on the coarse mesh Q1 0033 % lintc = lista of internal nodes of the coarse mesh 0034 % ldirc = lista of boudanry nodes of the coarse mesh 0035 % 0036 % References: CHQZ2 = C. Canuto, M.Y. Hussaini, A. Quarteroni, T.A. Zang, 0037 % "Spectral Methods. Fundamentals in Single Domains" 0038 % Springer Verlag, Berlin Heidelberg New York, 2006. 0039 % CHQZ3 = C. Canuto, M.Y. Hussaini, A. Quarteroni, T.A. Zang, 0040 % "Spectral Methods. Evolution to Complex Geometries 0041 % and Applications to Fluid DynamicsSpectral Methods" 0042 % Springer Verlag, Berlin Heidelberg New York, 2007. 0043 0044 % Written by Paola Gervasio 0045 % $Date: 2007/04/01$ 0046 0047 npdx=nx+1;npdy=ny+1;ne=nex*ney; 0048 npxc=nex+1;npyc=ney+1;noec=npxc*npyc; 0049 mn=npdx*npdy;nm1=mn-npdx+1; 0050 0051 % costruisco la mesh coarse 0052 [lista_coarse,novc,novcg,jacxe,jacye]=meshq1_coarse(npdx,npdy,nov,xy); 0053 noec=length(lista_coarse); 0054 0055 xyc=zeros(noec,2); xyc(:,1)=xy(lista_coarse,1);xyc(:,2)=xy(lista_coarse,2); 0056 ifroc=zeros(noec,1); ifroc=ifro(lista_coarse); 0057 0058 x=[-1;1]; wx=[1;1]; dx=[-0.5,0.5;-0.5,0.5]; 0059 y=x;wy=wx;dy=dx; 0060 iparl=zeros(10,1); 0061 iparl(7)=noec;iparl(6)=ne;iparl(1)=2;iparl(2)=2;iparl(3)=4; 0062 0063 0064 0065 % costruisco la matrice Q1 del laplaciano sull'elemento esteso. 0066 [Ac,wwc]=stiffq1_se(iparl,ifroc,novc,wx,dx,jacxe,wy,dy,jacye); 0067 0068 % filtro la matrice sui nodi interni e sui nodi di bordo 0069 % e costruisco la fattorizzazione di Cholesky di Al 0070 0071 [ldirc,lintc]=liste1(ifroc); 0072 Acb=Ac(lintc,ldirc); 0073 Ac=Ac(lintc,lintc); 0074 wwc=wwc(lintc); 0075 Ac=chol(Ac); 0076 0077 return