Home > Src > Elliptic_1d > precofem_1d_se.m

precofem_1d_se

PURPOSE ^

PRECOFEM_1D_SE P1 matrices (stiffness, mass, discrete mass) on macro spectral elements

SYNOPSIS ^

function [AFE,MFE,MFEd]=precofem_1d_se(AFE,MFE,MFEd,ne,xy,nov,nu,beta,gam,param);

DESCRIPTION ^

 PRECOFEM_1D_SE   P1 matrices (stiffness, mass, discrete mass) on macro spectral elements

 [AFE,MFE,MFEd]=precofem_1d_se(AFE,MFE,MFEd,ne,xy,nov,nu,beta,gam,param);

 Input:  AFE = stiffness P1 matrix
         MFE = mass P1 matrix
         MFEd = mass P1 matrix with numerical integration
         ne =  number of P1 element inside spectral element Omega_ie
         xy = mesh on the spectral element Omega_ie
         nov = map mesh Q_n in Omega_ie ---> global mesh in Omega 
         nu   = viscosity (constant>0)
         beta  = coefficient of first order term (constant)
         gam   = coefficient of zeroth order term (constant>=0)

 Output: AFE = stiffness P1 matrix 
         MFE = mass P1 matrix
         MFEd = mass P1 matrix with numerical integration


 Reference: CHQZ2 = C. Canuto, M.Y. Hussaini, A. Quarteroni, T.A. Zang,
                    "Spectral Methods. Fundamentals in Single Domains"
                    Springer Verlag, Berlin Heidelberg New York, 2006.

CROSS-REFERENCE INFORMATION ^

This function calls: This function is called by:

SOURCE CODE ^

0001 function [AFE,MFE,MFEd]=precofem_1d_se(AFE,MFE,MFEd,ne,xy,nov,nu,beta,gam,param);
0002 % PRECOFEM_1D_SE   P1 matrices (stiffness, mass, discrete mass) on macro spectral elements
0003 %
0004 % [AFE,MFE,MFEd]=precofem_1d_se(AFE,MFE,MFEd,ne,xy,nov,nu,beta,gam,param);
0005 %
0006 % Input:  AFE = stiffness P1 matrix
0007 %         MFE = mass P1 matrix
0008 %         MFEd = mass P1 matrix with numerical integration
0009 %         ne =  number of P1 element inside spectral element Omega_ie
0010 %         xy = mesh on the spectral element Omega_ie
0011 %         nov = map mesh Q_n in Omega_ie ---> global mesh in Omega
0012 %         nu   = viscosity (constant>0)
0013 %         beta  = coefficient of first order term (constant)
0014 %         gam   = coefficient of zeroth order term (constant>=0)
0015 %
0016 % Output: AFE = stiffness P1 matrix
0017 %         MFE = mass P1 matrix
0018 %         MFEd = mass P1 matrix with numerical integration
0019 %
0020 %
0021 % Reference: 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 
0025 %   Written by Paola Gervasio
0026 %   $Date: 2007/04/01$
0027 
0028 %    nov_q1 map: mesh_Q1 -----> mesh Q_N in Omega_ie
0029 %    nov_q1_g map: mesh_Q1 -----> global mesh in Omega
0030 %
0031     [nov_p1,nov_p1_g]=nov_p1_1d(nov);
0032     
0033 for ie_p1=1:ne
0034     jacx=5.d-1*(xy(nov_p1(2,ie_p1))-xy(nov_p1(1,ie_p1)));
0035     [Al,Ml,Mld]=matricesp1_1d(nu,beta,gam,jacx,param);
0036     AFE(nov_p1_g(1:2,ie_p1),nov_p1_g(1:2,ie_p1))=AFE(nov_p1_g(1:2,ie_p1),nov_p1_g(1:2,ie_p1))+Al;
0037     MFE(nov_p1_g(1:2,ie_p1),nov_p1_g(1:2,ie_p1))=MFE(nov_p1_g(1:2,ie_p1),nov_p1_g(1:2,ie_p1))+Ml;
0038     MFEd(nov_p1_g(1:2,ie_p1),nov_p1_g(1:2,ie_p1))=MFEd(nov_p1_g(1:2,ie_p1),nov_p1_g(1:2,ie_p1))+Mld;
0039 end
0040 return

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