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_renumber (d06cc)

 Contents

    1  Purpose
    2  Syntax
    7  Accuracy
    9  Example

Purpose

nag_mesh_2d_renumber (d06cc) renumbers the vertices of a given mesh using a Gibbs method, in order the reduce the bandwidth of Finite Element matrices associated with that mesh.

Syntax

[nz, coor, edge, conn, irow, icol, ifail] = d06cc(nnzmax, coor, edge, conn, itrace, 'nv', nv, 'nelt', nelt, 'nedge', nedge)
[nz, coor, edge, conn, irow, icol, ifail] = nag_mesh_2d_renumber(nnzmax, coor, edge, conn, itrace, 'nv', nv, 'nelt', nelt, 'nedge', nedge)

Description

nag_mesh_2d_renumber (d06cc) uses a Gibbs method to renumber the vertices of a given mesh in order to reduce the bandwidth of the associated finite element matrix A. This matrix has elements aij such that:
aij0i​ and ​j​ are vertices belonging to the same triangle.  
This function reduces the bandwidth m, which is the smallest integer such that aij0 whenever i-j>m (see Gibbs et al. (1976) for details about that method). nag_mesh_2d_renumber (d06cc) also returns the sparsity structure of the matrix associated with the renumbered mesh.
This function is derived from material in the MODULEF package from INRIA (Institut National de Recherche en Informatique et Automatique).

References

Gibbs N E, Poole W G Jr and Stockmeyer P K (1976) An algorithm for reducing the bandwidth and profile of a sparse matrix SIAM J. Numer. Anal. 13 236–250

Parameters

Compulsory Input Parameters

1:     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_renumber (d06cc) is called.
Constraint: 4×nelt+nvnnzmaxnv2.
2:     coor2nv – double array
coor1i contains the x coordinate of the ith input mesh vertex, for i=1,2,,nv; while coor2i contains the corresponding y coordinate.
3:     edge3nedge int64int32nag_int array
The specification of the boundary or interface edges. edge1j and edge2j contain the vertex numbers of the two end points of the jth boundary edge. edge3j is a user-supplied tag for the jth boundary or interface edge: edge3j=0 for an interior edge and has a nonzero tag otherwise.
Constraint: 1edgeijnv and edge1jedge2j, for i=1,2 and j=1,2,,nedge.
4:     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.
5:     itrace int64int32nag_int scalar
The level of trace information required from nag_mesh_2d_renumber (d06cc).
itrace0
No output is generated.
itrace=1
Information about the effect of the renumbering on the finite element matrix are output. This information includes the half bandwidth and the sparsity structure of this matrix before and after renumbering.
itrace>1
The output is similar to that produced when itrace=1 but the sparsities (for each row of the matrix, indices of nonzero entries) of the matrix before and after renumbering are also output.

Optional Input Parameters

1:     nv int64int32nag_int scalar
Default: the dimension of the array coor.
The total number of vertices in the input mesh.
Constraint: nv3.
2:     nelt int64int32nag_int scalar
Default: the dimension of the array conn.
The number of triangles in the input mesh.
Constraint: nelt2×nv-1.
3:     nedge int64int32nag_int scalar
Default: the dimension of the array edge.
The number of boundary edges in the input mesh.
Constraint: nedge1.

Output Parameters

1:     nz int64int32nag_int scalar
The number of nonzero entries in the matrix based on the input mesh.
2:     coor2nv – double array
coor1i will contain the x coordinate of the ith renumbered mesh vertex, for i=1,2,,nv; while coor2i will contain the corresponding y coordinate.
3:     edge3nedge int64int32nag_int array
The renumbered specification of the boundary or interface edges.
4:     conn3nelt int64int32nag_int array
The renumbered connectivity of the mesh between triangles and vertices.
5:     irownnzmax int64int32nag_int array
6:     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.
7:     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,
ornedge<1,
ornnzmax<4×nelt+nv or nnzmax>nv2
orconnij<1 or connij>nv for some i=1,2,3 and j=1,2,,nelt,
orconn1j=conn2j or conn1j=conn3j or
conn2j=conn3j for some j=1,2,,nelt,
oredgeij<1 or edgeij>nv for some i=1,2 and j=1,2,,nedge,
oredge1j=edge2j for some j=1,2,,nedge,
orliwork<maxnnzmax,20×nv,
orlrwork<nv.
   ifail=2
A serious error has occurred during the computation of the compact sparsity of the finite element matrix or in an internal call to the renumbering function. Check the input mesh, especially the connectivity between triangles and vertices (the argument conn). 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

In this example, a geometry with two holes (two interior circles inside an exterior one) is considered. The geometry has been meshed using the simple incremental method (nag_mesh_2d_gen_inc (d06aa)) and it has 250 vertices and 402 triangles. The function nag_mesh_2d_gen_boundary (d06ba) is used to renumber the vertices, and one can see the benefit in terms of the sparsity of the finite element matrix based on the renumbered mesh.
function d06cc_example


fprintf('d06cc 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');



d06cc example results


Number of non-zero entries in input mesh:  1556
Number of non-zero entries in output mesh: 1556
d06cc_fig1.png
d06cc_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