Home > Src > Elliptic_2d > Schwarz > matr0t.m

matr0t

PURPOSE ^

MATR0T Constructs matrix R_H^T used in (6.3.21), pag. 373 CHQZ3

SYNOPSIS ^

function r0t=matr0t(nx,ny,xy,nov, novc,noec,lista_coarse);

DESCRIPTION ^

 MATR0T Constructs matrix R_H^T used in (6.3.21), pag. 373 CHQZ3
 
 Construction of matrix R_H^T used in (6.3.21), pag. 373 CHQZ3
 
 The Q1 solution on the coarse grid is interpolated on the fine grid

 called by schwarz_2d.m  and by eig_schwarz_2d.m

 Input: nx = polynomial degree in each element (the same in each element)
               along x-direction
        ny = polynomial degree in each element (the same in each element)
               along y-direction
        xy = 2-indexes array wiht coordinates of 2D LGL mesh
        nov = local -global map, previously generated by cosnov_2d
        novc(4,ne) = it maps each macro-Q1 element on the coarse mesh Q1
        noec = number of nodes of the coarse mesh
        lista_coarse = list of nodes of Omega wich are nodes of teh coarse mesh

 Output: r0t =  matrix R_H^T

 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.

CROSS-REFERENCE INFORMATION ^

This function calls: This function is called by:

SOURCE CODE ^

0001 function r0t=matr0t(nx,ny,xy,nov, novc,noec,lista_coarse);
0002 % MATR0T Constructs matrix R_H^T used in (6.3.21), pag. 373 CHQZ3
0003 %
0004 % Construction of matrix R_H^T used in (6.3.21), pag. 373 CHQZ3
0005 %
0006 % The Q1 solution on the coarse grid is interpolated on the fine grid
0007 %
0008 % called by schwarz_2d.m  and by eig_schwarz_2d.m
0009 %
0010 % Input: nx = polynomial degree in each element (the same in each element)
0011 %               along x-direction
0012 %        ny = polynomial degree in each element (the same in each element)
0013 %               along y-direction
0014 %        xy = 2-indexes array wiht coordinates of 2D LGL mesh
0015 %        nov = local -global map, previously generated by cosnov_2d
0016 %        novc(4,ne) = it maps each macro-Q1 element on the coarse mesh Q1
0017 %        noec = number of nodes of the coarse mesh
0018 %        lista_coarse = list of nodes of Omega wich are nodes of teh coarse mesh
0019 %
0020 % Output: r0t =  matrix R_H^T
0021 %
0022 % References: CHQZ2 = C. Canuto, M.Y. Hussaini, A. Quarteroni, T.A. Zang,
0023 %                    "Spectral Methods. Fundamentals in Single Domains"
0024 %                    Springer Verlag, Berlin Heidelberg New York, 2006.
0025 %             CHQZ3 = C. Canuto, M.Y. Hussaini, A. Quarteroni, T.A. Zang,
0026 %                    "Spectral Methods. Evolution to Complex Geometries
0027 %                     and Applications to Fluid DynamicsSpectral Methods"
0028 %                    Springer Verlag, Berlin Heidelberg New York, 2007.
0029 
0030 %   Written by Paola Gervasio
0031 %   $Date: 2007/04/01$
0032 
0033 noe=length(xy);
0034 npdx=nx+1; npdy=ny+1;
0035 [ldnov,ne]=size(nov);
0036 mn=npdx*npdy;
0037 
0038 x=xwlgl(npdx); y=x;
0039 [px,py]=meshgrid(1-x,1-y); phi1=px.*py;phi1=phi1';phi1=phi1(:); clear px; clear py;
0040 [px,py]=meshgrid(1+x,1-y); phi2=px.*py;phi2=phi2';phi2=phi2(:); clear px; clear py;
0041 [px,py]=meshgrid(1-x,1+y); phi3=px.*py;phi3=phi3';phi3=phi3(:); clear px; clear py;
0042 [px,py]=meshgrid(1+x,1+y); phi4=px.*py;phi4=phi4';phi4=phi4(:); clear px; clear py;
0043 r0t=sparse(noe,noec);
0044 for ie=1:ne
0045     ifine=nov(1:mn,ie);
0046     icoarse=novc(1:4,ie);
0047     r0t(ifine,icoarse)=r0t(ifine,icoarse)+0.25*[phi1,phi2,phi3,phi4];
0048 end
0049     
0050 return

Generated on Fri 21-Sep-2007 10:07:00 by m2html © 2003