Home > Src > Elliptic_2d > Schur > cosrgam.m

cosrgam

PURPOSE ^

COSRGAM Computes matrix R_Gamma for Schur complement preconditioners

SYNOPSIS ^

function [Rgamma]=cosrgam(novg,LGG,ne,ngamma);

DESCRIPTION ^

 COSRGAM  Computes  matrix R_Gamma for Schur complement preconditioners

   [Rgamma]=cosrgam(novg,LGG,ne,ngamma);

 Computes  matrix R_Gamma of size (ne, ngamma), defined as follows:

 If \Gamma is the interface (uninion of interfaces between subdomains)
 and \Gamma_i :=\partial \Omega_i \cap \Gamma, then
 (R_\Gamma)_ij: =  1/n_j     if  x_j \in \Gamma_i
                =  0         otherwise

   where n_j = is the number of subdomains x_j belongs to

 (see CHQZ3, pag. 399)

 Matrix R_gamma is needed to construct the coarse operator A_H for 
 Balancing Neumann Neumann preconditioner (Schur complement matrix)
 and also to set the unity partition matrix D used in schur_preconn 
 and in preco_bnn

 Input : novg = restriction map from (\cup_i \Omega_i) to \Gamma
                 set in cosnovg.m
         LGG   = structure with maps from \Gamma_i to \Gamma
                 set in schur_local.m
         ne = number of subdomains
         ngamma = number of nodes on Gamma

 Output : Rgamma


 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 [Rgamma]=cosrgam(novg,LGG,ne,ngamma);
0002 % COSRGAM  Computes  matrix R_Gamma for Schur complement preconditioners
0003 %
0004 %   [Rgamma]=cosrgam(novg,LGG,ne,ngamma);
0005 %
0006 % Computes  matrix R_Gamma of size (ne, ngamma), defined as follows:
0007 %
0008 % If \Gamma is the interface (uninion of interfaces between subdomains)
0009 % and \Gamma_i :=\partial \Omega_i \cap \Gamma, then
0010 % (R_\Gamma)_ij: =  1/n_j     if  x_j \in \Gamma_i
0011 %                =  0         otherwise
0012 %
0013 %   where n_j = is the number of subdomains x_j belongs to
0014 %
0015 % (see CHQZ3, pag. 399)
0016 %
0017 % Matrix R_gamma is needed to construct the coarse operator A_H for
0018 % Balancing Neumann Neumann preconditioner (Schur complement matrix)
0019 % and also to set the unity partition matrix D used in schur_preconn
0020 % and in preco_bnn
0021 %
0022 % Input : novg = restriction map from (\cup_i \Omega_i) to \Gamma
0023 %                 set in cosnovg.m
0024 %         LGG   = structure with maps from \Gamma_i to \Gamma
0025 %                 set in schur_local.m
0026 %         ne = number of subdomains
0027 %         ngamma = number of nodes on Gamma
0028 %
0029 % Output : Rgamma
0030 %
0031 %
0032 % References: CHQZ2 = C. Canuto, M.Y. Hussaini, A. Quarteroni, T.A. Zang,
0033 %                    "Spectral Methods. Fundamentals in Single Domains"
0034 %                    Springer Verlag, Berlin Heidelberg New York, 2006.
0035 %             CHQZ3 = C. Canuto, M.Y. Hussaini, A. Quarteroni, T.A. Zang,
0036 %                    "Spectral Methods. Evolution to Complex Geometries
0037 %                     and Applications to Fluid DynamicsSpectral Methods"
0038 %                    Springer Verlag, Berlin Heidelberg New York, 2007.
0039 
0040 %   Written by Paola Gervasio
0041 %   $Date: 2007/04/01$
0042 
0043 
0044 Rgamma=zeros(ne,ngamma);
0045 
0046 for ie=1:ne
0047 lgamma=LGG{ie};
0048 Rgamma(ie,novg(lgamma,ie))=Rgamma(ie,novg(lgamma,ie))+1;
0049 end
0050 
0051 for i=1:ngamma
0052 Rgamma(:,i)=Rgamma(:,i)/sum(Rgamma(:,i));
0053 end
0054 
0055 return

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