MESHQ1_COARSE Construction of structures for the coarse mesh [lista_coarse,novc,novcg,jacx,jacy]=meshq1_coarse(nov,xy,ipar); Construction of structures for the coarse mesh Input: npdx = number of nodes in each original spectral element (along x direction) npdy = number of nodes in each original spectral element (along y direction) nov = local -global map, previously generated by cosnov_2d xy = column array with global mesh, length: noe=nov(npdx*npdy,ne) Output: lista_coarse = list of nodes of Omega wich are nodes of teh coarse mesh novc(4,ne) = it maps each macro-Q1 element on the coarse mesh Q1 novcg(4,ne) = it maps each macro-Q1 element on the fine mesh Qn jacxe, jacye = jacobians of coarse (extended) Q1 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 [lista_coarse,novc,novcg,jacxe,jacye]=meshq1_coarse(npdx,npdy,nov,xy); 0002 % MESHQ1_COARSE Construction of structures for the coarse mesh 0003 % 0004 % [lista_coarse,novc,novcg,jacx,jacy]=meshq1_coarse(nov,xy,ipar); 0005 % 0006 % Construction of structures for the coarse mesh 0007 % Input: 0008 % npdx = number of nodes in each original spectral element (along x 0009 % direction) 0010 % npdy = number of nodes in each original spectral element (along y 0011 % direction) 0012 % nov = local -global map, previously generated by cosnov_2d 0013 % xy = column array with global mesh, length: noe=nov(npdx*npdy,ne) 0014 % 0015 % Output: 0016 % lista_coarse = list of nodes of Omega wich are nodes of teh coarse mesh 0017 % novc(4,ne) = it maps each macro-Q1 element on the coarse mesh Q1 0018 % novcg(4,ne) = it maps each macro-Q1 element on the fine mesh Qn 0019 % jacxe, jacye = jacobians of coarse (extended) Q1 mesh 0020 % 0021 % References: CHQZ2 = C. Canuto, M.Y. Hussaini, A. Quarteroni, T.A. Zang, 0022 % "Spectral Methods. Fundamentals in Single Domains" 0023 % Springer Verlag, Berlin Heidelberg New York, 2006. 0024 % CHQZ3 = C. Canuto, M.Y. Hussaini, A. Quarteroni, T.A. Zang, 0025 % "Spectral Methods. Evolution to Complex Geometries 0026 % and Applications to Fluid DynamicsSpectral Methods" 0027 % Springer Verlag, Berlin Heidelberg New York, 2007. 0028 0029 % Written by Paola Gervasio 0030 % $Date: 2007/04/01$ 0031 0032 [ldnov,ne]=size(nov); 0033 noe=length(xy); 0034 mn=npdx*npdy; 0035 %ipar(7);ne=ipar(6);npdx=ipar(1);npdy=ipar(2);mn=ipar(3); 0036 novc=zeros(4,ne);novcg=zeros(4,ne); 0037 0038 jacxe=zeros(ne,1); jacye=zeros(ne,1); 0039 lista_coarse=[]; 0040 for ie=1:ne 0041 i1=nov(1,ie);i2=nov(npdx,ie);i3=nov(mn-npdx+1,ie);i4=nov(mn,ie); 0042 novcg(1:4,ie)=[i1;i2;i3;i4]; 0043 lista_coarse=[lista_coarse;i1;i2;i3;i4]; 0044 jacxe(ie)=(xy(i2,1)-xy(i1,1))*0.5; 0045 jacye(ie)=(xy(i3,2)-xy(i1,2))*0.5; 0046 end 0047 lista_coarse=unique(lista_coarse); 0048 noec=length(lista_coarse); 0049 for ie=1:ne 0050 for i=1:4 0051 knovcg=novcg(i,ie); 0052 trov=0;j=1; 0053 while j<=noec &trov==0 0054 if knovcg==lista_coarse(j) 0055 novc(i,ie)=j; 0056 trov=1; 0057 end 0058 j=j+1; 0059 end 0060 end 0061 end 0062 0063 0064 0065 0066 0067