Home > Src > Elliptic_2d > Schwarz > stiffq1.m

stiffq1

PURPOSE ^

STIFFQ1 Constructs local stiffness Q1 matrices on extended elements for Schwarz

SYNOPSIS ^

function [Aq1,Abq1,wwq1,linte,ldire,nove]=stiffq1(ifro,nov,xy,nove,nvle);

DESCRIPTION ^

 STIFFQ1 Constructs local stiffness Q1 matrices on extended elements for Schwarz

  [Aq1,Abq1,wwq1,linte,ldire,nove]=stiffq1(ifro,nov,xy,nove,nvle);

 called by schwarz_2d.m and by eig_schwarz_2d.m

 Input: 
        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,
        nov = local -global map, previously generated by cosnov_2d
        xy = column array with global mesh, length: noe=nov(npdx*npdy,ne)
        nove = 2-index array of "extended local" to global map, 
        nvle = 2-index array:
                 nvle (:,1)=number of nodes of the extendend elements
                 nvle (:,2)=number of Q1 elements of the extendend elements
                            along x-direction
                 nvle (:,3)=number of Q1 elements of the extendend elements
                            along y-direction

 Output: Aq1 = cell array with Q1 stiffness matrices (internal/internal nodes) 
               on extended elements
         Abq1 = cell array with Q1 stiffness matrices (internal/boundary nodes) 
               on extended elements
         wwq1 = cell array with Q1 mass matrices (internal nodes) 
                on extended elements
         linte = cell array with list of internal nodes of extended
                 elements
         ldire = cell array with list of boundary nodes of extended
                 elements
         nove = 2-index array of "extended local" to global map, permuted.

 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 [Aq1,Abq1,wwq1,linte,ldire,nove]=stiffq1(ifro,nov,xy,nove,nvle);
0002 % STIFFQ1 Constructs local stiffness Q1 matrices on extended elements for Schwarz
0003 %
0004 %  [Aq1,Abq1,wwq1,linte,ldire,nove]=stiffq1(ifro,nov,xy,nove,nvle);
0005 %
0006 % called by schwarz_2d.m and by eig_schwarz_2d.m
0007 %
0008 % Input:
0009 %        ifro = column array of length noe=nov(npdx*npdy,ne):
0010 %            if (x_i,y_i) is internal to Omega then ifro(i)=0,
0011 %            if (x_i,y_i) is on \partial\Omega then ifro(i)=1,
0012 %        nov = local -global map, previously generated by cosnov_2d
0013 %        xy = column array with global mesh, length: noe=nov(npdx*npdy,ne)
0014 %        nove = 2-index array of "extended local" to global map,
0015 %        nvle = 2-index array:
0016 %                 nvle (:,1)=number of nodes of the extendend elements
0017 %                 nvle (:,2)=number of Q1 elements of the extendend elements
0018 %                            along x-direction
0019 %                 nvle (:,3)=number of Q1 elements of the extendend elements
0020 %                            along y-direction
0021 %
0022 % Output: Aq1 = cell array with Q1 stiffness matrices (internal/internal nodes)
0023 %               on extended elements
0024 %         Abq1 = cell array with Q1 stiffness matrices (internal/boundary nodes)
0025 %               on extended elements
0026 %         wwq1 = cell array with Q1 mass matrices (internal nodes)
0027 %                on extended elements
0028 %         linte = cell array with list of internal nodes of extended
0029 %                 elements
0030 %         ldire = cell array with list of boundary nodes of extended
0031 %                 elements
0032 %         nove = 2-index array of "extended local" to global map, permuted.
0033 %
0034 % References: CHQZ2 = C. Canuto, M.Y. Hussaini, A. Quarteroni, T.A. Zang,
0035 %                    "Spectral Methods. Fundamentals in Single Domains"
0036 %                    Springer Verlag, Berlin Heidelberg New York, 2006.
0037 %             CHQZ3 = C. Canuto, M.Y. Hussaini, A. Quarteroni, T.A. Zang,
0038 %                    "Spectral Methods. Evolution to Complex Geometries
0039 %                     and Applications to Fluid DynamicsSpectral Methods"
0040 %                    Springer Verlag, Berlin Heidelberg New York, 2007.
0041 
0042 %   Written by Paola Gervasio
0043 %   $Date: 2007/04/01$
0044 
0045 [ldnov,ne]=size(nov);
0046 
0047 Aq1=cell(ne,1); Abq1=cell(ne,1);
0048 wwq1=cell(ne,1);
0049 linte=cell(ne,1);
0050 ldire=cell(ne,1);
0051 permu=cell(ne,1);
0052 
0053 % loop on extended spectral elements
0054 for ie=1:ne
0055 mne=nvle(ie,1); nexl=nvle(ie,2);neyl=nvle(ie,3); nel=nexl*neyl;
0056 npxl=nexl+1; npyl=neyl+1;
0057 xyl=zeros(mne,2);
0058 xyl(:,1)=xy(nove(1:mne,ie),1);
0059 xyl(:,2)=xy(nove(1:mne,ie),2);
0060 
0061 % xyl and nove are sorted following lexicographic ordering
0062 
0063 [p]=reorder(xyl,nexl,neyl);
0064 xyl(:,1)=xyl(p,1);xyl(:,2)=xyl(p,2);
0065 nove(1:mne,ie)=nove(p,ie);
0066 %
0067 ifroe=zeros(mne,1);
0068 ifroe(1:npxl)=1;ifroe(1:npxl:mne)=1;ifroe(npxl:npxl:mne)=1;
0069 ifroe(mne-npxl+1:mne)=1;
0070 
0071 permu(ie)={p};
0072 
0073 % Structures for Q1 stiffness matrices
0074 
0075 x=[-1;1]; wx=[1;1]; dx=[-0.5,0.5;-0.5,0.5];
0076 y=x;wy=wx;dy=dx;
0077 iparl=zeros(10,1);
0078 iparl(7)=mne;iparl(6)=nel;iparl(1)=2;iparl(2)=2;iparl(3)=4;
0079 iparl(4)=nexl;
0080 [novl,jaclx,jacly]=meshq1_ie(nove(:,ie),nvle(ie,:),xyl,iparl);
0081 
0082 % stiffness Q1 matrix on the extended element
0083 [Al,wwl]=stiffq1_se(iparl,ifroe,novl,wx,dx,jaclx,wy,dy,jacly);
0084 
0085 % filtering  of both stiffness and mass matrices
0086 
0087 [listaint,listadir]=liste2(ifroe);
0088 Alb=Al(listaint,listadir);
0089 Al=Al(listaint,listaint);
0090 wwl=wwl(listaint);
0091 
0092 Aq1(ie)={chol(Al)}; Abq1(ie)={Alb}; wwq1(ie)={wwl};
0093 linte(ie)={listaint};
0094 ldire(ie)={listadir};
0095 clear Al;
0096 end
0097

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