hide long namesshow long names
hide short namesshow short names
Integer type:  int32  int64  nag_int  show int32  show int32  show int64  show int64  show nag_int  show nag_int

PDF version (NAG web site, 64-bit version, 64-bit version)
Chapter Contents
Chapter Introduction
NAG Toolbox

NAG Toolbox: nag_mesh_2d_sparsity (d06cb)

 Contents

    1  Purpose
    2  Syntax
    7  Accuracy
    9  Example

Purpose

nag_mesh_2d_sparsity (d06cb) generates the sparsity pattern of a finite element matrix associated with a given mesh.

Syntax

[nz, irow, icol, ifail] = d06cb(nv, nnzmax, conn, 'nelt', nelt)
[nz, irow, icol, ifail] = nag_mesh_2d_sparsity(nv, nnzmax, conn, 'nelt', nelt)

Description

nag_mesh_2d_sparsity (d06cb) generates the sparsity pattern of a finite element matrix associated with a given mesh. The sparsity pattern is returned in a coordinate storage format consistent with the sparse linear algebra functions in Chapter F11. More precisely nag_mesh_2d_sparsity (d06cb) returns the number of nonzero elements in the associated sparse matrix, and their row and column indices. This is designed to assist you in applying finite element discretization to meshes from the D06 Chapter Introduction and in solving the resulting sparse linear system using functions from Chapter F11.
The output sparsity pattern is based on the fact that finite element matrix A has elements aij satisfying:
aij0 i​ and ​j ​ are vertices belonging to the same triangle.  

References

None.

Parameters

Compulsory Input Parameters

1:     nv int64int32nag_int scalar
The total number of vertices in the input mesh.
Constraint: nv3.
2:     nnzmax int64int32nag_int scalar
The maximum number of nonzero entries in the matrix based on the input mesh. It is the dimension of the arrays irow and icol as declared in the function from which nag_mesh_2d_sparsity (d06cb) is called.
Constraint: 4×nelt+nvnnzmaxnv2.
3:     conn3nelt int64int32nag_int array
The connectivity of the mesh between triangles and vertices. For each triangle j, connij gives the indices of its three vertices (in anticlockwise order), for i=1,2,3 and j=1,2,,nelt.
Constraint: 1connijnv and conn1jconn2j and conn1jconn3j and conn2jconn3j, for i=1,2,3 and j=1,2,,nelt.

Optional Input Parameters

1:     nelt int64int32nag_int scalar
Default: the dimension of the array conn.
The number of triangles in the input mesh.
Constraint: nelt2×nv-1.

Output Parameters

1:     nz int64int32nag_int scalar
The number of nonzero entries in the matrix associated with the input mesh.
2:     irownnzmax int64int32nag_int array
3:     icolnnzmax int64int32nag_int array
The first nz elements contain the row and column indices of the nonzero elements supplied in the finite element matrix A.
4:     ifail int64int32nag_int scalar
ifail=0 unless the function detects an error (see Error Indicators and Warnings).

Error Indicators and Warnings

Errors or warnings detected by the function:
   ifail=1
On entry,nv<3,
ornelt>2×nv-1,
ornnzmax<4×nelt+nv or nnzmax>nv2
orconnij<1 or connij>nv for some i=1,3 and j, 1jnelt,
orconn1j=conn2j or conn1j=conn3j or
conn2j=conn3j for some j=1,2,,nelt.
   ifail=2
A serious error has occurred in an internal call to an auxiliary function. Check the input mesh, especially the connectivity between triangles and vertices (the argument conn). Array dimensions should be checked as well. If the problem persists, contact NAG.
   ifail=-99
An unexpected error has been triggered by this routine. Please contact NAG.
   ifail=-399
Your licence key may have expired or may not have been installed correctly.
   ifail=-999
Dynamic memory allocation failed.

Accuracy

Not applicable.

Further Comments

None.

Example

See Example in nag_mesh_2d_renumber (d06cc).
function d06cb_example


fprintf('d06cb example results\n\n');


edge = zeros(3, 100, 'int64');
coor = zeros(2, 250);

% Define boundaries
ncirc     = 3; % 3 circles
nvertices = [40, 30, 30];
radii     = [1, 0.49, 0.15];
centres   = [0, 0; -0.5, 0; -0.5, 0.65];

% First circle is outer circle
csign = 1;
i1 = 0;
nvb = 0;
for icirc = 1:ncirc
   for i = 0:nvertices(icirc)-1
      i1 = i1+1;
      theta = 2*pi*i/nvertices(icirc);
      coor(1,i1) = radii(icirc)*cos(theta) + centres(icirc, 1);
      coor(2,i1) = csign*radii(icirc)*sin(theta) +  centres(icirc, 2);
      edge(1,i1) = i1;
      edge(2,i1) = i1 + 1;
      edge(3,i1) = 1;
   end
   edge(2,i1) = nvb + 1;
   nvb = nvb + nvertices(icirc);
   % Subsequent circles are inner circles
   csign = -1;
end
nedge = nvb;

% Initialise mesh control parameters
bspace = zeros(1, 100);
bspace(1:nvb) = 0.05;
smooth = true;
itrace = int64(0);

nnzmax = int64(3000);

% Mesh geometry
[nv, nelt, coor, conn, ifail] = d06aa(edge, coor, bspace, smooth, itrace);

% Compute the sparsity of the FE matrix from the input geometry
[nz, irow, icol, ifail] = d06cb(nv, nnzmax, conn, 'nelt', nelt);

fprintf('\nNumber of non-zero entries in input mesh:  %d\n', nz);

% Plot sparsity of input mesh
fig1 = figure;
plot(irow(1:double(nz)), icol(1:double(nz)), '.');
title ('Input Mesh');
set(gca,'YDir','reverse');


% Call the renumbering routine and get the new sparsity
[nz, coor, edge, conn, irow, icol, ifail] = ...
    d06cc(nnzmax, coor, edge, conn, itrace, 'nelt', nelt);

fprintf('Number of non-zero entries in output mesh: %d\n', nz);
% Plot smoothed mesh
fig2 = figure;
plot(irow(1:double(nz)), icol(1:double(nz)), '.');
title ('Output Mesh');
set(gca,'YDir','reverse');


d06cb example results


Number of non-zero entries in input mesh:  1556
Number of non-zero entries in output mesh: 1556
d06cb_fig1.png
d06cb_fig2.png

PDF version (NAG web site, 64-bit version, 64-bit version)
Chapter Contents
Chapter Introduction
NAG Toolbox

© The Numerical Algorithms Group Ltd, Oxford, UK. 2009–2015