function [ a, f ] = assemble_poisson_sparse ( node_num, node_xy, ... element_num, element_node, quad_num, nz_num ) %*****************************************************************************80 % %% ASSEMBLE_POISSON_SPARSE assembles the system for the Poisson equation. % % Discussion: % % The matrix is stored in MATLAB sparse matrix format. % % Licensing: % % This code is distributed under the GNU LGPL license. % % Modified: % % 14 February 2003 % % Author: % % John Burkardt % % Parameters: % % Input, integer NODE_NUM, the number of nodes. % % Input, real NODE_XY(2,NODE_NUM), the % coordinates of nodes. % % Input, integer ELEMENT_NUM, the number of elements. % % Input, integer ELEMENT_NODE(3,ELEMENT_NUM); % element_NODE(I,J) is the global index of local node I in element J. % % Input, integer QUAD_NUM, the number of quadrature points used in assembly. % % Input, integer NZ_NUM, the (maximum) number of nonzeros in the matrix. % If set to 0 on input, we hope MATLAB's sparse utility will be able % to take over the task of reallocating space as necessary. % % Output, real sparse A, the coefficient matrix, stored in MATLAB % sparse matrix format. % % Output, real F(NODE_NUM), the right hand side. % % Local parameters: % % Local, real BI, DBIDX, DBIDY, the value of some basis function % and its first derivatives at a quadrature point. % % Local, real BJ, DBJDX, DBJDY, the value of another basis % function and its first derivatives at a quadrature point. % % % Initialize the arrays to zero. % f(1:node_num) = 0.0; fprintf ( 1, '\n' ); fprintf ( 1, 'ASSEMBLE_POISSON_SPARSE:\n' ); fprintf ( 1, ' Setting up sparse Poisson matrix with NZ_NUM = %d\n', nz_num ); a = sparse ( [], [], [], node_num, node_num, nz_num ); % % Get the quadrature weights and nodes. % [ quad_w, quad_xy ] = quad_rule ( quad_num ); % % Add up all quantities associated with the element-th element. % for element = 1 : element_num % % Make a copy of the element. % t3(1:2,1:3) = node_xy(1:2,element_node(1:3,element)); % % Map the quadrature points QUAD_XY to points XY in the physical element. % xy(1:2,1:quad_num) = reference_to_physical_t3 ( t3, quad_num, quad_xy ); area = abs ( triangle_area_2d ( t3 ) ); w(1:quad_num) = quad_w(1:quad_num) * area; quad_rhs = rhs ( quad_num, quad_xy ); quad_k = k_coef ( quad_num, quad_xy ); % % Consider the QUAD-th quadrature point.. % for quad = 1 : quad_num % % Consider the TEST-th test function. % % We generate an integral for every node associated with an unknown. % But if a node is associated with a boundary condition, we do nothing. % for test = 1 : 3 i = element_node(test,element); [ bi, dbidx, dbidy ] = basis_11_t3 ( t3, test, xy(1:2,quad) ); f(i) = f(i) + w(quad) * quad_rhs(quad) * bi; % % Consider the BASIS-th basis function, which is used to form the % value of the solution function. % for basis = 1 : 3 j = element_node(basis,element); [ bj, dbjdx, dbjdy ] = basis_11_t3 ( t3, basis, xy(1:2,quad) ); a(i,j) = a(i,j) + w(quad) * ( ... dbidx * dbjdx + dbidy * dbjdy + quad_k(quad) * bj * bi ); end end end end return end